1 | MODULE 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 | !!---------------------------------------------------------------------- |
---|
36 | CONTAINS |
---|
37 | |
---|
38 | SUBROUTINE bdy_dyn2d( kt, pua2d, pva2d, pub2d, pvb2d, phur, phvr, pssh & |
---|
39 | & , kdbi, kdei, kdbj, kdej, ldcomall, pumask, pvmask, khlcom ) |
---|
40 | !!---------------------------------------------------------------------- |
---|
41 | !! *** SUBROUTINE bdy_dyn2d *** |
---|
42 | !! |
---|
43 | !! ** Purpose : - Apply open boundary conditions for barotropic variables |
---|
44 | !! |
---|
45 | !!---------------------------------------------------------------------- |
---|
46 | INTEGER, INTENT(in ) :: kt ! Main time step counter |
---|
47 | REAL(wp), DIMENSION(kdbi:kdei,kdbj:kdej), INTENT(inout) :: pua2d, pva2d |
---|
48 | REAL(wp), DIMENSION(kdbi:kdei,kdbj:kdej), INTENT(in ) :: pub2d, pvb2d |
---|
49 | REAL(wp), DIMENSION(kdbi:kdei,kdbj:kdej), INTENT(in ) :: phur, phvr |
---|
50 | REAL(wp), DIMENSION(kdbi:kdei,kdbj:kdej), INTENT(in ) :: pssh |
---|
51 | INTEGER , INTENT(in ) :: kdbi, kdei, kdbj, kdej ! size of array |
---|
52 | LOGICAL , OPTIONAL, INTENT(in ) :: ldcomall ! communicate with all neighbours |
---|
53 | REAL(wp), OPTIONAL, DIMENSION(kdbi:kdei,kdbj:kdej), INTENT(in ) :: pumask ! optional mask for extended domain |
---|
54 | REAL(wp), OPTIONAL, DIMENSION(kdbi:kdei,kdbj:kdej), INTENT(in ) :: pvmask ! - - |
---|
55 | INTEGER , OPTIONAL, INTENT(in ) :: khlcom ! number of halos to communicate |
---|
56 | !! |
---|
57 | INTEGER :: ib_bdy, ir ! BDY set index, rim index |
---|
58 | LOGICAL :: llrim0 ! indicate if rim 0 is treated |
---|
59 | LOGICAL, DIMENSION(4) :: llsend2, llrecv2, llsend3, llrecv3 ! indicate how communications are to be carried out |
---|
60 | |
---|
61 | llsend2(:) = .false. ; llrecv2(:) = .false. |
---|
62 | llsend3(:) = .false. ; llrecv3(:) = .false. |
---|
63 | DO ir = 1, 0, -1 ! treat rim 1 before rim 0 |
---|
64 | IF( ir == 0 ) THEN ; llrim0 = .TRUE. |
---|
65 | ELSE ; llrim0 = .FALSE. |
---|
66 | END IF |
---|
67 | DO ib_bdy=1, nb_bdy |
---|
68 | SELECT CASE( cn_dyn2d(ib_bdy) ) |
---|
69 | CASE('none') |
---|
70 | CYCLE |
---|
71 | CASE('frs') ! treat the whole boundary at once |
---|
72 | IF( llrim0 ) CALL bdy_dyn2d_frs( idx_bdy(ib_bdy), dta_bdy(ib_bdy), ib_bdy, pua2d, pva2d, & |
---|
73 | & kdbi, kdei, kdbj, kdej, pumask=pumask, pvmask=pvmask ) |
---|
74 | CASE('flather') |
---|
75 | CALL bdy_dyn2d_fla( idx_bdy(ib_bdy), dta_bdy(ib_bdy), ib_bdy, pua2d, pva2d, pssh, phur, phvr, llrim0, & |
---|
76 | & kdbi, kdei, kdbj, kdej ) |
---|
77 | CASE('orlanski') |
---|
78 | CALL bdy_dyn2d_orlanski( idx_bdy(ib_bdy), dta_bdy(ib_bdy), ib_bdy, & |
---|
79 | & pua2d, pva2d, pub2d, pvb2d, .false., llrim0, kdbi, kdei, kdbj, kdej ) |
---|
80 | CASE('orlanski_npo') |
---|
81 | CALL bdy_dyn2d_orlanski( idx_bdy(ib_bdy), dta_bdy(ib_bdy), ib_bdy, & |
---|
82 | & pua2d, pva2d, pub2d, pvb2d, .true. , llrim0, kdbi, kdei, kdbj, kdej ) |
---|
83 | CASE DEFAULT |
---|
84 | CALL ctl_stop( 'bdy_dyn2d : unrecognised option for open boundaries for barotropic variables' ) |
---|
85 | END SELECT |
---|
86 | END DO |
---|
87 | ! |
---|
88 | IF( nn_hls > 1 .AND. ir == 1 ) CYCLE ! at least 2 halos will be corrected -> no need to correct rim 1 before rim 0 |
---|
89 | IF( nn_hls == 1 ) THEN |
---|
90 | llsend2(:) = .false. ; llrecv2(:) = .false. |
---|
91 | llsend3(:) = .false. ; llrecv3(:) = .false. |
---|
92 | END IF |
---|
93 | DO ib_bdy=1, nb_bdy |
---|
94 | SELECT CASE( cn_dyn2d(ib_bdy) ) |
---|
95 | CASE('flather') |
---|
96 | llsend2(1:2) = llsend2(1:2) .OR. lsend_bdyint(ib_bdy,2,1:2,ir) ! west/east, U points |
---|
97 | llsend2(1) = llsend2(1) .OR. lsend_bdyext(ib_bdy,2,1,ir) ! neighbour might search point towards its east bdy |
---|
98 | llrecv2(1:2) = llrecv2(1:2) .OR. lrecv_bdyint(ib_bdy,2,1:2,ir) ! west/east, U points |
---|
99 | llrecv2(2) = llrecv2(2) .OR. lrecv_bdyext(ib_bdy,2,2,ir) ! might search point towards bdy on the east |
---|
100 | llsend3(3:4) = llsend3(3:4) .OR. lsend_bdyint(ib_bdy,3,3:4,ir) ! north/south, V points |
---|
101 | llsend3(3) = llsend3(3) .OR. lsend_bdyext(ib_bdy,3,3,ir) ! neighbour might search point towards its north bdy |
---|
102 | llrecv3(3:4) = llrecv3(3:4) .OR. lrecv_bdyint(ib_bdy,3,3:4,ir) ! north/south, V points |
---|
103 | llrecv3(4) = llrecv3(4) .OR. lrecv_bdyext(ib_bdy,3,4,ir) ! might search point towards bdy on the north |
---|
104 | CASE('orlanski', 'orlanski_npo') |
---|
105 | llsend2(:) = llsend2(:) .OR. lsend_bdy(ib_bdy,2,:,ir) ! possibly every direction, U points |
---|
106 | llrecv2(:) = llrecv2(:) .OR. lrecv_bdy(ib_bdy,2,:,ir) ! possibly every direction, U points |
---|
107 | llsend3(:) = llsend3(:) .OR. lsend_bdy(ib_bdy,3,:,ir) ! possibly every direction, V points |
---|
108 | llrecv3(:) = llrecv3(:) .OR. lrecv_bdy(ib_bdy,3,:,ir) ! possibly every direction, V points |
---|
109 | END SELECT |
---|
110 | END DO |
---|
111 | IF( PRESENT(ldcomall) ) THEN |
---|
112 | IF( ldcomall ) THEN ! if ldcomall is present and true then communicate with all neighbours |
---|
113 | CALL lbc_lnk_multi( 'bdydyn2d', pua2d, 'U', 1., pva2d, 'V', 1., kfillmode=jpfillnothing, khlcom=khlcom ) |
---|
114 | CYCLE |
---|
115 | END IF |
---|
116 | END IF |
---|
117 | IF( ANY(llsend2) .OR. ANY(llrecv2) ) THEN ! if need to send/recv in at least one direction |
---|
118 | CALL lbc_lnk( 'bdydyn2d', pua2d, 'U', -1., kfillmode=jpfillnothing ,ldsend=llsend2, ldrecv=llrecv2, khlcom=khlcom ) |
---|
119 | END IF |
---|
120 | IF( ANY(llsend3) .OR. ANY(llrecv3) ) THEN ! if need to send/recv in at least one direction |
---|
121 | CALL lbc_lnk( 'bdydyn2d', pva2d, 'V', -1., kfillmode=jpfillnothing ,ldsend=llsend3, ldrecv=llrecv3, khlcom=khlcom ) |
---|
122 | END IF |
---|
123 | ! |
---|
124 | END DO ! ir |
---|
125 | ! |
---|
126 | END SUBROUTINE bdy_dyn2d |
---|
127 | |
---|
128 | SUBROUTINE bdy_dyn2d_frs( idx, dta, ib_bdy, pua2d, pva2d, kdbi, kdei, kdbj, kdej, pumask, pvmask ) |
---|
129 | !!---------------------------------------------------------------------- |
---|
130 | !! *** SUBROUTINE bdy_dyn2d_frs *** |
---|
131 | !! |
---|
132 | !! ** Purpose : - Apply the Flow Relaxation Scheme for barotropic velocities |
---|
133 | !! at open boundaries. |
---|
134 | !! |
---|
135 | !! References :- Engedahl H., 1995: Use of the flow relaxation scheme in |
---|
136 | !! a three-dimensional baroclinic ocean model with realistic |
---|
137 | !! topography. Tellus, 365-382. |
---|
138 | !!---------------------------------------------------------------------- |
---|
139 | TYPE(OBC_INDEX), INTENT(in ) :: idx ! OBC indices |
---|
140 | TYPE(OBC_DATA), INTENT(in ) :: dta ! OBC external data |
---|
141 | INTEGER, INTENT(in ) :: ib_bdy ! BDY set index |
---|
142 | REAL(wp), DIMENSION(kdbi:kdei,kdbj:kdej), INTENT(inout) :: pua2d, pva2d |
---|
143 | INTEGER , INTENT(in ) :: kdbi, kdei, kdbj, kdej ! size of array |
---|
144 | REAL(wp), OPTIONAL, TARGET, DIMENSION(kdbi:kdei,kdbj:kdej,1), INTENT(in ) :: pumask ! optional mask for extended domain |
---|
145 | REAL(wp), OPTIONAL, TARGET, DIMENSION(kdbi:kdei,kdbj:kdej,1), INTENT(in ) :: pvmask ! - - |
---|
146 | !! |
---|
147 | INTEGER :: jb ! dummy loop indices |
---|
148 | INTEGER :: ii, ij, igrd ! local integers |
---|
149 | REAL(wp) :: zwgt ! boundary weight |
---|
150 | REAL(wp), POINTER, DIMENSION(:,:,:) :: pmask ! land/sea mask for field |
---|
151 | !!---------------------------------------------------------------------- |
---|
152 | ! |
---|
153 | igrd = 2 ! Relaxation of zonal velocity |
---|
154 | IF( PRESENT(pumask) ) THEN ; pmask => pumask |
---|
155 | ELSE ; pmask => umask |
---|
156 | END IF |
---|
157 | DO jb = 1, idx%nblen(igrd) |
---|
158 | ii = idx%nbi(jb,igrd) |
---|
159 | ij = idx%nbj(jb,igrd) |
---|
160 | zwgt = idx%nbw(jb,igrd) |
---|
161 | pua2d(ii,ij) = ( pua2d(ii,ij) + zwgt * ( dta%u2d(jb) - pua2d(ii,ij) ) ) * pmask(ii,ij,1) |
---|
162 | END DO |
---|
163 | ! |
---|
164 | igrd = 3 ! Relaxation of meridional velocity |
---|
165 | IF( PRESENT(pvmask) ) THEN ; pmask => pvmask |
---|
166 | ELSE ; pmask => vmask |
---|
167 | END IF |
---|
168 | DO jb = 1, idx%nblen(igrd) |
---|
169 | ii = idx%nbi(jb,igrd) |
---|
170 | ij = idx%nbj(jb,igrd) |
---|
171 | zwgt = idx%nbw(jb,igrd) |
---|
172 | pva2d(ii,ij) = ( pva2d(ii,ij) + zwgt * ( dta%v2d(jb) - pva2d(ii,ij) ) ) * pmask(ii,ij,1) |
---|
173 | END DO |
---|
174 | ! |
---|
175 | END SUBROUTINE bdy_dyn2d_frs |
---|
176 | |
---|
177 | |
---|
178 | SUBROUTINE bdy_dyn2d_fla( idx, dta, ib_bdy, pua2d, pva2d, pssh, phur, phvr, ldrim0 & |
---|
179 | & , kdbi, kdei, kdbj, kdej ) |
---|
180 | !!---------------------------------------------------------------------- |
---|
181 | !! *** SUBROUTINE bdy_dyn2d_fla *** |
---|
182 | !! |
---|
183 | !! - Apply Flather boundary conditions on normal barotropic velocities |
---|
184 | !! |
---|
185 | !! ** WARNINGS about FLATHER implementation: |
---|
186 | !!1. According to Palma and Matano, 1998 "after ssh" is used. |
---|
187 | !! In ROMS and POM implementations, it is "now ssh". In the current |
---|
188 | !! implementation (tested only in the EEL-R5 conf.), both cases were unstable. |
---|
189 | !! So I use "before ssh" in the following. |
---|
190 | !! |
---|
191 | !!2. We assume that the normal ssh gradient at the bdy is zero. As a matter of |
---|
192 | !! fact, the model ssh just inside the dynamical boundary is used (the outside |
---|
193 | !! ssh in the code is not updated). |
---|
194 | !! |
---|
195 | !! References: Flather, R. A., 1976: A tidal model of the northwest European |
---|
196 | !! continental shelf. Mem. Soc. R. Sci. Liege, Ser. 6,10, 141-164. |
---|
197 | !!---------------------------------------------------------------------- |
---|
198 | TYPE(OBC_INDEX), INTENT(in ) :: idx ! OBC indices |
---|
199 | TYPE(OBC_DATA), INTENT(in ) :: dta ! OBC external data |
---|
200 | INTEGER, INTENT(in ) :: ib_bdy ! BDY set index |
---|
201 | REAL(wp), DIMENSION(kdbi:kdei,kdbj:kdej), INTENT(inout) :: pua2d, pva2d |
---|
202 | REAL(wp), DIMENSION(kdbi:kdei,kdbj:kdej), INTENT(in ) :: pssh, phur, phvr |
---|
203 | LOGICAL , INTENT(in ) :: ldrim0 ! indicate if rim 0 is treated |
---|
204 | INTEGER , INTENT(in ) :: kdbi, kdei, kdbj, kdej ! size of array |
---|
205 | !! |
---|
206 | INTEGER :: ibeg, iend ! length of rim to be treated (rim 0 or rim 1) |
---|
207 | INTEGER :: jb, igrd ! dummy loop indices |
---|
208 | INTEGER :: ii, ij ! 2D addresses |
---|
209 | INTEGER :: iiTrim, ijTrim ! T pts i/j-indice on the rim |
---|
210 | INTEGER :: iiToce, ijToce, iiUoce, ijVoce ! T, U and V pts i/j-indice of the ocean next to the rim |
---|
211 | REAL(wp) :: flagu, flagv ! short cuts |
---|
212 | REAL(wp) :: zfla ! Flather correction |
---|
213 | REAL(wp) :: z1_2 ! |
---|
214 | REAL(wp), DIMENSION(kdbi:kdei,kdbj:kdej) :: sshdta ! 2D version of dta%ssh |
---|
215 | !!---------------------------------------------------------------------- |
---|
216 | |
---|
217 | z1_2 = 0.5_wp |
---|
218 | |
---|
219 | ! ---------------------------------! |
---|
220 | ! Flather boundary conditions :! |
---|
221 | ! ---------------------------------! |
---|
222 | |
---|
223 | ! Fill temporary array with ssh data (here we use spgu with the alias sshdta): |
---|
224 | igrd = 1 |
---|
225 | IF( ldrim0 ) THEN ; ibeg = 1 ; iend = idx%nblenrim0(igrd) |
---|
226 | ELSE ; ibeg = idx%nblenrim0(igrd)+1 ; iend = idx%nblenrim(igrd) |
---|
227 | END IF |
---|
228 | ! |
---|
229 | DO jb = ibeg, iend |
---|
230 | ii = idx%nbi(jb,igrd) |
---|
231 | ij = idx%nbj(jb,igrd) |
---|
232 | IF( ll_wd ) THEN ; sshdta(ii, ij) = dta%ssh(jb) - ssh_ref |
---|
233 | ELSE ; sshdta(ii, ij) = dta%ssh(jb) |
---|
234 | ENDIF |
---|
235 | END DO |
---|
236 | ! |
---|
237 | igrd = 2 ! Flather bc on u-velocity |
---|
238 | ! ! remember that flagu=-1 if normal velocity direction is outward |
---|
239 | ! ! I think we should rather use after ssh ? |
---|
240 | IF( ldrim0 ) THEN ; ibeg = 1 ; iend = idx%nblenrim0(igrd) |
---|
241 | ELSE ; ibeg = idx%nblenrim0(igrd)+1 ; iend = idx%nblenrim(igrd) |
---|
242 | END IF |
---|
243 | DO jb = ibeg, iend |
---|
244 | ii = idx%nbi(jb,igrd) |
---|
245 | ij = idx%nbj(jb,igrd) |
---|
246 | flagu = idx%flagu(jb,igrd) |
---|
247 | IF( flagu == 0. ) THEN |
---|
248 | pua2d(ii,ij) = dta%u2d(jb) |
---|
249 | ELSE ! T pts j-indice on the rim on the ocean next to the rim on T and U points |
---|
250 | IF( flagu == 1. ) THEN ; iiTrim = ii ; iiToce = ii+1 ; iiUoce = ii+1 ; ENDIF |
---|
251 | IF( flagu == -1. ) THEN ; iiTrim = ii+1 ; iiToce = ii ; iiUoce = ii-1 ; ENDIF |
---|
252 | ! |
---|
253 | ! Rare case : rim is parallel to the mpi subdomain border and on the halo : point will be received |
---|
254 | IF( iiTrim > kdei .OR. iiToce > kdei .OR. iiUoce > kdei .OR. iiUoce < kdbi ) CYCLE |
---|
255 | ! |
---|
256 | zfla = dta%u2d(jb) - flagu * SQRT( grav * phur(ii, ij) ) * ( pssh(iiToce,ij) - sshdta(iiTrim,ij) ) |
---|
257 | ! |
---|
258 | ! jchanut tschanges, use characteristics method (Blayo et Debreu, 2005) : |
---|
259 | ! mix Flather scheme with velocity of the ocean next to the rim |
---|
260 | pua2d(ii,ij) = z1_2 * ( pua2d(iiUoce,ij) + zfla ) |
---|
261 | END IF |
---|
262 | END DO |
---|
263 | ! |
---|
264 | igrd = 3 ! Flather bc on v-velocity |
---|
265 | ! ! remember that flagv=-1 if normal velocity direction is outward |
---|
266 | IF( ldrim0 ) THEN ; ibeg = 1 ; iend = idx%nblenrim0(igrd) |
---|
267 | ELSE ; ibeg = idx%nblenrim0(igrd)+1 ; iend = idx%nblenrim(igrd) |
---|
268 | END IF |
---|
269 | DO jb = ibeg, iend |
---|
270 | ii = idx%nbi(jb,igrd) |
---|
271 | ij = idx%nbj(jb,igrd) |
---|
272 | flagv = idx%flagv(jb,igrd) |
---|
273 | IF( flagv == 0. ) THEN |
---|
274 | pva2d(ii,ij) = dta%v2d(jb) |
---|
275 | ELSE ! T pts j-indice on the rim on the ocean next to the rim on T and V points |
---|
276 | IF( flagv == 1. ) THEN ; ijTrim = ij ; ijToce = ij+1 ; ijVoce = ij+1 ; ENDIF |
---|
277 | IF( flagv == -1. ) THEN ; ijTrim = ij+1 ; ijToce = ij ; ijVoce = ij-1 ; ENDIF |
---|
278 | ! |
---|
279 | ! Rare case : rim is parallel to the mpi subdomain border and on the halo : point will be received |
---|
280 | IF( ijTrim > kdej .OR. ijToce > kdej .OR. ijVoce > kdej .OR. ijVoce < kdbj ) CYCLE |
---|
281 | ! |
---|
282 | zfla = dta%v2d(jb) - flagv * SQRT( grav * phvr(ii, ij) ) * ( pssh(ii,ijToce) - sshdta(ii,ijTrim) ) |
---|
283 | ! |
---|
284 | ! jchanut tschanges, use characteristics method (Blayo et Debreu, 2005) : |
---|
285 | ! mix Flather scheme with velocity of the ocean next to the rim |
---|
286 | pva2d(ii,ij) = z1_2 * ( pva2d(ii,ijVoce) + zfla ) |
---|
287 | END IF |
---|
288 | END DO |
---|
289 | ! |
---|
290 | END SUBROUTINE bdy_dyn2d_fla |
---|
291 | |
---|
292 | |
---|
293 | SUBROUTINE bdy_dyn2d_orlanski( idx, dta, ib_bdy, pua2d, pva2d, pub2d, pvb2d, ld_npo, ldrim0, kdbi, kdei, kdbj, kdej ) |
---|
294 | !!---------------------------------------------------------------------- |
---|
295 | !! *** SUBROUTINE bdy_dyn2d_orlanski *** |
---|
296 | !! |
---|
297 | !! - Apply Orlanski radiation condition adaptively: |
---|
298 | !! - radiation plus weak nudging at outflow points |
---|
299 | !! - no radiation and strong nudging at inflow points |
---|
300 | !! |
---|
301 | !! |
---|
302 | !! References: Marchesiello, McWilliams and Shchepetkin, Ocean Modelling vol. 3 (2001) |
---|
303 | !!---------------------------------------------------------------------- |
---|
304 | TYPE(OBC_INDEX), INTENT(in ) :: idx ! OBC indices |
---|
305 | TYPE(OBC_DATA), INTENT(in ) :: dta ! OBC external data |
---|
306 | INTEGER, INTENT(in ) :: ib_bdy ! number of current open boundary set |
---|
307 | REAL(wp), DIMENSION(kdbi:kdei,kdbj:kdej), INTENT(inout) :: pua2d, pva2d |
---|
308 | REAL(wp), DIMENSION(kdbi:kdei,kdbj:kdej), INTENT(in ) :: pub2d, pvb2d |
---|
309 | LOGICAL, INTENT(in ) :: ld_npo ! flag for NPO version |
---|
310 | LOGICAL, INTENT(in ) :: ldrim0 ! indicate if rim 0 is treated |
---|
311 | INTEGER , INTENT(in ) :: kdbi, kdei, kdbj, kdej ! size of array |
---|
312 | INTEGER :: ib, igrd ! dummy loop indices |
---|
313 | INTEGER :: ii, ij, iibm1, ijbm1 ! indices |
---|
314 | !!---------------------------------------------------------------------- |
---|
315 | ! |
---|
316 | igrd = 2 ! Orlanski bc on u-velocity; |
---|
317 | ! |
---|
318 | CALL bdy_orlanski_2d( idx, igrd, pub2d, pua2d, dta%u2d, ldrim0, ld_npo ) |
---|
319 | |
---|
320 | igrd = 3 ! Orlanski bc on v-velocity |
---|
321 | ! |
---|
322 | CALL bdy_orlanski_2d( idx, igrd, pvb2d, pva2d, dta%v2d, ldrim0, ld_npo ) |
---|
323 | ! |
---|
324 | END SUBROUTINE bdy_dyn2d_orlanski |
---|
325 | |
---|
326 | |
---|
327 | SUBROUTINE bdy_ssh( zssh, kdbi, kdei, kdbj, kdej, ldcomall, pmask, khlcom ) |
---|
328 | !!---------------------------------------------------------------------- |
---|
329 | !! *** SUBROUTINE bdy_ssh *** |
---|
330 | !! |
---|
331 | !! ** Purpose : Duplicate sea level across open boundaries |
---|
332 | !! |
---|
333 | !!---------------------------------------------------------------------- |
---|
334 | REAL(wp), DIMENSION(kdbi:kdei,kdbj:kdej,1), INTENT(inout) :: zssh ! Sea level, need 3 dimensions to be used by bdy_nmn |
---|
335 | INTEGER , INTENT(in ) :: kdbi, kdei, kdbj, kdej ! size of ssh array |
---|
336 | LOGICAL , OPTIONAL, INTENT(in ) :: ldcomall ! communicate with all neighbours |
---|
337 | REAL(wp), OPTIONAL, DIMENSION(kdbi:kdei,kdbj:kdej,1), INTENT(in) :: pmask ! optional mask for extended domain |
---|
338 | INTEGER , OPTIONAL, INTENT(in) :: khlcom ! number of halos to communicate |
---|
339 | !! |
---|
340 | INTEGER :: ib_bdy, ir ! bdy index, rim index |
---|
341 | INTEGER :: ibeg, iend ! length of rim to be treated (rim 0 or rim 1) |
---|
342 | INTEGER :: ihl ! thickness of halo |
---|
343 | LOGICAL :: llrim0 ! indicate if rim 0 is treated |
---|
344 | LOGICAL, DIMENSION(4) :: llsend1, llrecv1 ! indicate how communications are to be carried out |
---|
345 | !!---------------------------------------------------------------------- |
---|
346 | llsend1(:) = .false. ; llrecv1(:) = .false. |
---|
347 | ! |
---|
348 | DO ir = 1, 0, -1 ! treat rim 1 before rim 0 |
---|
349 | IF( nn_hls == 1 ) THEN ; llsend1(:) = .false. ; llrecv1(:) = .false. ; END IF |
---|
350 | IF( ir == 0 ) THEN ; llrim0 = .TRUE. |
---|
351 | ELSE ; llrim0 = .FALSE. |
---|
352 | END IF |
---|
353 | DO ib_bdy = 1, nb_bdy |
---|
354 | CALL bdy_nmn( idx_bdy(ib_bdy), 1, kdbi, kdei, kdbj, kdej, 1, zssh, ldrim0=llrim0, pmask=pmask ) ! zssh is masked |
---|
355 | llsend1(:) = llsend1(:) .OR. lsend_bdyint(ib_bdy,1,:,ir) ! possibly every direction, T points |
---|
356 | llrecv1(:) = llrecv1(:) .OR. lrecv_bdyint(ib_bdy,1,:,ir) ! possibly every direction, T points |
---|
357 | END DO |
---|
358 | IF( nn_hls > 1 .AND. ir == 1 ) CYCLE ! at least 2 halos will be corrected -> no need to correct rim 1 before rim 0 |
---|
359 | ! N.B. ihl>1 is not enough as values are usually wrong on extended domain |
---|
360 | IF( PRESENT(ldcomall) ) THEN |
---|
361 | IF( ldcomall ) THEN ! if ldcomall is present and true |
---|
362 | CALL lbc_lnk( 'bdydyn2d', zssh(:,:,1), 'T', 1., kfillmode=jpfillnothing, khlcom=khlcom ) ! com with all neighbours |
---|
363 | END IF |
---|
364 | ELSEIF( ANY(llsend1) .OR. ANY(llrecv1) ) THEN ! if need to send/recv in at least one direction |
---|
365 | CALL lbc_lnk( 'bdydyn2d', zssh(:,:,1), 'T', 1., kfillmode=jpfillnothing, ldsend=llsend1, ldrecv=llrecv1, khlcom=khlcom ) |
---|
366 | END IF |
---|
367 | END DO |
---|
368 | ! |
---|
369 | END SUBROUTINE bdy_ssh |
---|
370 | |
---|
371 | !!====================================================================== |
---|
372 | END MODULE bdydyn2d |
---|
373 | |
---|