ROMS
Loading...
Searching...
No Matches
wetdry.F
Go to the documentation of this file.
1#include "cppdefs.h"
2 MODULE wetdry_mod
3#ifdef WET_DRY
4!
5!git $Id$
6!=======================================================================
7! Copyright (c) 2002-2025 The ROMS Group !
8! Licensed under a MIT/X style license !
9! See License_ROMS.md Hernan G. Arango !
10!==================================================== John C. Warner ===
11! !
12! This routine computes the wet/dry masking arrays. !
13! !
14!=======================================================================
15!
16 implicit none
17!
18 CONTAINS
19!
20!***********************************************************************
21 SUBROUTINE wetdry (ng, tile, Tindex, Linitialize)
22!***********************************************************************
23!
24 USE mod_param
25 USE mod_coupling
26 USE mod_grid
27 USE mod_ocean
28!
29! Imported variable declarations.
30!
31 logical, intent(in) :: Linitialize
32 integer, intent(in) :: ng, tile, Tindex
33!
34! Local variable declarations.
35!
36# include "tile.h"
37!
38 IF (linitialize) THEN
39 CALL wetdry_ini_tile (ng, tile, &
40 & lbi, ubi, lbj, ubj, &
41 & imins, imaxs, jmins, jmaxs, &
42# ifdef MASKING
43 & grid(ng)%pmask, &
44 & grid(ng)%rmask, &
45 & grid(ng)%umask, &
46 & grid(ng)%vmask, &
47# endif
48 & grid(ng)%h, &
49 & ocean(ng)%zeta(:,:,tindex), &
50# ifdef SOLVE3D
51 & ocean(ng)%ubar(:,:,tindex), &
52 & ocean(ng)%vbar(:,:,tindex), &
53# endif
54 & grid(ng)%pmask_wet, &
55 & grid(ng)%pmask_full, &
56 & grid(ng)%rmask_wet, &
57 & grid(ng)%rmask_full, &
58 & grid(ng)%umask_wet, &
59 & grid(ng)%umask_full, &
60 & grid(ng)%vmask_wet, &
61 & grid(ng)%vmask_full)
62 ELSE
63 CALL wetdry_tile (ng, tile, &
64 & lbi, ubi, lbj, ubj, &
65 & imins, imaxs, jmins, jmaxs, &
66# ifdef MASKING
67 & grid(ng)%pmask, &
68 & grid(ng)%rmask, &
69 & grid(ng)%umask, &
70 & grid(ng)%vmask, &
71# endif
72 & grid(ng)%h, &
73 & ocean(ng)%zeta(:,:,tindex), &
74# ifdef SOLVE3D
75 & coupling(ng)%DU_avg1, &
76 & coupling(ng)%DV_avg1, &
77 & grid(ng)%rmask_wet_avg, &
78# endif
79 & grid(ng)%pmask_wet, &
80 & grid(ng)%pmask_full, &
81 & grid(ng)%rmask_wet, &
82 & grid(ng)%rmask_full, &
83 & grid(ng)%umask_wet, &
84 & grid(ng)%umask_full, &
85 & grid(ng)%vmask_wet, &
86 & grid(ng)%vmask_full)
87 END IF
88
89 RETURN
90 END SUBROUTINE wetdry
91!
92!***********************************************************************
93 SUBROUTINE wetdry_tile (ng, tile, &
94 & LBi, UBi, LBj, UBj, &
95 & IminS, ImaxS, JminS, JmaxS, &
96# ifdef MASKING
97 & pmask, rmask, umask, vmask, &
98# endif
99 & h, zeta, &
100# ifdef SOLVE3D
101 & DU_avg1, DV_avg1, &
102 & rmask_wet_avg, &
103# endif
104 & pmask_wet, pmask_full, &
105 & rmask_wet, rmask_full, &
106 & umask_wet, umask_full, &
107 & vmask_wet, vmask_full)
108!***********************************************************************
109!
110 USE mod_param
111 USE mod_ncparam
112 USE mod_scalars
113 USE mod_sources
114!
116# ifdef DISTRIBUTE
118# endif
119!
120! Imported variable declarations.
121!
122 integer, intent(in) :: ng, tile
123 integer, intent(in) :: LBi, UBi, LBj, UBj
124 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
125!
126# ifdef ASSUMED_SHAPE
127 real(r8), intent(in) :: h(LBi:,LBj:)
128# ifdef MASKING
129 real(r8), intent(in) :: pmask(LBi:,LBj:)
130 real(r8), intent(in) :: rmask(LBi:,LBj:)
131 real(r8), intent(in) :: umask(LBi:,LBj:)
132 real(r8), intent(in) :: vmask(LBi:,LBj:)
133# endif
134 real(r8), intent(in) :: zeta(LBi:,LBj:)
135# ifdef SOLVE3D
136 real(r8), intent(in) :: DU_avg1(LBi:,LBj:)
137 real(r8), intent(in) :: DV_avg1(LBi:,LBj:)
138 real(r8), intent(inout) :: rmask_wet_avg(LBi:,LBj:)
139# endif
140 real(r8), intent(out) :: pmask_full(LBi:,LBj:)
141 real(r8), intent(out) :: rmask_full(LBi:,LBj:)
142 real(r8), intent(out) :: umask_full(LBi:,LBj:)
143 real(r8), intent(out) :: vmask_full(LBi:,LBj:)
144
145 real(r8), intent(out) :: pmask_wet(LBi:,LBj:)
146 real(r8), intent(out) :: rmask_wet(LBi:,LBj:)
147 real(r8), intent(out) :: umask_wet(LBi:,LBj:)
148 real(r8), intent(out) :: vmask_wet(LBi:,LBj:)
149
150# else
151
152 real(r8), intent(in) :: h(LBi:UBi,LBj:UBj)
153# ifdef MASKING
154 real(r8), intent(in) :: pmask(LBi:UBi,LBj:UBj)
155 real(r8), intent(in) :: rmask(LBi:UBi,LBj:UBj)
156 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
157 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
158# endif
159 real(r8), intent(in) :: zeta(LBi:UBi,LBj:UBj)
160# ifdef SOLVE3D
161 real(r8), intent(in) :: DU_avg1(LBi:UBi,LBj:UBj)
162 real(r8), intent(in) :: DV_avg1(LBi:UBi,LBj:UBj)
163 real(r8), intent(inout) :: rmask_wet_avg(LBi:UBi,LBj:UBj)
164# endif
165 real(r8), intent(out) :: pmask_full(LBi:UBi,LBj:UBj)
166 real(r8), intent(out) :: rmask_full(LBi:UBi,LBj:UBj)
167 real(r8), intent(out) :: umask_full(LBi:UBi,LBj:UBj)
168 real(r8), intent(out) :: vmask_full(LBi:UBi,LBj:UBj)
169
170 real(r8), intent(out) :: pmask_wet(LBi:UBi,LBj:UBj)
171 real(r8), intent(out) :: rmask_wet(LBi:UBi,LBj:UBj)
172 real(r8), intent(out) :: umask_wet(LBi:UBi,LBj:UBj)
173 real(r8), intent(out) :: vmask_wet(LBi:UBi,LBj:UBj)
174# endif
175!
176! Local variable declarations.
177!
178 integer :: i, is, j
179
180 real(r8) :: cff
181 real(r8), parameter :: eps = 1.0e-10_r8
182
183 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: wetdry
184
185# include "set_bounds.h"
186!
187!-----------------------------------------------------------------------
188! If wet/drying, compute new masks for cells with depth < Dcrit.
189!-----------------------------------------------------------------------
190!
191! Compute local mask at RHO-points in terms of free-surface. The local
192! array allows computations in the same parallel region.
193!
194 DO j=jstr-1,jendr
195 DO i=istr-1,iendr
196 wetdry(i,j)=1.0_r8
197# ifdef MASKING
198 wetdry(i,j)=wetdry(i,j)*rmask(i,j)
199# endif
200 IF ((zeta(i,j)+h(i,j)).le.(dcrit(ng)+eps)) THEN
201 wetdry(i,j)=0.0_r8
202 END IF
203 END DO
204 END DO
205!
206! Compute wet/dry mask arrays.
207!
208 IF (iif(ng).le.nfast(ng)) THEN
209 CALL wetdry_mask_tile (ng, tile, &
210 & lbi, ubi, lbj, ubj, &
211 & imins, imaxs, jmins, jmaxs, &
212 & wetdry, &
213 & pmask_wet, rmask_wet, &
214 & umask_wet, vmask_wet)
215 END IF
216
217# ifdef SOLVE3D
218!
219! Wet/dry mask at RHO-points, averaged over all fast time-steps.
220!
221 IF (iif(ng).le.nfast(ng)) THEN
222 IF (predictor_2d_step(ng).and.(first_2d_step)) THEN
223 DO j=jstrr,jendr
224 DO i=istrr,iendr
225 rmask_wet_avg(i,j)=wetdry(i,j)
226 END DO
227 END DO
228 ELSE
229 DO j=jstrr,jendr
230 DO i=istrr,iendr
231 rmask_wet_avg(i,j)=rmask_wet_avg(i,j)+wetdry(i,j)
232 END DO
233 END DO
234 END IF
235!
236 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
237 CALL exchange_r2d_tile (ng, tile, &
238 & lbi, ubi, lbj, ubj, &
239 & rmask_wet_avg)
240 END IF
241
242# ifdef DISTRIBUTE
243!
244 CALL mp_exchange2d (ng, tile, inlm, 1, &
245 & lbi, ubi, lbj, ubj, &
246 & nghostpoints, &
247 & ewperiodic(ng), nsperiodic(ng), &
248 & rmask_wet_avg)
249# endif
250!
251! If done witn 2D fast time-stepping, scale time-averaged mask by
252! the inverse of fast time-step (2*nfast).
253!
254 ELSE
255 cff=1.0_r8/real(2*nfast(ng),r8)
256 DO j=jstr-1,jendr ! executed in a different
257 DO i=istr-1,iendr ! paralllel region
258 wetdry(i,j)=aint(rmask_wet_avg(i,j)*cff)
259 END DO
260 END DO
261!
262! Compute time averaged wet/dry mask arrays.
263!
264 CALL wetdry_avg_mask_tile (ng, tile, &
265 & lbi, ubi, lbj, ubj, &
266 & imins, imaxs, jmins, jmaxs, &
267 & du_avg1, dv_avg1, &
268 & wetdry, &
269 & pmask_wet, rmask_wet, &
270 & umask_wet, vmask_wet)
271 END IF
272# endif
273!
274!-----------------------------------------------------------------------
275! Set masks full time-dependent masks.
276!-----------------------------------------------------------------------
277!
278# ifdef SOLVE3D
279 IF (iif(ng).gt.nfast(ng)) THEN
280# else
281 IF (iif(ng).eq.nfast(ng)) THEN
282# endif
283 DO j=jstrr,jendr
284 DO i=istrr,iendr
285 rmask_full(i,j)=rmask_wet(i,j)*rmask(i,j)
286 END DO
287 END DO
288 DO j=jstr,jendr
289 DO i=istr,iendr
290 pmask_full(i,j)=max(pmask_wet(i,j)*pmask(i,j), 2.0_r8)
291 END DO
292 END DO
293 DO j=jstrr,jendr
294 DO i=istr,iendr
295 umask_full(i,j)=umask_wet(i,j)*umask(i,j)
296 END DO
297 END DO
298 DO j=jstr,jendr
299 DO i=istrr,iendr
300 vmask_full(i,j)=vmask_wet(i,j)*vmask(i,j)
301 END DO
302 END DO
303!
304! Insure that masks at mass point source locations are set to water
305! to avoid writting output with FillValue at those locations.
306!
307 IF (luvsrc(ng)) THEN
308 DO is=1,nsrc(ng)
309 i=sources(ng)%Isrc(is)
310 j=sources(ng)%Jsrc(is)
311 IF (((istrr.le.i).and.(i.le.iendr)).and. &
312 & ((jstrr.le.j).and.(j.le.jendr))) THEN
313 IF (int(sources(ng)%Dsrc(is)).eq.0) THEN
314 umask_full(i,j)=1.0_r8
315 ELSE IF (int(sources(ng)%Dsrc(is)).eq.1) THEN
316 vmask_full(i,j)=1.0_r8
317 END IF
318 END IF
319 END DO
320 END IF
321!
322! Exchange boundary data.
323!
324 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
325 CALL exchange_p2d_tile (ng, tile, &
326 & lbi, ubi, lbj, ubj, &
327 & pmask_full)
328 CALL exchange_r2d_tile (ng, tile, &
329 & lbi, ubi, lbj, ubj, &
330 & rmask_full)
331 CALL exchange_u2d_tile (ng, tile, &
332 & lbi, ubi, lbj, ubj, &
333 & umask_full)
334 CALL exchange_v2d_tile (ng, tile, &
335 & lbi, ubi, lbj, ubj, &
336 & vmask_full)
337 END IF
338
339# ifdef DISTRIBUTE
340!
341 CALL mp_exchange2d (ng, tile, inlm, 4, &
342 & lbi, ubi, lbj, ubj, &
343 & nghostpoints, &
344 & ewperiodic(ng), nsperiodic(ng), &
345 & pmask_full, rmask_full, &
346 & umask_full, vmask_full)
347# endif
348 END IF
349
350 RETURN
351 END SUBROUTINE wetdry_tile
352!
353!***********************************************************************
354 SUBROUTINE wetdry_ini_tile (ng, tile, &
355 & LBi, UBi, LBj, UBj, &
356 & IminS, ImaxS, JminS, JmaxS, &
357# ifdef MASKING
358 & pmask, rmask, umask, vmask, &
359# endif
360 & h, zeta, &
361# ifdef SOLVE3D
362 & ubar, vbar, &
363# endif
364 & pmask_wet, pmask_full, &
365 & rmask_wet, rmask_full, &
366 & umask_wet, umask_full, &
367 & vmask_wet, vmask_full)
368!***********************************************************************
369!
370 USE mod_param
371 USE mod_ncparam
372 USE mod_scalars
373 USE mod_sources
374!
376# ifdef DISTRIBUTE
378# endif
379!
380! Imported variable declarations.
381!
382 integer, intent(in) :: ng, tile
383 integer, intent(in) :: LBi, UBi, LBj, UBj
384 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
385!
386# ifdef ASSUMED_SHAPE
387 real(r8), intent(in) :: h(LBi:,LBj:)
388# ifdef MASKING
389 real(r8), intent(in) :: pmask(LBi:,LBj:)
390 real(r8), intent(in) :: rmask(LBi:,LBj:)
391 real(r8), intent(in) :: umask(LBi:,LBj:)
392 real(r8), intent(in) :: vmask(LBi:,LBj:)
393# endif
394 real(r8), intent(in) :: zeta(LBi:,LBj:)
395# ifdef SOLVE3D
396 real(r8), intent(in) :: ubar(LBi:,LBj:)
397 real(r8), intent(in) :: vbar(LBi:,LBj:)
398# endif
399 real(r8), intent(out) :: pmask_full(LBi:,LBj:)
400 real(r8), intent(out) :: rmask_full(LBi:,LBj:)
401 real(r8), intent(out) :: umask_full(LBi:,LBj:)
402 real(r8), intent(out) :: vmask_full(LBi:,LBj:)
403
404 real(r8), intent(out) :: pmask_wet(LBi:,LBj:)
405 real(r8), intent(out) :: rmask_wet(LBi:,LBj:)
406 real(r8), intent(out) :: umask_wet(LBi:,LBj:)
407 real(r8), intent(out) :: vmask_wet(LBi:,LBj:)
408
409# else
410
411 real(r8), intent(in) :: h(LBi:UBi,LBj:UBj)
412# ifdef MASKING
413 real(r8), intent(in) :: pmask(LBi:UBi,LBj:UBj)
414 real(r8), intent(in) :: rmask(LBi:UBi,LBj:UBj)
415 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
416 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
417# endif
418 real(r8), intent(in) :: zeta(LBi:UBi,LBj:UBj)
419# ifdef SOLVE3D
420 real(r8), intent(in) :: ubar(LBi:UBi,LBj:UBj)
421 real(r8), intent(in) :: vbar(LBi:UBi,LBj:UBj)
422# endif
423 real(r8), intent(out) :: pmask_full(LBi:UBi,LBj:UBj)
424 real(r8), intent(out) :: rmask_full(LBi:UBi,LBj:UBj)
425 real(r8), intent(out) :: umask_full(LBi:UBi,LBj:UBj)
426 real(r8), intent(out) :: vmask_full(LBi:UBi,LBj:UBj)
427
428 real(r8), intent(out) :: pmask_wet(LBi:UBi,LBj:UBj)
429 real(r8), intent(out) :: rmask_wet(LBi:UBi,LBj:UBj)
430 real(r8), intent(out) :: umask_wet(LBi:UBi,LBj:UBj)
431 real(r8), intent(out) :: vmask_wet(LBi:UBi,LBj:UBj)
432# endif
433!
434! Local variable declarations.
435!
436 integer :: i, is, j
437
438 real(r8) :: cff
439 real(r8), parameter :: eps = 1.0e-10_r8
440
441 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: wetdry
442
443# include "set_bounds.h"
444!
445!-----------------------------------------------------------------------
446! If wet/drying, compute new masks for cells with depth < Dcrit.
447!-----------------------------------------------------------------------
448!
449! Compute local mask at RHO-points in terms of free-surface. The local
450! array allows computations in the same parallel region.
451!
452 DO j=jstr-1,jendr
453 DO i=istr-1,iendr
454 wetdry(i,j)=1.0_r8
455# ifdef MASKING
456 wetdry(i,j)=wetdry(i,j)*rmask(i,j)
457# endif
458 IF ((zeta(i,j)+h(i,j)).le.(dcrit(ng)+eps)) THEN
459 wetdry(i,j)=0.0_r8
460 END IF
461 END DO
462 END DO
463!
464! Compute initial wet/dry mask arrays.
465!
466# ifdef SOLVE3D
467 CALL wetdry_avg_mask_tile (ng, tile, &
468 & lbi, ubi, lbj, ubj, &
469 & imins, imaxs, jmins, jmaxs, &
470 & ubar, vbar, &
471 & wetdry, &
472 & pmask_wet, rmask_wet, &
473 & umask_wet, vmask_wet)
474# else
475 CALL wetdry_mask_tile (ng, tile, &
476 & lbi, ubi, lbj, ubj, &
477 & imins, imaxs, jmins, jmaxs, &
478 & wetdry, &
479 & pmask_wet, rmask_wet, &
480 & umask_wet, vmask_wet)
481# endif
482!
483!-----------------------------------------------------------------------
484! Set masks full time-dependent masks.
485!-----------------------------------------------------------------------
486!
487 DO j=jstrr,jendr
488 DO i=istrr,iendr
489 rmask_full(i,j)=rmask_wet(i,j)*rmask(i,j)
490 END DO
491 END DO
492 DO j=jstr,jendr
493 DO i=istr,iendr
494 pmask_full(i,j)=max(pmask_wet(i,j)*pmask(i,j), 2.0_r8)
495 END DO
496 END DO
497 DO j=jstrr,jendr
498 DO i=istr,iendr
499 umask_full(i,j)=umask_wet(i,j)*umask(i,j)
500 END DO
501 END DO
502 DO j=jstr,jendr
503 DO i=istrr,iendr
504 vmask_full(i,j)=vmask_wet(i,j)*vmask(i,j)
505 END DO
506 END DO
507!
508! Insure that masks at mass point source locations are set to water
509! to avoid writting output with FillValue at those locations.
510!
511 IF (luvsrc(ng)) THEN
512 DO is=1,nsrc(ng)
513 i=sources(ng)%Isrc(is)
514 j=sources(ng)%Jsrc(is)
515 IF (((istrr.le.i).and.(i.le.iendr)).and. &
516 & ((jstrr.le.j).and.(j.le.jendr))) THEN
517 IF (int(sources(ng)%Dsrc(is)).eq.0) THEN
518 umask_full(i,j)=1.0_r8
519 ELSE IF (int(sources(ng)%Dsrc(is)).eq.1) THEN
520 vmask_full(i,j)=1.0_r8
521 END IF
522 END IF
523 END DO
524 END IF
525!
526! Exchange boundary data.
527!
528 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
529 CALL exchange_p2d_tile (ng, tile, &
530 & lbi, ubi, lbj, ubj, &
531 & pmask_full)
532 CALL exchange_r2d_tile (ng, tile, &
533 & lbi, ubi, lbj, ubj, &
534 & rmask_full)
535 CALL exchange_u2d_tile (ng, tile, &
536 & lbi, ubi, lbj, ubj, &
537 & umask_full)
538 CALL exchange_v2d_tile (ng, tile, &
539 & lbi, ubi, lbj, ubj, &
540 & vmask_full)
541 END IF
542
543# ifdef DISTRIBUTE
544!
545 CALL mp_exchange2d (ng, tile, inlm, 4, &
546 & lbi, ubi, lbj, ubj, &
547 & nghostpoints, &
548 & ewperiodic(ng), nsperiodic(ng), &
549 & pmask_full, rmask_full, &
550 & umask_full, vmask_full)
551# endif
552
553 RETURN
554 END SUBROUTINE wetdry_ini_tile
555!
556!***********************************************************************
557 SUBROUTINE wetdry_mask_tile (ng, tile, &
558 & LBi, UBi, LBj, UBj, &
559 & IminS, ImaxS, JminS, JmaxS, &
560 & wetdry, &
561 & pmask_wet, rmask_wet, &
562 & umask_wet, vmask_wet)
563!***********************************************************************
564!
565 USE mod_param
566 USE mod_scalars
567!
569# ifdef DISTRIBUTE
571# endif
572!
573! Imported variable declarations.
574!
575 integer, intent(in) :: ng, tile
576 integer, intent(in) :: LBi, UBi, LBj, UBj
577 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
578!
579# ifdef ASSUMED_SHAPE
580 real(r8), intent(in) :: wetdry(IminS:,JminS:)
581 real(r8), intent(out) :: pmask_wet(LBi:,LBj:)
582 real(r8), intent(out) :: rmask_wet(LBi:,LBj:)
583 real(r8), intent(out) :: umask_wet(LBi:,LBj:)
584 real(r8), intent(out) :: vmask_wet(LBi:,LBj:)
585# else
586 real(r8), intent(in) :: wetdry(IminS:ImaxS,JminS:JmaxS)
587 real(r8), intent(out) :: pmask_wet(LBi:UBi,LBj:UBj)
588 real(r8), intent(out) :: rmask_wet(LBi:UBi,LBj:UBj)
589 real(r8), intent(out) :: umask_wet(LBi:UBi,LBj:UBj)
590 real(r8), intent(out) :: vmask_wet(LBi:UBi,LBj:UBj)
591# endif
592!
593! Local variable declarations.
594!
595 integer :: i, j
596
597 real(r8) :: cff1, cff2
598
599# include "set_bounds.h"
600!
601!-----------------------------------------------------------------------
602! Compute wet/dry masks in terms of RHO-points mask (wetdry).
603!-----------------------------------------------------------------------
604!
605! Wet/dry mask at RHO-points.
606!
607 DO j=jstrr,jendr
608 DO i=istrr,iendr
609 rmask_wet(i,j)=wetdry(i,j)
610 END DO
611 END DO
612!
613! Wet/dry mask at U-points.
614!
615 DO j=jstrr,jendr
616 DO i=istr,iendr
617 umask_wet(i,j)=wetdry(i-1,j)+wetdry(i,j)
618 IF (umask_wet(i,j).eq.1.0_r8) THEN
619 umask_wet(i,j)=wetdry(i-1,j)-wetdry(i,j)
620 END IF
621 END DO
622 END DO
623!
624! Wet/dry mask at V-points.
625!
626 DO j=jstr,jendr
627 DO i=istrr,iendr
628 vmask_wet(i,j)=wetdry(i,j-1)+wetdry(i,j)
629 IF (vmask_wet(i,j).eq.1.0_r8) THEN
630 vmask_wet(i,j)=wetdry(i,j-1)-wetdry(i,j)
631 END IF
632 END DO
633 END DO
634!
635! Wet/dry mask at PSI-points.
636!
637 cff1=1.0_r8
638 cff2=2.0_r8
639 DO j=jstr,jendr
640 DO i=istr,iendr
641 IF ((wetdry(i-1,j ).gt.0.5_r8).and. &
642 & (wetdry(i ,j ).gt.0.5_r8).and. &
643 & (wetdry(i-1,j-1).gt.0.5_r8).and. &
644 & (wetdry(i ,j-1).gt.0.5_r8)) THEN
645 pmask_wet(i,j)=1.0_r8
646 ELSE IF ((wetdry(i-1,j ).lt.0.5_r8).and. &
647 & (wetdry(i ,j ).gt.0.5_r8).and. &
648 & (wetdry(i-1,j-1).gt.0.5_r8).and. &
649 & (wetdry(i ,j-1).gt.0.5_r8)) THEN
650 pmask_wet(i,j)=cff1
651 ELSE IF ((wetdry(i-1,j ).gt.0.5_r8).and. &
652 & (wetdry(i ,j ).lt.0.5_r8).and. &
653 & (wetdry(i-1,j-1).gt.0.5_r8).and. &
654 & (wetdry(i ,j-1).gt.0.5_r8)) THEN
655 pmask_wet(i,j)=cff1
656 ELSE IF ((wetdry(i-1,j ).gt.0.5_r8).and. &
657 & (wetdry(i ,j ).gt.0.5_r8).and. &
658 & (wetdry(i-1,j-1).lt.0.5_r8).and. &
659 & (wetdry(i ,j-1).gt.0.5_r8)) THEN
660 pmask_wet(i,j)=cff1
661 ELSE IF ((wetdry(i-1,j ).gt.0.5_r8).and. &
662 & (wetdry(i ,j ).gt.0.5_r8).and. &
663 & (wetdry(i-1,j-1).gt.0.5_r8).and. &
664 & (wetdry(i ,j-1).lt.0.5_r8)) THEN
665 pmask_wet(i,j)=cff1
666 ELSE IF ((wetdry(i-1,j ).gt.0.5_r8).and. &
667 & (wetdry(i ,j ).lt.0.5_r8).and. &
668 & (wetdry(i-1,j-1).gt.0.5_r8).and. &
669 & (wetdry(i ,j-1).lt.0.5_r8)) THEN
670 pmask_wet(i,j)=cff2
671 ELSE IF ((wetdry(i-1,j ).lt.0.5_r8).and. &
672 & (wetdry(i ,j ).gt.0.5_r8).and. &
673 & (wetdry(i-1,j-1).lt.0.5_r8).and. &
674 & (wetdry(i ,j-1).gt.0.5_r8)) THEN
675 pmask_wet(i,j)=cff2
676 ELSE IF ((wetdry(i-1,j ).gt.0.5_r8).and. &
677 & (wetdry(i ,j ).gt.0.5_r8).and. &
678 & (wetdry(i-1,j-1).lt.0.5_r8).and. &
679 & (wetdry(i ,j-1).lt.0.5_r8)) THEN
680 pmask_wet(i,j)=cff2
681 ELSE IF ((wetdry(i-1,j ).lt.0.5_r8).and. &
682 & (wetdry(i ,j ).lt.0.5_r8).and. &
683 & (wetdry(i-1,j-1).gt.0.5_r8).and. &
684 & (wetdry(i ,j-1).gt.0.5_r8)) THEN
685 pmask_wet(i,j)=cff2
686 ELSE
687 pmask_wet(i,j)=0.0_r8
688 END IF
689 END DO
690 END DO
691!
692! Exchange boundary data.
693!
694 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
695 CALL exchange_p2d_tile (ng, tile, &
696 & lbi, ubi, lbj, ubj, &
697 & pmask_wet)
698 CALL exchange_r2d_tile (ng, tile, &
699 & lbi, ubi, lbj, ubj, &
700 & rmask_wet)
701 CALL exchange_u2d_tile (ng, tile, &
702 & lbi, ubi, lbj, ubj, &
703 & umask_wet)
704 CALL exchange_v2d_tile (ng, tile, &
705 & lbi, ubi, lbj, ubj, &
706 & vmask_wet)
707 END IF
708
709# ifdef DISTRIBUTE
710!
711 CALL mp_exchange2d (ng, tile, inlm, 4, &
712 & lbi, ubi, lbj, ubj, &
713 & nghostpoints, &
714 & ewperiodic(ng), nsperiodic(ng), &
715 & pmask_wet, rmask_wet, umask_wet, vmask_wet)
716# endif
717
718 RETURN
719 END SUBROUTINE wetdry_mask_tile
720
721# ifdef SOLVE3D
722!
723!***********************************************************************
724 SUBROUTINE wetdry_avg_mask_tile (ng, tile, &
725 & LBi, UBi, LBj, UBj, &
726 & IminS, ImaxS, JminS, JmaxS, &
727 & DU_avg1, DV_avg1, &
728 & wetdry, &
729 & pmask_wet, rmask_wet, &
730 & umask_wet, vmask_wet)
731!***********************************************************************
732!
733 USE mod_param
734 USE mod_scalars
735!
737# ifdef DISTRIBUTE
739# endif
740!
741! Imported variable declarations.
742!
743 integer, intent(in) :: ng, tile
744 integer, intent(in) :: LBi, UBi, LBj, UBj
745 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
746!
747# ifdef ASSUMED_SHAPE
748 real(r8), intent(in) :: wetdry(IminS:,JminS:)
749 real(r8), intent(in) :: DU_avg1(LBi:,LBj:)
750 real(r8), intent(in) :: DV_avg1(LBi:,LBj:)
751 real(r8), intent(out) :: pmask_wet(LBi:,LBj:)
752 real(r8), intent(out) :: rmask_wet(LBi:,LBj:)
753 real(r8), intent(out) :: umask_wet(LBi:,LBj:)
754 real(r8), intent(out) :: vmask_wet(LBi:,LBj:)
755# else
756 real(r8), intent(in) :: wetdry(IminS:ImaxS,JminS:JmaxS)
757 real(r8), intent(in) :: DU_avg1(LBi:UBi,LBj:UBj)
758 real(r8), intent(in) :: DV_avg1(LBi:UBi,LBj:UBj)
759 real(r8), intent(out) :: pmask_wet(LBi:UBi,LBj:UBj)
760 real(r8), intent(out) :: rmask_wet(LBi:UBi,LBj:UBj)
761 real(r8), intent(out) :: umask_wet(LBi:UBi,LBj:UBj)
762 real(r8), intent(out) :: vmask_wet(LBi:UBi,LBj:UBj)
763# endif
764!
765! Local variable declarations.
766!
767 integer :: i, j
768
769 real(r8) :: cff1, cff2, cff5, cff6
770
771# include "set_bounds.h"
772!
773!-----------------------------------------------------------------------
774! Compute time-averaged wet/dry masks in terms of RHO-points mask
775! (wetdry).
776!-----------------------------------------------------------------------
777!
778! Wet/dry mask at RHO-points, averaged over all fast time-steps.
779!
780 DO j=jstrr,jendr
781 DO i=istrr,iendr
782 rmask_wet(i,j)=wetdry(i,j)
783 END DO
784 END DO
785!
786! Wet/dry mask at U-points, averaged over all fast time-steps.
787!
788 DO j=jstrr,jendr
789 DO i=istr,iendr
790 cff1=wetdry(i-1,j)+wetdry(i,j)
791 IF (cff1.eq.1.0_r8) THEN
792 cff1=wetdry(i-1,j)-wetdry(i,j)
793 END IF
794 cff5=abs(abs(cff1)-1.0_r8)
795 cff6=0.5_r8+dsign(0.5_r8,du_avg1(i,j))*cff1
796 umask_wet(i,j)=0.5_r8*cff1*cff5+cff6*(1.0_r8-cff5)
797! catch lone ponds
798 IF (du_avg1(i,j).eq.0.0_r8) THEN
799 IF ((wetdry(i-1,j)+wetdry(i,j)).le.1.0_r8) THEN
800 umask_wet(i,j)=0.0_r8
801 END IF
802 END IF
803 END DO
804 END DO
805!
806! Wet/dry mask at V-points, averaged over all fast time-steps.
807!
808 DO j=jstr,jendr
809 DO i=istrr,iendr
810 cff1=wetdry(i,j-1)+wetdry(i,j)
811 IF (cff1.eq.1.0_r8) THEN
812 cff1=wetdry(i,j-1)-wetdry(i,j)
813 END IF
814 cff5=abs(abs(cff1)-1.0_r8)
815 cff6=0.5_r8+dsign(0.5_r8,dv_avg1(i,j))*cff1
816 vmask_wet(i,j)=0.5_r8*cff1*cff5+cff6*(1.0_r8-cff5)
817! catch lone ponds
818 IF (dv_avg1(i,j).eq.0.0_r8) THEN
819 IF ((wetdry(i,j-1)+wetdry(i,j)).le.1.0_r8) THEN
820 vmask_wet(i,j)=0.0_r8
821 END IF
822 END IF
823 END DO
824 END DO
825!
826! Wet/dry mask at PSI-points.
827!
828 cff1=1.0_r8
829 cff2=2.0_r8
830 DO j=jstr,jendr
831 DO i=istr,iendr
832 IF ((wetdry(i-1,j ).gt.0.5_r8).and. &
833 & (wetdry(i ,j ).gt.0.5_r8).and. &
834 & (wetdry(i-1,j-1).gt.0.5_r8).and. &
835 & (wetdry(i ,j-1).gt.0.5_r8)) THEN
836 pmask_wet(i,j)=1.0_r8
837 ELSE IF ((wetdry(i-1,j ).lt.0.5_r8).and. &
838 & (wetdry(i ,j ).gt.0.5_r8).and. &
839 & (wetdry(i-1,j-1).gt.0.5_r8).and. &
840 & (wetdry(i ,j-1).gt.0.5_r8)) THEN
841 pmask_wet(i,j)=cff1
842 ELSE IF ((wetdry(i-1,j ).gt.0.5_r8).and. &
843 & (wetdry(i ,j ).lt.0.5_r8).and. &
844 & (wetdry(i-1,j-1).gt.0.5_r8).and. &
845 & (wetdry(i ,j-1).gt.0.5_r8)) THEN
846 pmask_wet(i,j)=cff1
847 ELSE IF ((wetdry(i-1,j ).gt.0.5_r8).and. &
848 & (wetdry(i ,j ).gt.0.5_r8).and. &
849 & (wetdry(i-1,j-1).lt.0.5_r8).and. &
850 & (wetdry(i ,j-1).gt.0.5_r8)) THEN
851 pmask_wet(i,j)=cff1
852 ELSE IF ((wetdry(i-1,j ).gt.0.5_r8).and. &
853 & (wetdry(i ,j ).gt.0.5_r8).and. &
854 & (wetdry(i-1,j-1).gt.0.5_r8).and. &
855 & (wetdry(i ,j-1).lt.0.5_r8)) THEN
856 pmask_wet(i,j)=cff1
857 ELSE IF ((wetdry(i-1,j ).gt.0.5_r8).and. &
858 & (wetdry(i ,j ).lt.0.5_r8).and. &
859 & (wetdry(i-1,j-1).gt.0.5_r8).and. &
860 & (wetdry(i ,j-1).lt.0.5_r8)) THEN
861 pmask_wet(i,j)=cff2
862 ELSE IF ((wetdry(i-1,j ).lt.0.5_r8).and. &
863 & (wetdry(i ,j ).gt.0.5_r8).and. &
864 & (wetdry(i-1,j-1).lt.0.5_r8).and. &
865 & (wetdry(i ,j-1).gt.0.5_r8)) THEN
866 pmask_wet(i,j)=cff2
867 ELSE IF ((wetdry(i-1,j ).gt.0.5_r8).and. &
868 & (wetdry(i ,j ).gt.0.5_r8).and. &
869 & (wetdry(i-1,j-1).lt.0.5_r8).and. &
870 & (wetdry(i ,j-1).lt.0.5_r8)) THEN
871 pmask_wet(i,j)=cff2
872 ELSE IF ((wetdry(i-1,j ).lt.0.5_r8).and. &
873 & (wetdry(i ,j ).lt.0.5_r8).and. &
874 & (wetdry(i-1,j-1).gt.0.5_r8).and. &
875 & (wetdry(i ,j-1).gt.0.5_r8)) THEN
876 pmask_wet(i,j)=cff2
877 ELSE
878 pmask_wet(i,j)=0.0_r8
879 END IF
880 END DO
881 END DO
882!
883! Exchange boundary data.
884!
885 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
886 CALL exchange_p2d_tile (ng, tile, &
887 & lbi, ubi, lbj, ubj, &
888 & pmask_wet)
889 CALL exchange_r2d_tile (ng, tile, &
890 & lbi, ubi, lbj, ubj, &
891 & rmask_wet)
892 CALL exchange_u2d_tile (ng, tile, &
893 & lbi, ubi, lbj, ubj, &
894 & umask_wet)
895 CALL exchange_v2d_tile (ng, tile, &
896 & lbi, ubi, lbj, ubj, &
897 & vmask_wet)
898 END IF
899
900# ifdef DISTRIBUTE
901!
902 CALL mp_exchange2d (ng, tile, inlm, 4, &
903 & lbi, ubi, lbj, ubj, &
904 & nghostpoints, &
905 & ewperiodic(ng), nsperiodic(ng), &
906 & pmask_wet, rmask_wet, umask_wet, vmask_wet)
907# endif
908
909 RETURN
910 END SUBROUTINE wetdry_avg_mask_tile
911# endif
912#endif
913 END MODULE wetdry_mod
subroutine exchange_r2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
subroutine exchange_u2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
subroutine exchange_p2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
Definition exchange_2d.F:66
subroutine exchange_v2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
type(t_coupling), dimension(:), allocatable coupling
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
integer, parameter inlm
Definition mod_param.F:662
integer nghostpoints
Definition mod_param.F:710
logical, dimension(:), allocatable luvsrc
real(r8), dimension(:), allocatable dcrit
logical, dimension(:), allocatable ewperiodic
logical, dimension(:), allocatable nsperiodic
logical, dimension(:), allocatable predictor_2d_step
integer, dimension(:), allocatable nfast
integer, dimension(:), allocatable iif
type(t_sources), dimension(:), allocatable sources
Definition mod_sources.F:90
integer, dimension(:), allocatable nsrc
Definition mod_sources.F:97
subroutine mp_exchange2d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, nghost, ew_periodic, ns_periodic, a, b, c, d)
subroutine wetdry(ng, tile, tindex, linitialize)
Definition wetdry.F:22
subroutine wetdry_ini_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, pmask, rmask, umask, vmask, h, zeta, ubar, vbar, pmask_wet, pmask_full, rmask_wet, rmask_full, umask_wet, umask_full, vmask_wet, vmask_full)
Definition wetdry.F:368
subroutine wetdry_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, pmask, rmask, umask, vmask, h, zeta, du_avg1, dv_avg1, rmask_wet_avg, pmask_wet, pmask_full, rmask_wet, rmask_full, umask_wet, umask_full, vmask_wet, vmask_full)
Definition wetdry.F:108
subroutine wetdry_mask_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, wetdry, pmask_wet, rmask_wet, umask_wet, vmask_wet)
Definition wetdry.F:563
subroutine wetdry_avg_mask_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, du_avg1, dv_avg1, wetdry, pmask_wet, rmask_wet, umask_wet, vmask_wet)
Definition wetdry.F:731