New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
bdydyn3d.F90 in NEMO/branches/2018/dev_r10164_HPC09_ESIWACE_PREP_MERGE/src/OCE/BDY – NEMO

source: NEMO/branches/2018/dev_r10164_HPC09_ESIWACE_PREP_MERGE/src/OCE/BDY/bdydyn3d.F90 @ 10314

Last change on this file since 10314 was 10314, checked in by smasson, 5 years ago

dev_r10164_HPC09_ESIWACE_PREP_MERGE: action 2: add generic glob_min/max/sum and locmin/max, complete timing and report (including bdy and icb), see #2133

  • Property svn:keywords set to Id
File size: 15.6 KB
Line 
1MODULE bdydyn3d
2   !!======================================================================
3   !!                       ***  MODULE  bdydyn3d  ***
4   !! Unstructured Open Boundary Cond. :   Flow relaxation scheme on baroclinic velocities
5   !!======================================================================
6   !! History :  3.4  !  2011     (D. Storkey) new module as part of BDY rewrite
7   !!            3.5  !  2012     (S. Mocavero, I. Epicoco) Optimization of BDY communications
8   !!----------------------------------------------------------------------
9   !!   bdy_dyn3d        : apply open boundary conditions to baroclinic velocities
10   !!   bdy_dyn3d_frs    : apply Flow Relaxation Scheme
11   !!----------------------------------------------------------------------
12   USE timing          ! Timing
13   USE oce             ! ocean dynamics and tracers
14   USE dom_oce         ! ocean space and time domain
15   USE bdy_oce         ! ocean open boundary conditions
16   USE bdylib          ! for orlanski library routines
17   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
18   USE in_out_manager  !
19   Use phycst
20
21   IMPLICIT NONE
22   PRIVATE
23
24   PUBLIC   bdy_dyn3d     ! routine called by bdy_dyn
25   PUBLIC   bdy_dyn3d_dmp ! routine called by step
26
27   !!----------------------------------------------------------------------
28   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
29   !! $Id$
30   !! Software governed by the CeCILL license (see ./LICENSE)
31   !!----------------------------------------------------------------------
32CONTAINS
33
34   SUBROUTINE bdy_dyn3d( kt )
35      !!----------------------------------------------------------------------
36      !!                  ***  SUBROUTINE bdy_dyn3d  ***
37      !!
38      !! ** Purpose : - Apply open boundary conditions for baroclinic velocities
39      !!
40      !!----------------------------------------------------------------------
41      INTEGER, INTENT(in) ::   kt   ! Main time step counter
42      !
43      INTEGER ::   ib_bdy   ! loop index
44      !!----------------------------------------------------------------------
45      !
46      DO ib_bdy=1, nb_bdy
47         !
48         SELECT CASE( cn_dyn3d(ib_bdy) )
49         CASE('none')        ;   CYCLE
50         CASE('frs' )        ;   CALL bdy_dyn3d_frs( idx_bdy(ib_bdy), dta_bdy(ib_bdy), kt, ib_bdy )
51         CASE('specified')   ;   CALL bdy_dyn3d_spe( idx_bdy(ib_bdy), dta_bdy(ib_bdy), kt, ib_bdy )
52         CASE('zero')        ;   CALL bdy_dyn3d_zro( idx_bdy(ib_bdy), dta_bdy(ib_bdy), kt, ib_bdy )
53         CASE('orlanski' )   ;   CALL bdy_dyn3d_orlanski( idx_bdy(ib_bdy), dta_bdy(ib_bdy), ib_bdy, ll_npo=.false. )
54         CASE('orlanski_npo');   CALL bdy_dyn3d_orlanski( idx_bdy(ib_bdy), dta_bdy(ib_bdy), ib_bdy, ll_npo=.true. )
55         CASE('zerograd')    ;   CALL bdy_dyn3d_zgrad( idx_bdy(ib_bdy), dta_bdy(ib_bdy), kt, ib_bdy )
56         CASE('neumann')     ;   CALL bdy_dyn3d_nmn( idx_bdy(ib_bdy), ib_bdy )
57         CASE DEFAULT        ;   CALL ctl_stop( 'bdy_dyn3d : unrecognised option for open boundaries for baroclinic velocities' )
58         END SELECT
59      END DO
60      !
61   END SUBROUTINE bdy_dyn3d
62
63
64   SUBROUTINE bdy_dyn3d_spe( idx, dta, kt , ib_bdy )
65      !!----------------------------------------------------------------------
66      !!                  ***  SUBROUTINE bdy_dyn3d_spe  ***
67      !!
68      !! ** Purpose : - Apply a specified value for baroclinic velocities
69      !!                at open boundaries.
70      !!
71      !!----------------------------------------------------------------------
72      INTEGER        , INTENT(in) ::   kt      ! time step index
73      TYPE(OBC_INDEX), INTENT(in) ::   idx     ! OBC indices
74      TYPE(OBC_DATA) , INTENT(in) ::   dta     ! OBC external data
75      INTEGER        , INTENT(in) ::   ib_bdy  ! BDY set index
76      !
77      INTEGER  ::   jb, jk         ! dummy loop indices
78      INTEGER  ::   ii, ij, igrd   ! local integers
79      REAL(wp) ::   zwgt           ! boundary weight
80      !!----------------------------------------------------------------------
81      !
82      igrd = 2                      ! Relaxation of zonal velocity
83      DO jb = 1, idx%nblenrim(igrd)
84         DO jk = 1, jpkm1
85            ii   = idx%nbi(jb,igrd)
86            ij   = idx%nbj(jb,igrd)
87            ua(ii,ij,jk) = dta%u3d(jb,jk) * umask(ii,ij,jk)
88         END DO
89      END DO
90      !
91      igrd = 3                      ! Relaxation of meridional velocity
92      DO jb = 1, idx%nblenrim(igrd)
93         DO jk = 1, jpkm1
94            ii   = idx%nbi(jb,igrd)
95            ij   = idx%nbj(jb,igrd)
96            va(ii,ij,jk) = dta%v3d(jb,jk) * vmask(ii,ij,jk)
97         END DO
98      END DO
99      CALL lbc_bdy_lnk( 'bdydyn3d', ua, 'U', -1., ib_bdy )   ! Boundary points should be updated 
100      CALL lbc_bdy_lnk( 'bdydyn3d', va, 'V', -1., ib_bdy )   
101      !
102      IF( kt == nit000 )   CLOSE( unit = 102 )
103      !
104   END SUBROUTINE bdy_dyn3d_spe
105
106
107   SUBROUTINE bdy_dyn3d_zgrad( idx, dta, kt , ib_bdy )
108      !!----------------------------------------------------------------------
109      !!                  ***  SUBROUTINE bdy_dyn3d_zgrad  ***
110      !!
111      !! ** Purpose : - Enforce a zero gradient of normal velocity
112      !!
113      !!----------------------------------------------------------------------
114      INTEGER                     ::   kt
115      TYPE(OBC_INDEX), INTENT(in) ::   idx  ! OBC indices
116      TYPE(OBC_DATA),  INTENT(in) ::   dta  ! OBC external data
117      INTEGER,         INTENT(in) ::   ib_bdy  ! BDY set index
118      !!
119      INTEGER  ::   jb, jk         ! dummy loop indices
120      INTEGER  ::   ii, ij, igrd   ! local integers
121      REAL(wp) ::   zwgt           ! boundary weight
122      INTEGER  ::   fu, fv
123      !!----------------------------------------------------------------------
124      !
125      igrd = 2                      ! Copying tangential velocity into bdy points
126      DO jb = 1, idx%nblenrim(igrd)
127         DO jk = 1, jpkm1
128            ii   = idx%nbi(jb,igrd)
129            ij   = idx%nbj(jb,igrd)
130            fu   = ABS( ABS (NINT( idx%flagu(jb,igrd) ) ) - 1 )
131            ua(ii,ij,jk) = ua(ii,ij,jk) * REAL( 1 - fu) + ( ua(ii,ij+fu,jk) * umask(ii,ij+fu,jk) &
132                        &+ ua(ii,ij-fu,jk) * umask(ii,ij-fu,jk) ) * umask(ii,ij,jk) * REAL( fu )
133         END DO
134      END DO
135      !
136      igrd = 3                      ! Copying tangential velocity into bdy points
137      DO jb = 1, idx%nblenrim(igrd)
138         DO jk = 1, jpkm1
139            ii   = idx%nbi(jb,igrd)
140            ij   = idx%nbj(jb,igrd)
141            fv   = ABS( ABS (NINT( idx%flagv(jb,igrd) ) ) - 1 )
142            va(ii,ij,jk) = va(ii,ij,jk) * REAL( 1 - fv ) + ( va(ii+fv,ij,jk) * vmask(ii+fv,ij,jk) &
143                        &+ va(ii-fv,ij,jk) * vmask(ii-fv,ij,jk) ) * vmask(ii,ij,jk) * REAL( fv )
144         END DO
145      END DO
146      CALL lbc_bdy_lnk( 'bdydyn3d', ua, 'U', -1., ib_bdy )   ! Boundary points should be updated 
147      CALL lbc_bdy_lnk( 'bdydyn3d', va, 'V', -1., ib_bdy )   
148      !
149      IF( kt == nit000 )   CLOSE( unit = 102 )
150      !
151   END SUBROUTINE bdy_dyn3d_zgrad
152
153
154   SUBROUTINE bdy_dyn3d_zro( idx, dta, kt, ib_bdy )
155      !!----------------------------------------------------------------------
156      !!                  ***  SUBROUTINE bdy_dyn3d_zro  ***
157      !!
158      !! ** Purpose : - baroclinic velocities = 0. at open boundaries.
159      !!
160      !!----------------------------------------------------------------------
161      INTEGER        , INTENT(in) ::   kt      ! time step index
162      TYPE(OBC_INDEX), INTENT(in) ::   idx     ! OBC indices
163      TYPE(OBC_DATA) , INTENT(in) ::   dta     ! OBC external data
164      INTEGER,         INTENT(in) ::   ib_bdy  ! BDY set index
165      !
166      INTEGER  ::   ib, ik         ! dummy loop indices
167      INTEGER  ::   ii, ij, igrd   ! local integers
168      REAL(wp) ::   zwgt           ! boundary weight
169      !!----------------------------------------------------------------------
170      !
171      igrd = 2                       ! Everything is at T-points here
172      DO ib = 1, idx%nblenrim(igrd)
173         ii = idx%nbi(ib,igrd)
174         ij = idx%nbj(ib,igrd)
175         DO ik = 1, jpkm1
176            ua(ii,ij,ik) = 0._wp
177         END DO
178      END DO
179
180      igrd = 3                       ! Everything is at T-points here
181      DO ib = 1, idx%nblenrim(igrd)
182         ii = idx%nbi(ib,igrd)
183         ij = idx%nbj(ib,igrd)
184         DO ik = 1, jpkm1
185            va(ii,ij,ik) = 0._wp
186         END DO
187      END DO
188      !
189      CALL lbc_bdy_lnk( 'bdydyn3d', ua, 'U', -1., ib_bdy )   ;   CALL lbc_bdy_lnk( 'bdydyn3d', va, 'V', -1.,ib_bdy )   ! Boundary points should be updated
190      !
191      IF( kt == nit000 )   CLOSE( unit = 102 )
192      !
193   END SUBROUTINE bdy_dyn3d_zro
194
195
196   SUBROUTINE bdy_dyn3d_frs( idx, dta, kt, ib_bdy )
197      !!----------------------------------------------------------------------
198      !!                  ***  SUBROUTINE bdy_dyn3d_frs  ***
199      !!
200      !! ** Purpose : - Apply the Flow Relaxation Scheme for baroclinic velocities
201      !!                at open boundaries.
202      !!
203      !! References :- Engedahl H., 1995: Use of the flow relaxation scheme in
204      !!               a three-dimensional baroclinic ocean model with realistic
205      !!               topography. Tellus, 365-382.
206      !!----------------------------------------------------------------------
207      INTEGER        , INTENT(in) ::   kt      ! time step index
208      TYPE(OBC_INDEX), INTENT(in) ::   idx     ! OBC indices
209      TYPE(OBC_DATA) , INTENT(in) ::   dta     ! OBC external data
210      INTEGER,         INTENT(in) ::   ib_bdy  ! BDY set index
211      !
212      INTEGER  ::   jb, jk         ! dummy loop indices
213      INTEGER  ::   ii, ij, igrd   ! local integers
214      REAL(wp) ::   zwgt           ! boundary weight
215      !!----------------------------------------------------------------------
216      !
217      igrd = 2                      ! Relaxation of zonal velocity
218      DO jb = 1, idx%nblen(igrd)
219         DO jk = 1, jpkm1
220            ii   = idx%nbi(jb,igrd)
221            ij   = idx%nbj(jb,igrd)
222            zwgt = idx%nbw(jb,igrd)
223            ua(ii,ij,jk) = ( ua(ii,ij,jk) + zwgt * ( dta%u3d(jb,jk) - ua(ii,ij,jk) ) ) * umask(ii,ij,jk)
224         END DO
225      END DO
226      !
227      igrd = 3                      ! Relaxation of meridional velocity
228      DO jb = 1, idx%nblen(igrd)
229         DO jk = 1, jpkm1
230            ii   = idx%nbi(jb,igrd)
231            ij   = idx%nbj(jb,igrd)
232            zwgt = idx%nbw(jb,igrd)
233            va(ii,ij,jk) = ( va(ii,ij,jk) + zwgt * ( dta%v3d(jb,jk) - va(ii,ij,jk) ) ) * vmask(ii,ij,jk)
234         END DO
235      END DO
236      CALL lbc_bdy_lnk( 'bdydyn3d', ua, 'U', -1., ib_bdy )    ! Boundary points should be updated
237      CALL lbc_bdy_lnk( 'bdydyn3d', va, 'V', -1., ib_bdy )   
238      !
239      IF( kt == nit000 )   CLOSE( unit = 102 )
240      !
241   END SUBROUTINE bdy_dyn3d_frs
242
243
244   SUBROUTINE bdy_dyn3d_orlanski( idx, dta, ib_bdy, ll_npo )
245      !!----------------------------------------------------------------------
246      !!                 ***  SUBROUTINE bdy_dyn3d_orlanski  ***
247      !!             
248      !!              - Apply Orlanski radiation to baroclinic velocities.
249      !!              - Wrapper routine for bdy_orlanski_3d
250      !!
251      !!
252      !! References:  Marchesiello, McWilliams and Shchepetkin, Ocean Modelling vol. 3 (2001)   
253      !!----------------------------------------------------------------------
254      TYPE(OBC_INDEX),              INTENT(in) ::   idx  ! OBC indices
255      TYPE(OBC_DATA),               INTENT(in) ::   dta  ! OBC external data
256      INTEGER,                      INTENT(in) ::   ib_bdy  ! BDY set index
257      LOGICAL,                      INTENT(in) ::   ll_npo  ! switch for NPO version
258
259      INTEGER  ::   jb, igrd                               ! dummy loop indices
260      !!----------------------------------------------------------------------
261      !
262      !! Note that at this stage the ub and ua arrays contain the baroclinic velocities.
263      !
264      igrd = 2      ! Orlanski bc on u-velocity;
265      !           
266      CALL bdy_orlanski_3d( idx, igrd, ub, ua, dta%u3d, ll_npo )
267
268      igrd = 3      ! Orlanski bc on v-velocity
269     
270      CALL bdy_orlanski_3d( idx, igrd, vb, va, dta%v3d, ll_npo )
271      !
272      CALL lbc_bdy_lnk( 'bdydyn3d', ua, 'U', -1., ib_bdy )    ! Boundary points should be updated
273      CALL lbc_bdy_lnk( 'bdydyn3d', va, 'V', -1., ib_bdy )   
274      !
275   END SUBROUTINE bdy_dyn3d_orlanski
276
277
278   SUBROUTINE bdy_dyn3d_dmp( kt )
279      !!----------------------------------------------------------------------
280      !!                  ***  SUBROUTINE bdy_dyn3d_dmp  ***
281      !!
282      !! ** Purpose : Apply damping for baroclinic velocities at open boundaries.
283      !!
284      !!----------------------------------------------------------------------
285      INTEGER, INTENT(in) ::   kt   ! time step index
286      !
287      INTEGER  ::   jb, jk         ! dummy loop indices
288      INTEGER  ::   ib_bdy         ! loop index
289      INTEGER  ::   ii, ij, igrd   ! local integers
290      REAL(wp) ::   zwgt           ! boundary weight
291      !!----------------------------------------------------------------------
292      !
293      IF( ln_timing )   CALL timing_start('bdy_dyn3d_dmp')
294      !
295      DO ib_bdy=1, nb_bdy
296         IF ( ln_dyn3d_dmp(ib_bdy) .and. cn_dyn3d(ib_bdy) /= 'none' ) THEN
297            igrd = 2                      ! Relaxation of zonal velocity
298            DO jb = 1, idx_bdy(ib_bdy)%nblen(igrd)
299               ii   = idx_bdy(ib_bdy)%nbi(jb,igrd)
300               ij   = idx_bdy(ib_bdy)%nbj(jb,igrd)
301               zwgt = idx_bdy(ib_bdy)%nbd(jb,igrd)
302               DO jk = 1, jpkm1
303                  ua(ii,ij,jk) = ( ua(ii,ij,jk) + zwgt * ( dta_bdy(ib_bdy)%u3d(jb,jk) - &
304                                   ub(ii,ij,jk) + ub_b(ii,ij)) ) * umask(ii,ij,jk)
305               END DO
306            END DO
307            !
308            igrd = 3                      ! Relaxation of meridional velocity
309            DO jb = 1, idx_bdy(ib_bdy)%nblen(igrd)
310               ii   = idx_bdy(ib_bdy)%nbi(jb,igrd)
311               ij   = idx_bdy(ib_bdy)%nbj(jb,igrd)
312               zwgt = idx_bdy(ib_bdy)%nbd(jb,igrd)
313               DO jk = 1, jpkm1
314                  va(ii,ij,jk) = ( va(ii,ij,jk) + zwgt * ( dta_bdy(ib_bdy)%v3d(jb,jk) -  &
315                                   vb(ii,ij,jk) + vb_b(ii,ij)) ) * vmask(ii,ij,jk)
316               END DO
317            END DO
318         ENDIF
319      END DO
320      !
321      CALL lbc_lnk_multi( 'bdydyn3d', ua, 'U', -1.,  va, 'V', -1. )   ! Boundary points should be updated
322      !
323      IF( ln_timing )   CALL timing_stop('bdy_dyn3d_dmp')
324      !
325   END SUBROUTINE bdy_dyn3d_dmp
326
327
328   SUBROUTINE bdy_dyn3d_nmn( idx, ib_bdy )
329      !!----------------------------------------------------------------------
330      !!                 ***  SUBROUTINE bdy_dyn3d_nmn  ***
331      !!             
332      !!              - Apply Neumann condition to baroclinic velocities.
333      !!              - Wrapper routine for bdy_nmn
334      !!
335      !!
336      !!----------------------------------------------------------------------
337      TYPE(OBC_INDEX),              INTENT(in) ::   idx  ! OBC indices
338      INTEGER,                      INTENT(in) ::   ib_bdy  ! BDY set index
339
340      INTEGER  ::   jb, igrd                               ! dummy loop indices
341      !!----------------------------------------------------------------------
342      !
343      !! Note that at this stage the ub and ua arrays contain the baroclinic velocities.
344      !
345      igrd = 2      ! Neumann bc on u-velocity;
346      !           
347      CALL bdy_nmn( idx, igrd, ua )
348
349      igrd = 3      ! Neumann bc on v-velocity
350     
351      CALL bdy_nmn( idx, igrd, va )
352      !
353      CALL lbc_bdy_lnk( 'bdydyn3d', ua, 'U', -1., ib_bdy )    ! Boundary points should be updated
354      CALL lbc_bdy_lnk( 'bdydyn3d', va, 'V', -1., ib_bdy )
355      !
356   END SUBROUTINE bdy_dyn3d_nmn
357
358   !!======================================================================
359END MODULE bdydyn3d
Note: See TracBrowser for help on using the repository browser.