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/2019/dev_r10984_HPC-13_IRRMANN_BDY_optimization/src/OCE/BDY – NEMO

source: NEMO/branches/2019/dev_r10984_HPC-13_IRRMANN_BDY_optimization/src/OCE/BDY/bdydyn3d.F90 @ 11049

Last change on this file since 11049 was 11049, checked in by girrmann, 5 years ago

dev_r10984_HPC-13 : CYCLE instruction is not systematic anymore, computation is done on the halo whenever possible and overwritten by lbc_bdy instruction, see #2285

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