ROMS
Loading...
Searching...
No Matches
set_masks_mod Module Reference

Functions/Subroutines

subroutine, public set_masks (ng, tile, model)
 
subroutine set_masks_tile (ng, tile, model, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, pmask, rmask, umask, vmask, pmask_avg, rmask_avg, umask_avg, vmask_avg, pmask_dia, rmask_dia, umask_dia, vmask_dia, pmask_full, rmask_full, umask_full, vmask_full)
 
subroutine, public set_avg_masks (ng, tile, model, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, pmask_avg, rmask_avg, umask_avg, vmask_avg)
 

Function/Subroutine Documentation

◆ set_avg_masks()

subroutine, public set_masks_mod::set_avg_masks ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) model,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs,
real(r8), dimension(lbi:,lbj:), intent(inout) pmask_avg,
real(r8), dimension(lbi:,lbj:), intent(inout) rmask_avg,
real(r8), dimension(lbi:,lbj:), intent(inout) umask_avg,
real(r8), dimension(lbi:,lbj:), intent(inout) vmask_avg )

Definition at line 412 of file set_masks.F.

417!***********************************************************************
418!
419 USE mod_param
420 USE mod_scalars
421!
423# ifdef DISTRIBUTE
425# endif
426!
427! Imported variable declarations.
428!
429 integer, intent(in) :: ng, tile, model
430 integer, intent(in) :: LBi, UBi, LBj, UBj
431 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
432!
433# ifdef ASSUMED_SHAPE
434 real(r8), intent(inout) :: pmask_avg(LBi:,LBj:)
435 real(r8), intent(inout) :: rmask_avg(LBi:,LBj:)
436 real(r8), intent(inout) :: umask_avg(LBi:,LBj:)
437 real(r8), intent(inout) :: vmask_avg(LBi:,LBj:)
438# else
439 real(r8), intent(inout) :: pmask_avg(LBi:UBi,LBj:UBj)
440 real(r8), intent(inout) :: rmask_avg(LBi:UBi,LBj:UBj)
441 real(r8), intent(inout) :: umask_avg(LBi:UBi,LBj:UBj)
442 real(r8), intent(inout) :: vmask_avg(LBi:UBi,LBj:UBj)
443# endif
444!
445!
446! Local variable declarations.
447!
448 integer :: i, j
449
450# include "set_bounds.h"
451!
452!-----------------------------------------------------------------------
453! Return if time-averaging window is zero.
454!-----------------------------------------------------------------------
455!
456 IF (navg(ng).eq.0) RETURN
457!
458!-----------------------------------------------------------------------
459! If last time-step of average window, convert time dependent counters
460! for wet points to time-averaged Land/Sea masks (dry=0, wet=1) for
461! the current average window period. Notice that a grid point is wet
462! if the count is greater than zero for the current time average
463! window.
464!-----------------------------------------------------------------------
465!
466 IF ((iic(ng).gt.ntsavg(ng)).and. &
467 & (mod(iic(ng)-1,navg(ng)).eq.0).and. &
468 & ((iic(ng).ne.ntstart(ng)).or.(nrrec(ng).eq.0))) THEN
469
470 DO j=jstrp,jendp
471 DO i=istrp,iendp
472 pmask_avg(i,j)=min(1.0_r8, pmask_avg(i,j))
473 END DO
474 END DO
475 DO j=jstrt,jendt
476 DO i=istrt,iendt
477 rmask_avg(i,j)=min(1.0_r8, rmask_avg(i,j))
478 END DO
479 END DO
480 DO j=jstrt,jendt
481 DO i=istrp,iendt
482 umask_avg(i,j)=min(1.0_r8, umask_avg(i,j))
483 END DO
484 END DO
485 DO j=jstrp,jendt
486 DO i=istrt,iendt
487 vmask_avg(i,j)=min(1.0_r8, vmask_avg(i,j))
488 END DO
489 END DO
490!
491 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
492 CALL exchange_p2d_tile (ng, tile, &
493 & lbi, ubi, lbj, ubj, &
494 & pmask_avg)
495 CALL exchange_r2d_tile (ng, tile, &
496 & lbi, ubi, lbj, ubj, &
497 & rmask_avg)
498 CALL exchange_u2d_tile (ng, tile, &
499 & lbi, ubi, lbj, ubj, &
500 & umask_avg)
501 CALL exchange_v2d_tile (ng, tile, &
502 & lbi, ubi, lbj, ubj, &
503 & vmask_avg)
504 END IF
505
506# ifdef DISTRIBUTE
507 CALL mp_exchange2d (ng, tile, model, 4, &
508 & lbi, ubi, lbj, ubj, &
509 & nghostpoints, &
510 & ewperiodic(ng), nsperiodic(ng), &
511 & pmask_avg, rmask_avg, umask_avg, vmask_avg)
512# endif
513
514 END IF
515!
516 RETURN
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)
integer nghostpoints
Definition mod_param.F:710
integer, dimension(:), allocatable nrrec
integer, dimension(:), allocatable iic
logical, dimension(:), allocatable ewperiodic
logical, dimension(:), allocatable nsperiodic
integer, dimension(:), allocatable navg
integer, dimension(:), allocatable ntstart
integer, dimension(:), allocatable ntsavg
subroutine mp_exchange2d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, nghost, ew_periodic, ns_periodic, a, b, c, d)

References mod_scalars::ewperiodic, exchange_2d_mod::exchange_p2d_tile(), exchange_2d_mod::exchange_r2d_tile(), exchange_2d_mod::exchange_u2d_tile(), exchange_2d_mod::exchange_v2d_tile(), mod_scalars::iic, mp_exchange_mod::mp_exchange2d(), mod_scalars::navg, mod_param::nghostpoints, mod_scalars::nrrec, mod_scalars::nsperiodic, mod_scalars::ntsavg, and mod_scalars::ntstart.

Referenced by set_avg_mod::set_avg().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_masks()

subroutine, public set_masks_mod::set_masks ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) model )

Definition at line 43 of file set_masks.F.

44!***********************************************************************
45!
46 USE mod_param
47 USE mod_grid
48!
49! Imported variable declarations.
50!
51 integer, intent(in) :: ng, tile, model
52!
53! Local variable declarations.
54!
55 character (len=*), parameter :: MyFile = &
56 & __FILE__
57!
58# include "tile.h"
59!
60# ifdef PROFILE
61 CALL wclock_on (ng, model, 2, __line__, myfile)
62# endif
63 CALL set_masks_tile (ng, tile, model, &
64 & lbi, ubi, lbj, ubj, &
65 & imins, imaxs, jmins, jmaxs, &
66 & grid(ng) % pmask, &
67 & grid(ng) % rmask, &
68 & grid(ng) % umask, &
69 & grid(ng) % vmask, &
70# if defined AVERAGES || \
71 (defined ad_averages && defined adjoint) || \
72 (defined rp_averages && defined tl_ioms) || \
73 (defined tl_averages && defined tangent)
74 & grid(ng) % pmask_avg, &
75 & grid(ng) % rmask_avg, &
76 & grid(ng) % umask_avg, &
77 & grid(ng) % vmask_avg, &
78# endif
79# ifdef DIAGNOSTICS
80 & grid(ng) % pmask_dia, &
81 & grid(ng) % rmask_dia, &
82 & grid(ng) % umask_dia, &
83 & grid(ng) % vmask_dia, &
84# endif
85 & grid(ng) % pmask_full, &
86 & grid(ng) % rmask_full, &
87 & grid(ng) % umask_full, &
88 & grid(ng) % vmask_full)
89# ifdef PROFILE
90 CALL wclock_off (ng, model, 2, __line__, myfile)
91# endif
92!
93 RETURN
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
recursive subroutine wclock_off(ng, model, region, line, routine)
Definition timers.F:148
recursive subroutine wclock_on(ng, model, region, line, routine)
Definition timers.F:3

References mod_grid::grid, set_masks_tile(), wclock_off(), and wclock_on().

Referenced by ad_initial(), roms_kernel_mod::adm_initial(), i4dvar_mod::analysis(), initial(), roms_kernel_mod::nlm_initial(), rp_initial(), tl_initial(), and roms_kernel_mod::tlm_initial().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_masks_tile()

subroutine set_masks_mod::set_masks_tile ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) model,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs,
real(r8), dimension(lbi:,lbj:), intent(in) pmask,
real(r8), dimension(lbi:,lbj:), intent(in) rmask,
real(r8), dimension(lbi:,lbj:), intent(in) umask,
real(r8), dimension(lbi:,lbj:), intent(in) vmask,
real(r8), dimension(lbi:,lbj:), intent(inout) pmask_avg,
real(r8), dimension(lbi:,lbj:), intent(inout) rmask_avg,
real(r8), dimension(lbi:,lbj:), intent(inout) umask_avg,
real(r8), dimension(lbi:,lbj:), intent(inout) vmask_avg,
real(r8), dimension(lbi:,lbj:), intent(inout) pmask_dia,
real(r8), dimension(lbi:,lbj:), intent(inout) rmask_dia,
real(r8), dimension(lbi:,lbj:), intent(inout) umask_dia,
real(r8), dimension(lbi:,lbj:), intent(inout) vmask_dia,
real(r8), dimension(lbi:,lbj:), intent(inout) pmask_full,
real(r8), dimension(lbi:,lbj:), intent(inout) rmask_full,
real(r8), dimension(lbi:,lbj:), intent(inout) umask_full,
real(r8), dimension(lbi:,lbj:), intent(inout) vmask_full )
private

Definition at line 97 of file set_masks.F.

115!***********************************************************************
116!
117 USE mod_param
118 USE mod_scalars
119 USE mod_sources
120!
122# ifdef DISTRIBUTE
124# endif
125!
126! Imported variable declarations.
127!
128 integer, intent(in) :: ng, tile, model
129 integer, intent(in) :: LBi, UBi, LBj, UBj
130 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
131!
132# ifdef ASSUMED_SHAPE
133 real(r8), intent(in) :: pmask(LBi:,LBj:)
134 real(r8), intent(in) :: rmask(LBi:,LBj:)
135 real(r8), intent(in) :: umask(LBi:,LBj:)
136 real(r8), intent(in) :: vmask(LBi:,LBj:)
137# if defined AVERAGES || \
138 (defined ad_averages && defined adjoint) || \
139 (defined rp_averages && defined tl_ioms) || \
140 (defined tl_averages && defined tangent)
141 real(r8), intent(inout) :: pmask_avg(LBi:,LBj:)
142 real(r8), intent(inout) :: rmask_avg(LBi:,LBj:)
143 real(r8), intent(inout) :: umask_avg(LBi:,LBj:)
144 real(r8), intent(inout) :: vmask_avg(LBi:,LBj:)
145# endif
146# ifdef DIAGNOSTICS
147 real(r8), intent(inout) :: pmask_dia(LBi:,LBj:)
148 real(r8), intent(inout) :: rmask_dia(LBi:,LBj:)
149 real(r8), intent(inout) :: umask_dia(LBi:,LBj:)
150 real(r8), intent(inout) :: vmask_dia(LBi:,LBj:)
151# endif
152 real(r8), intent(inout) :: pmask_full(LBi:,LBj:)
153 real(r8), intent(inout) :: rmask_full(LBi:,LBj:)
154 real(r8), intent(inout) :: umask_full(LBi:,LBj:)
155 real(r8), intent(inout) :: vmask_full(LBi:,LBj:)
156# else
157 real(r8), intent(in) :: pmask(LBi:UBi,LBj:UBj)
158 real(r8), intent(in) :: rmask(LBi:UBi,LBj:UBj)
159 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
160 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
161# if defined AVERAGES || \
162 (defined ad_averages && defined adjoint) || \
163 (defined rp_averages && defined tl_ioms) || \
164 (defined tl_averages && defined tangent)
165 real(r8), intent(inout) :: pmask_avg(LBi:UBi,LBj:UBj)
166 real(r8), intent(inout) :: rmask_avg(LBi:UBi,LBj:UBj)
167 real(r8), intent(inout) :: umask_avg(LBi:UBi,LBj:UBj)
168 real(r8), intent(inout) :: vmask_avg(LBi:UBi,LBj:UBj)
169# endif
170# ifdef DIAGNOSTICS
171 real(r8), intent(inout) :: pmask_dia(LBi:UBi,LBj:UBj)
172 real(r8), intent(inout) :: rmask_dia(LBi:UBi,LBj:UBj)
173 real(r8), intent(inout) :: umask_dia(LBi:UBi,LBj:UBj)
174 real(r8), intent(inout) :: vmask_dia(LBi:UBi,LBj:UBj)
175# endif
176 real(r8), intent(inout) :: pmask_full(LBi:UBi,LBj:UBj)
177 real(r8), intent(inout) :: rmask_full(LBi:UBi,LBj:UBj)
178 real(r8), intent(inout) :: umask_full(LBi:UBi,LBj:UBj)
179 real(r8), intent(inout) :: vmask_full(LBi:UBi,LBj:UBj)
180# endif
181!
182! Local variable declarations.
183!
184 integer :: i, is, j
185
186# include "set_bounds.h"
187!
188!-----------------------------------------------------------------------
189! Initialize internal full Land/Sea masks with its respective
190! application time-indrpendent values.
191!-----------------------------------------------------------------------
192!
193! The full mask values are updated with time-dependent values in
194! file "wetdry.F" if wetting and drying is activated.
195!
196 DO j=jstrp,jendp
197 DO i=istrp,iendp
198 pmask_full(i,j)=pmask(i,j)
199 END DO
200 END DO
201 DO j=jstrt,jendt
202 DO i=istrt,iendt
203 rmask_full(i,j)=rmask(i,j)
204 END DO
205 END DO
206 DO j=jstrt,jendt
207 DO i=istrp,iendt
208 umask_full(i,j)=umask(i,j)
209 END DO
210 END DO
211 DO j=jstrp,jendt
212 DO i=istrt,iendt
213 vmask_full(i,j)=vmask(i,j)
214 END DO
215 END DO
216!
217! Insure that masks at mass point source locations are set to water
218! to avoid masking with _FillValue at those locations.
219!
220 IF (luvsrc(ng)) THEN
221 DO is=1,nsrc(ng)
222 i=sources(ng)%Isrc(is)
223 j=sources(ng)%Jsrc(is)
224 IF (((istrt.le.i).and.(i.le.iendt)).and. &
225 & ((jstrt.le.j).and.(j.le.jendt))) THEN
226 IF (int(sources(ng)%Dsrc(is)).eq.0) THEN
227 umask_full(i,j)=1.0_r8
228 ELSE
229 vmask_full(i,j)=1.0_r8
230 END IF
231 END IF
232 END DO
233 END IF
234!
235 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
236 CALL exchange_p2d_tile (ng, tile, &
237 & lbi, ubi, lbj, ubj, &
238 & pmask_full)
239 CALL exchange_r2d_tile (ng, tile, &
240 & lbi, ubi, lbj, ubj, &
241 & rmask_full)
242 CALL exchange_u2d_tile (ng, tile, &
243 & lbi, ubi, lbj, ubj, &
244 & umask_full)
245 CALL exchange_v2d_tile (ng, tile, &
246 & lbi, ubi, lbj, ubj, &
247 & vmask_full)
248 END IF
249
250# ifdef DISTRIBUTE
251 CALL mp_exchange2d (ng, tile, model, 4, &
252 & lbi, ubi, lbj, ubj, &
253 & nghostpoints, &
254 & ewperiodic(ng), nsperiodic(ng), &
255 & pmask_full, rmask_full, umask_full, vmask_full)
256# endif
257
258# if defined AVERAGES || \
259 (defined ad_averages && defined adjoint) || \
260 (defined rp_averages && defined tl_ioms) || \
261 (defined tl_averages && defined tangent)
262!
263!-----------------------------------------------------------------------
264! Initialize average file Land/Sea masks for time-averaged fields.
265!-----------------------------------------------------------------------
266!
267 DO j=jstrp,jendp
268 DO i=istrp,iendp
269# ifdef WET_DRY
270 pmask_avg(i,j)=0.0_r8
271# else
272 pmask_avg(i,j)=pmask_full(i,j)
273# endif
274 END DO
275 END DO
276
277 DO j=jstrt,jendt
278 DO i=istrt,iendt
279# ifdef WET_DRY
280 rmask_avg(i,j)=0.0_r8
281# else
282 rmask_avg(i,j)=rmask_full(i,j)
283# endif
284 END DO
285 END DO
286
287 DO j=jstrt,jendt
288 DO i=istrp,iendt
289# ifdef WET_DRY
290 umask_avg(i,j)=0.0_r8
291# else
292 umask_avg(i,j)=umask_full(i,j)
293# endif
294 END DO
295 END DO
296
297 DO j=jstrp,jendt
298 DO i=istrt,iendt
299# ifdef WET_DRY
300 vmask_avg(i,j)=0.0_r8
301# else
302 vmask_avg(i,j)=vmask_full(i,j)
303# endif
304 END DO
305 END DO
306!
307 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
308 CALL exchange_p2d_tile (ng, tile, &
309 & lbi, ubi, lbj, ubj, &
310 & pmask_avg)
311 CALL exchange_r2d_tile (ng, tile, &
312 & lbi, ubi, lbj, ubj, &
313 & rmask_avg)
314 CALL exchange_u2d_tile (ng, tile, &
315 & lbi, ubi, lbj, ubj, &
316 & umask_avg)
317 CALL exchange_v2d_tile (ng, tile, &
318 & lbi, ubi, lbj, ubj, &
319 & vmask_avg)
320 END IF
321
322# ifdef DISTRIBUTE
323 CALL mp_exchange2d (ng, tile, model, 4, &
324 & lbi, ubi, lbj, ubj, &
325 & nghostpoints, &
326 & ewperiodic(ng), nsperiodic(ng), &
327 & pmask_avg, rmask_avg, umask_avg, vmask_avg)
328# endif
329# endif
330
331# ifdef DIAGNOSTICS
332!
333!-----------------------------------------------------------------------
334! Initialize diagnostic file Land/Sea masks for time-averaged fields.
335!-----------------------------------------------------------------------
336!
337 DO j=jstrp,jendp
338 DO i=istrp,iendp
339# ifdef WET_DRY
340 pmask_dia(i,j)=0.0_r8
341# else
342 pmask_dia(i,j)=pmask_full(i,j)
343# endif
344 END DO
345 END DO
346
347 DO j=jstrt,jendt
348 DO i=istrt,iendt
349# ifdef WET_DRY
350 rmask_dia(i,j)=0.0_r8
351# else
352 rmask_dia(i,j)=rmask_full(i,j)
353# endif
354 END DO
355 END DO
356
357 DO j=jstrt,jendt
358 DO i=istrp,iendt
359# ifdef WET_DRY
360 umask_dia(i,j)=0.0_r8
361# else
362 umask_dia(i,j)=umask_full(i,j)
363# endif
364 END DO
365 END DO
366
367 DO j=jstrp,jendt
368 DO i=istrt,iendt
369# ifdef WET_DRY
370 vmask_dia(i,j)=0.0_r8
371# else
372 vmask_dia(i,j)=vmask_full(i,j)
373# endif
374 END DO
375 END DO
376!
377 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
378 CALL exchange_p2d_tile (ng, tile, &
379 & lbi, ubi, lbj, ubj, &
380 & pmask_dia)
381 CALL exchange_r2d_tile (ng, tile, &
382 & lbi, ubi, lbj, ubj, &
383 & rmask_dia)
384 CALL exchange_u2d_tile (ng, tile, &
385 & lbi, ubi, lbj, ubj, &
386 & umask_dia)
387 CALL exchange_v2d_tile (ng, tile, &
388 & lbi, ubi, lbj, ubj, &
389 & vmask_dia)
390 END IF
391
392# ifdef DISTRIBUTE
393 CALL mp_exchange2d (ng, tile, model, 4, &
394 & lbi, ubi, lbj, ubj, &
395 & nghostpoints, &
396 & ewperiodic(ng), nsperiodic(ng), &
397 & pmask_dia, rmask_dia, umask_dia, vmask_dia)
398# endif
399
400# endif
401
402 RETURN
logical, dimension(:), allocatable luvsrc
type(t_sources), dimension(:), allocatable sources
Definition mod_sources.F:90
integer, dimension(:), allocatable nsrc
Definition mod_sources.F:97

References mod_scalars::ewperiodic, exchange_2d_mod::exchange_p2d_tile(), exchange_2d_mod::exchange_r2d_tile(), exchange_2d_mod::exchange_u2d_tile(), exchange_2d_mod::exchange_v2d_tile(), mod_scalars::luvsrc, mp_exchange_mod::mp_exchange2d(), mod_param::nghostpoints, mod_scalars::nsperiodic, mod_sources::nsrc, and mod_sources::sources.

Referenced by set_masks().

Here is the call graph for this function:
Here is the caller graph for this function: