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.
bdydyn2d.F90 in NEMO/branches/2020/dev_r13296_HPC-07_mocavero_mpi3/src/OCE/BDY – NEMO

source: NEMO/branches/2020/dev_r13296_HPC-07_mocavero_mpi3/src/OCE/BDY/bdydyn2d.F90 @ 13630

Last change on this file since 13630 was 13630, checked in by mocavero, 4 years ago

Add neighborhood collectives calls in the NEMO src - ticket #2496

  • Property svn:keywords set to Id
File size: 18.4 KB
Line 
1MODULE bdydyn2d
2   !!======================================================================
3   !!                       ***  MODULE  bdydyn  ***
4   !! Unstructured Open Boundary Cond. :   Apply boundary conditions to barotropic solution
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   !!            3.5  !  2013-07  (J. Chanut) Compliant with time splitting changes
9   !!----------------------------------------------------------------------
10   !!   bdy_dyn2d          : Apply open boundary conditions to barotropic variables.
11   !!   bdy_dyn2d_frs      : Apply Flow Relaxation Scheme
12   !!   bdy_dyn2d_fla      : Apply Flather condition
13   !!   bdy_dyn2d_orlanski : Orlanski Radiation
14   !!   bdy_ssh            : Duplicate sea level across open boundaries
15   !!----------------------------------------------------------------------
16   USE dom_oce         ! ocean space and time domain
17   USE bdy_oce         ! ocean open boundary conditions
18   USE bdylib          ! BDY library routines
19   USE phycst          ! physical constants
20   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
21   USE wet_dry         ! Use wet dry to get reference ssh level
22   USE in_out_manager  !
23   USE lib_mpp, ONLY: ctl_stop
24
25   IMPLICIT NONE
26   PRIVATE
27
28   PUBLIC   bdy_dyn2d   ! routine called in dynspg_ts and bdy_dyn
29   PUBLIC   bdy_ssh       ! routine called in dynspg_ts or sshwzv
30
31   !!----------------------------------------------------------------------
32   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
33   !! $Id$
34   !! Software governed by the CeCILL license (see ./LICENSE)
35   !!----------------------------------------------------------------------
36CONTAINS
37
38   SUBROUTINE bdy_dyn2d( kt, pua2d, pva2d, pub2d, pvb2d, phur, phvr, pssh  )
39      !!----------------------------------------------------------------------
40      !!                  ***  SUBROUTINE bdy_dyn2d  ***
41      !!
42      !! ** Purpose : - Apply open boundary conditions for barotropic variables
43      !!
44      !!----------------------------------------------------------------------
45      INTEGER,                      INTENT(in) ::   kt   ! Main time step counter
46      REAL(wp), DIMENSION(jpi,jpj), INTENT(inout) :: pua2d, pva2d 
47      REAL(wp), DIMENSION(jpi,jpj), INTENT(in   ) :: pub2d, pvb2d
48      REAL(wp), DIMENSION(jpi,jpj), INTENT(in   ) :: phur, phvr
49      REAL(wp), DIMENSION(jpi,jpj), INTENT(in   ) :: pssh
50      !!
51      INTEGER  ::   ib_bdy, ir     ! BDY set index, rim index
52      LOGICAL  ::   llrim0         ! indicate if rim 0 is treated
53      LOGICAL, DIMENSION(4) :: llsend2, llrecv2, llsend3, llrecv3  ! indicate how communications are to be carried out
54     
55      llsend2(:) = .false.   ;   llrecv2(:) = .false.
56      llsend3(:) = .false.   ;   llrecv3(:) = .false.
57      DO ir = 1, 0, -1   ! treat rim 1 before rim 0
58         IF( ir == 0 ) THEN   ;   llrim0 = .TRUE.
59         ELSE                 ;   llrim0 = .FALSE.
60         END IF
61         DO ib_bdy=1, nb_bdy
62            SELECT CASE( cn_dyn2d(ib_bdy) )
63            CASE('none')
64               CYCLE
65            CASE('frs')   ! treat the whole boundary at once
66               IF( llrim0 )   CALL bdy_dyn2d_frs( idx_bdy(ib_bdy), dta_bdy(ib_bdy), ib_bdy, pua2d, pva2d )
67            CASE('flather')
68               CALL bdy_dyn2d_fla( idx_bdy(ib_bdy), dta_bdy(ib_bdy), ib_bdy, pua2d, pva2d, pssh, phur, phvr, llrim0 )
69            CASE('orlanski')
70               CALL bdy_dyn2d_orlanski( idx_bdy(ib_bdy), dta_bdy(ib_bdy), ib_bdy, &
71                    & pua2d, pva2d, pub2d, pvb2d, llrim0, ll_npo=.false. )
72            CASE('orlanski_npo')
73               CALL bdy_dyn2d_orlanski( idx_bdy(ib_bdy), dta_bdy(ib_bdy), ib_bdy, &
74                    & pua2d, pva2d, pub2d, pvb2d, llrim0, ll_npo=.true.  )
75            CASE DEFAULT
76               CALL ctl_stop( 'bdy_dyn2d : unrecognised option for open boundaries for barotropic variables' )
77            END SELECT
78         ENDDO
79         !
80         IF( nn_hls > 1 .AND. ir == 1 ) CYCLE   ! at least 2 halos will be corrected -> no need to correct rim 1 before rim 0
81         IF( nn_hls == 1 ) THEN
82            llsend2(:) = .false.   ;   llrecv2(:) = .false.
83            llsend3(:) = .false.   ;   llrecv3(:) = .false.
84         END IF
85         DO ib_bdy=1, nb_bdy
86            SELECT CASE( cn_dyn2d(ib_bdy) )
87            CASE('flather')
88               llsend2(1:2) = llsend2(1:2) .OR. lsend_bdyint(ib_bdy,2,1:2,ir)   ! west/east, U points
89               llsend2(1)   = llsend2(1)   .OR. lsend_bdyext(ib_bdy,2,1,ir)     ! neighbour might search point towards its east bdy
90               llrecv2(1:2) = llrecv2(1:2) .OR. lrecv_bdyint(ib_bdy,2,1:2,ir)   ! west/east, U points
91               llrecv2(2)   = llrecv2(2)   .OR. lrecv_bdyext(ib_bdy,2,2,ir)     ! might search point towards bdy on the east
92               llsend3(3:4) = llsend3(3:4) .OR. lsend_bdyint(ib_bdy,3,3:4,ir)   ! north/south, V points
93               llsend3(3)   = llsend3(3)   .OR. lsend_bdyext(ib_bdy,3,3,ir)     ! neighbour might search point towards its north bdy
94               llrecv3(3:4) = llrecv3(3:4) .OR. lrecv_bdyint(ib_bdy,3,3:4,ir)   ! north/south, V points
95               llrecv3(4)   = llrecv3(4)   .OR. lrecv_bdyext(ib_bdy,3,4,ir)     ! might search point towards bdy on the north
96            CASE('orlanski', 'orlanski_npo')
97               llsend2(:) = llsend2(:) .OR. lsend_bdy(ib_bdy,2,:,ir)   ! possibly every direction, U points
98               llrecv2(:) = llrecv2(:) .OR. lrecv_bdy(ib_bdy,2,:,ir)   ! possibly every direction, U points
99               llsend3(:) = llsend3(:) .OR. lsend_bdy(ib_bdy,3,:,ir)   ! possibly every direction, V points
100               llrecv3(:) = llrecv3(:) .OR. lrecv_bdy(ib_bdy,3,:,ir)   ! possibly every direction, V points
101            END SELECT
102         END DO
103         IF( ANY(llsend2) .OR. ANY(llrecv2) ) THEN   ! if need to send/recv in at least one direction
104#if defined key_mpi3
105            CALL lbc_lnk_nc_multi( 'bdydyn2d', pua2d, 'U', -1.0_wp, kfillmode=jpfillnothing ,lsend=llsend2, lrecv=llrecv2 )
106#else
107            CALL lbc_lnk( 'bdydyn2d', pua2d, 'U', -1.0_wp, kfillmode=jpfillnothing ,lsend=llsend2, lrecv=llrecv2 )
108#endif
109         END IF
110         IF( ANY(llsend3) .OR. ANY(llrecv3) ) THEN   ! if need to send/recv in at least one direction
111#if defined key_mpi3
112            CALL lbc_lnk_nc_multi( 'bdydyn2d', pva2d, 'V', -1.0_wp, kfillmode=jpfillnothing ,lsend=llsend3, lrecv=llrecv3 )
113#else
114            CALL lbc_lnk( 'bdydyn2d', pva2d, 'V', -1.0_wp, kfillmode=jpfillnothing ,lsend=llsend3, lrecv=llrecv3 )
115#endif
116         END IF
117         !
118      END DO   ! ir
119      !
120   END SUBROUTINE bdy_dyn2d
121
122   SUBROUTINE bdy_dyn2d_frs( idx, dta, ib_bdy, pua2d, pva2d )
123      !!----------------------------------------------------------------------
124      !!                  ***  SUBROUTINE bdy_dyn2d_frs  ***
125      !!
126      !! ** Purpose : - Apply the Flow Relaxation Scheme for barotropic velocities
127      !!                at open boundaries.
128      !!
129      !! References :- Engedahl H., 1995: Use of the flow relaxation scheme in
130      !!               a three-dimensional baroclinic ocean model with realistic
131      !!               topography. Tellus, 365-382.
132      !!----------------------------------------------------------------------
133      TYPE(OBC_INDEX), INTENT(in) ::   idx  ! OBC indices
134      TYPE(OBC_DATA),  INTENT(in) ::   dta  ! OBC external data
135      INTEGER,         INTENT(in) ::   ib_bdy  ! BDY set index
136      REAL(wp), DIMENSION(jpi,jpj), INTENT(inout) :: pua2d, pva2d 
137      !!
138      INTEGER  ::   jb             ! dummy loop indices
139      INTEGER  ::   ii, ij, igrd   ! local integers
140      REAL(wp) ::   zwgt           ! boundary weight
141      !!----------------------------------------------------------------------
142      !
143      igrd = 2                      ! Relaxation of zonal velocity
144      DO jb = 1, idx%nblen(igrd)
145         ii   = idx%nbi(jb,igrd)
146         ij   = idx%nbj(jb,igrd)
147         zwgt = idx%nbw(jb,igrd)
148         pua2d(ii,ij) = ( pua2d(ii,ij) + zwgt * ( dta%u2d(jb) - pua2d(ii,ij) ) ) * umask(ii,ij,1)
149      END DO
150      !
151      igrd = 3                      ! Relaxation of meridional velocity
152      DO jb = 1, idx%nblen(igrd)
153         ii   = idx%nbi(jb,igrd)
154         ij   = idx%nbj(jb,igrd)
155         zwgt = idx%nbw(jb,igrd)
156         pva2d(ii,ij) = ( pva2d(ii,ij) + zwgt * ( dta%v2d(jb) - pva2d(ii,ij) ) ) * vmask(ii,ij,1)
157      END DO 
158      !
159   END SUBROUTINE bdy_dyn2d_frs
160
161
162   SUBROUTINE bdy_dyn2d_fla( idx, dta, ib_bdy, pua2d, pva2d, pssh, phur, phvr, llrim0 )
163      !!----------------------------------------------------------------------
164      !!                 ***  SUBROUTINE bdy_dyn2d_fla  ***
165      !!             
166      !!              - Apply Flather boundary conditions on normal barotropic velocities
167      !!
168      !! ** WARNINGS about FLATHER implementation:
169      !!1. According to Palma and Matano, 1998 "after ssh" is used.
170      !!   In ROMS and POM implementations, it is "now ssh". In the current
171      !!   implementation (tested only in the EEL-R5 conf.), both cases were unstable.
172      !!   So I use "before ssh" in the following.
173      !!
174      !!2. We assume that the normal ssh gradient at the bdy is zero. As a matter of
175      !!   fact, the model ssh just inside the dynamical boundary is used (the outside 
176      !!   ssh in the code is not updated).
177      !!
178      !! References:  Flather, R. A., 1976: A tidal model of the northwest European
179      !!              continental shelf. Mem. Soc. R. Sci. Liege, Ser. 6,10, 141-164.     
180      !!----------------------------------------------------------------------
181      TYPE(OBC_INDEX),              INTENT(in) ::   idx  ! OBC indices
182      TYPE(OBC_DATA),               INTENT(in) ::   dta  ! OBC external data
183      INTEGER,                      INTENT(in) ::   ib_bdy  ! BDY set index
184      REAL(wp), DIMENSION(jpi,jpj), INTENT(inout) :: pua2d, pva2d
185      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   pssh, phur, phvr
186      LOGICAL                     , INTENT(in) ::   llrim0   ! indicate if rim 0 is treated
187      INTEGER  ::   ibeg, iend                       ! length of rim to be treated (rim 0 or rim 1)
188      INTEGER  ::   jb, igrd                         ! dummy loop indices
189      INTEGER  ::   ii, ij                           ! 2D addresses
190      INTEGER  ::   iiTrim, ijTrim                   ! T pts i/j-indice on the rim
191      INTEGER  ::   iiToce, ijToce, iiUoce, ijVoce   ! T, U and V pts i/j-indice of the ocean next to the rim
192      REAL(wp) ::   flagu, flagv                     ! short cuts
193      REAL(wp) ::   zfla                             ! Flather correction
194      REAL(wp) ::   z1_2                             !
195      REAL(wp), DIMENSION(jpi,jpj) ::   sshdta       ! 2D version of dta%ssh
196      !!----------------------------------------------------------------------
197
198      z1_2 = 0.5_wp
199
200      ! ---------------------------------!
201      ! Flather boundary conditions     :!
202      ! ---------------------------------!
203
204      ! Fill temporary array with ssh data (here we use spgu with the alias sshdta):
205      igrd = 1
206      IF( llrim0 ) THEN   ;   ibeg = 1                       ;   iend = idx%nblenrim0(igrd)
207      ELSE                ;   ibeg = idx%nblenrim0(igrd)+1   ;   iend = idx%nblenrim(igrd)
208      END IF
209      !
210      DO jb = ibeg, iend
211         ii = idx%nbi(jb,igrd)
212         ij = idx%nbj(jb,igrd)
213         IF( ll_wd ) THEN   ;   sshdta(ii, ij) = dta%ssh(jb) - ssh_ref 
214         ELSE               ;   sshdta(ii, ij) = dta%ssh(jb)
215         ENDIF
216      END DO
217      !
218      igrd = 2      ! Flather bc on u-velocity
219      !             ! remember that flagu=-1 if normal velocity direction is outward
220      !             ! I think we should rather use after ssh ?
221      IF( llrim0 ) THEN   ;   ibeg = 1                       ;   iend = idx%nblenrim0(igrd)
222      ELSE                ;   ibeg = idx%nblenrim0(igrd)+1   ;   iend = idx%nblenrim(igrd)
223      END IF
224      DO jb = ibeg, iend
225         ii    = idx%nbi(jb,igrd)
226         ij    = idx%nbj(jb,igrd)
227         flagu = idx%flagu(jb,igrd)
228         IF( flagu == 0. ) THEN
229            pua2d(ii,ij) = dta%u2d(jb)
230         ELSE      ! T pts j-indice       on the rim          on the ocean next to the rim on T and U points
231            IF( flagu ==  1. ) THEN   ;   iiTrim = ii     ;   iiToce = ii+1   ;   iiUoce = ii+1   ;   ENDIF
232            IF( flagu == -1. ) THEN   ;   iiTrim = ii+1   ;   iiToce = ii     ;   iiUoce = ii-1   ;   ENDIF
233            !
234            ! Rare case : rim is parallel to the mpi subdomain border and on the halo : point will be received
235            IF( iiTrim > jpi .OR. iiToce > jpi .OR. iiUoce > jpi .OR. iiUoce < 1 )   CYCLE   
236            !
237            zfla = dta%u2d(jb) - flagu * SQRT( grav * phur(ii, ij) ) * ( pssh(iiToce,ij) - sshdta(iiTrim,ij) )
238            !
239            ! jchanut tschanges, use characteristics method (Blayo et Debreu, 2005) :
240            ! mix Flather scheme with velocity of the ocean next to the rim
241            pua2d(ii,ij) =  z1_2 * ( pua2d(iiUoce,ij) + zfla )
242         END IF
243      END DO
244      !
245      igrd = 3      ! Flather bc on v-velocity
246      !             ! remember that flagv=-1 if normal velocity direction is outward
247      IF( llrim0 ) THEN   ;   ibeg = 1                       ;   iend = idx%nblenrim0(igrd)
248      ELSE                ;   ibeg = idx%nblenrim0(igrd)+1   ;   iend = idx%nblenrim(igrd)
249      END IF
250      DO jb = ibeg, iend
251         ii    = idx%nbi(jb,igrd)
252         ij    = idx%nbj(jb,igrd)
253         flagv = idx%flagv(jb,igrd)
254         IF( flagv == 0. ) THEN
255            pva2d(ii,ij) = dta%v2d(jb)
256         ELSE      ! T pts j-indice       on the rim          on the ocean next to the rim on T and V points
257            IF( flagv ==  1. ) THEN   ;   ijTrim = ij     ;   ijToce = ij+1   ;   ijVoce = ij+1   ;   ENDIF
258            IF( flagv == -1. ) THEN   ;   ijTrim = ij+1   ;   ijToce = ij     ;   ijVoce = ij-1   ;   ENDIF
259            !
260            ! Rare case : rim is parallel to the mpi subdomain border and on the halo : point will be received
261            IF( ijTrim > jpj .OR. ijToce > jpj .OR. ijVoce > jpj .OR. ijVoce < 1 )   CYCLE
262            !
263            zfla = dta%v2d(jb) - flagv * SQRT( grav * phvr(ii, ij) ) * ( pssh(ii,ijToce) - sshdta(ii,ijTrim) )
264            !
265            ! jchanut tschanges, use characteristics method (Blayo et Debreu, 2005) :
266            ! mix Flather scheme with velocity of the ocean next to the rim
267            pva2d(ii,ij) =  z1_2 * ( pva2d(ii,ijVoce) + zfla )
268         END IF
269      END DO
270      !
271   END SUBROUTINE bdy_dyn2d_fla
272
273
274   SUBROUTINE bdy_dyn2d_orlanski( idx, dta, ib_bdy, pua2d, pva2d, pub2d, pvb2d, llrim0, ll_npo )
275      !!----------------------------------------------------------------------
276      !!                 ***  SUBROUTINE bdy_dyn2d_orlanski  ***
277      !!             
278      !!              - Apply Orlanski radiation condition adaptively:
279      !!                  - radiation plus weak nudging at outflow points
280      !!                  - no radiation and strong nudging at inflow points
281      !!
282      !!
283      !! References:  Marchesiello, McWilliams and Shchepetkin, Ocean Modelling vol. 3 (2001)   
284      !!----------------------------------------------------------------------
285      TYPE(OBC_INDEX),              INTENT(in) ::   idx  ! OBC indices
286      TYPE(OBC_DATA),               INTENT(in) ::   dta  ! OBC external data
287      INTEGER,                      INTENT(in) ::   ib_bdy  ! number of current open boundary set
288      REAL(wp), DIMENSION(jpi,jpj), INTENT(inout) :: pua2d, pva2d
289      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: pub2d, pvb2d 
290      LOGICAL,                      INTENT(in) ::   ll_npo  ! flag for NPO version
291      LOGICAL,                      INTENT(in) ::   llrim0   ! indicate if rim 0 is treated
292      INTEGER  ::   ib, igrd                               ! dummy loop indices
293      INTEGER  ::   ii, ij, iibm1, ijbm1                   ! indices
294      !!----------------------------------------------------------------------
295      !
296      igrd = 2      ! Orlanski bc on u-velocity;
297      !           
298      CALL bdy_orlanski_2d( idx, igrd, pub2d, pua2d, dta%u2d, llrim0, ll_npo )
299
300      igrd = 3      ! Orlanski bc on v-velocity
301     
302      CALL bdy_orlanski_2d( idx, igrd, pvb2d, pva2d, dta%v2d, llrim0, ll_npo )
303      !
304   END SUBROUTINE bdy_dyn2d_orlanski
305
306
307   SUBROUTINE bdy_ssh( zssh )
308      !!----------------------------------------------------------------------
309      !!                  ***  SUBROUTINE bdy_ssh  ***
310      !!
311      !! ** Purpose : Duplicate sea level across open boundaries
312      !!
313      !!----------------------------------------------------------------------
314      REAL(wp), DIMENSION(jpi,jpj,1), INTENT(inout) ::   zssh ! Sea level, need 3 dimensions to be used by bdy_nmn
315      !!
316      INTEGER ::   ib_bdy, ir      ! bdy index, rim index
317      INTEGER ::   ibeg, iend      ! length of rim to be treated (rim 0 or rim 1)
318      LOGICAL ::   llrim0          ! indicate if rim 0 is treated
319      LOGICAL, DIMENSION(4) :: llsend1, llrecv1  ! indicate how communications are to be carried out
320      !!----------------------------------------------------------------------
321      llsend1(:) = .false.   ;   llrecv1(:) = .false.
322      DO ir = 1, 0, -1   ! treat rim 1 before rim 0
323         IF( nn_hls == 1 ) THEN   ;   llsend1(:) = .false.   ;   llrecv1(:) = .false.   ;   END IF
324         IF( ir == 0 ) THEN   ;   llrim0 = .TRUE.
325         ELSE                 ;   llrim0 = .FALSE.
326         END IF
327         DO ib_bdy = 1, nb_bdy
328            CALL bdy_nmn( idx_bdy(ib_bdy), 1, zssh, llrim0 )   ! zssh is masked
329            llsend1(:) = llsend1(:) .OR. lsend_bdyint(ib_bdy,1,:,ir)   ! possibly every direction, T points
330            llrecv1(:) = llrecv1(:) .OR. lrecv_bdyint(ib_bdy,1,:,ir)   ! possibly every direction, T points
331         END DO
332         IF( nn_hls > 1 .AND. ir == 1 ) CYCLE   ! at least 2 halos will be corrected -> no need to correct rim 1 before rim 0
333         IF( ANY(llsend1) .OR. ANY(llrecv1) ) THEN   ! if need to send/recv in at least one direction
334#if defined key_mpi3
335            CALL lbc_lnk_nc_multi( 'bdydyn2d', zssh(:,:,1), 'T',  1.0_wp, kfillmode=jpfillnothing ,lsend=llsend1, lrecv=llrecv1 )
336#else
337            CALL lbc_lnk( 'bdydyn2d', zssh(:,:,1), 'T',  1.0_wp, kfillmode=jpfillnothing ,lsend=llsend1, lrecv=llrecv1 )
338#endif
339         END IF
340      END DO
341      !
342   END SUBROUTINE bdy_ssh
343
344   !!======================================================================
345END MODULE bdydyn2d
346
Note: See TracBrowser for help on using the repository browser.