ROMS
Loading...
Searching...
No Matches
ana_mask.h
Go to the documentation of this file.
1!!
2 SUBROUTINE ana_mask (ng, tile, model)
3!
4!! git $Id$
5!!======================================================================
6!! Copyright (c) 2002-2025 The ROMS Group !
7!! Licensed under a MIT/X style license !
8!! See License_ROMS.md !
9!=======================================================================
10! !
11! This subroutine sets analytical Land/Sea masking. !
12! !
13!=======================================================================
14!
15 USE mod_param
16 USE mod_grid
17 USE mod_ncparam
18!
19! Imported variable declarations.
20!
21 integer, intent(in) :: ng, tile, model
22!
23! Local variable declarations.
24!
25 character (len=*), parameter :: MyFile = &
26 & __FILE__
27!
28#include "tile.h"
29!
30 CALL ana_mask_tile (ng, tile, model, &
31 & lbi, ubi, lbj, ubj, &
32 & imins, imaxs, jmins, jmaxs, &
33 & grid(ng) % pmask, &
34 & grid(ng) % rmask, &
35 & grid(ng) % umask, &
36 & grid(ng) % vmask)
37!
38! Set analytical header file name used.
39!
40#ifdef DISTRIBUTE
41 IF (lanafile) THEN
42#else
43 IF (lanafile.and.(tile.eq.0)) THEN
44#endif
45 ananame(15)=myfile
46 END IF
47!
48 RETURN
49 END SUBROUTINE ana_mask
50!
51!***********************************************************************
52 SUBROUTINE ana_mask_tile (ng, tile, model, &
53 & LBi, UBi, LBj, UBj, &
54 & IminS, ImaxS, JminS, JmaxS, &
55 & pmask, rmask, umask, vmask)
56!***********************************************************************
57!
58 USE mod_param
59 USE mod_parallel
60 USE mod_ncparam
61 USE mod_iounits
62 USE mod_scalars
63!
65#ifdef DISTRIBUTE
67#endif
68 USE stats_mod, ONLY : stats_2dfld
69!
70! Imported variable declarations.
71!
72 integer, intent(in) :: ng, tile, model
73 integer, intent(in) :: LBi, UBi, LBj, UBj
74 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
75!
76#ifdef ASSUMED_SHAPE
77 real(r8), intent(out) :: pmask(LBi:,LBj:)
78 real(r8), intent(out) :: rmask(LBi:,LBj:)
79 real(r8), intent(out) :: umask(LBi:,LBj:)
80 real(r8), intent(out) :: vmask(LBi:,LBj:)
81#else
82 real(r8), intent(out) :: pmask(LBi:UBi,LBj:UBj)
83 real(r8), intent(out) :: rmask(LBi:UBi,LBj:UBj)
84 real(r8), intent(out) :: umask(LBi:UBi,LBj:UBj)
85 real(r8), intent(out) :: vmask(LBi:UBi,LBj:UBj)
86#endif
87!
88! Local variable declarations.
89!
90 logical, save :: first = .true.
91!
92 integer :: Imin, Imax, Jmin, Jmax
93 integer :: i, j
94!
95 real(r8) :: mask(IminS:ImaxS,JminS:JmaxS)
96!
97 TYPE (T_STATS), save :: Stats(4)
98
99#include "set_bounds.h"
100!
101!-----------------------------------------------------------------------
102! Initialize field statictics structure.
103!-----------------------------------------------------------------------
104!
105 IF (first) THEN
106 first=.false.
107 DO i=1,SIZE(stats,1)
108 stats(i) % checksum=0_i8b
109 stats(i) % count=0
110 stats(i) % min=large
111 stats(i) % max=-large
112 stats(i) % avg=0.0_r8
113 stats(i) % rms=0.0_r8
114 END DO
115 END IF
116!
117!-----------------------------------------------------------------------
118! Set Land/Sea mask of RHO-points: Land=0, Sea=1.
119!-----------------------------------------------------------------------
120!
121! Notice that private scratch array "mask" is used to allow
122! computation within a parallel loop.
123!
124#ifdef DOUBLE_GYRE
125 imin=-2+(lm(ng)+1)/2
126 imax=imin+2
127 jmin=-2+(mm(ng)+1)/2
128 jmax=jmin+2
129 DO j=jstrm2,jendp2
130 DO i=istrm2,iendp2
131 mask(i,j)=1.0_r8
132 IF (((imin.le.i).and.(i.le.imax)).and. &
133 & ((jmin.le.j).and.(j.le.jmax))) THEN
134 mask(i,j)=0.0_r8
135 END IF
136 END DO
137 END DO
138#elif defined FLT_TEST
139 DO j=jstrm2,jendp2
140 DO i=istrm2,iendp2
141 mask(i,j)=1.0_r8
142 IF (j.eq.1 ) mask(i,j)=0.0_r8
143 IF (j.eq.mm(ng)) mask(i,j)=0.0_r8
144 IF ((i.ge.((lm(ng)+1)/2)).and. &
145 & (i.le.((lm(ng)+1)/2+1)).and. &
146 & (j.ge.((mm(ng)+1)/2)).and. &
147 & (j.le.((mm(ng)+1)/2+1))) mask(i,j)=0.0_r8
148 END DO
149 END DO
150#elif defined LAKE_SIGNELL
151 DO j=jstrm2,jendp2
152 DO i=istrm2,iendp2
153 mask(i,j)=1.0_r8
154 END DO
155 END DO
156 IF (domain(ng)%Western_Edge(tile)) THEN
157 DO j=jstrm1,jendp1
158 mask(istr-1,j)=0.0_r8
159 END DO
160 END IF
161 IF (domain(ng)%Eastern_Edge(tile)) THEN
162 DO j=jstrm1,jendp1
163 mask(iend+1,j)=0.0_r8
164 END DO
165 END IF
166 IF (domain(ng)%Southern_Edge(tile)) THEN
167 DO i=istrm1,iendp1
168 mask(i,jstr-1)=0.0_r8
169 END DO
170 END IF
171 IF (domain(ng)%Northern_Edge(tile)) THEN
172 DO i=istrm1,iendp1
173 mask(i,jend+1)=0.0_r8
174 END DO
175 END IF
176#elif defined RIVERPLUME1
177 DO j=jstrm2,jendp2
178 DO i=istrm2,iendp2
179 mask(i,j)=1.0_r8
180 END DO
181 END DO
182 DO i=istrm2,min(5,iendp2)
183 DO j=jstrm2,min(mm(ng)-18,jendp2)
184 mask(i,j)=0.0_r8
185 END DO
186 DO j=max(jstrm2,mm(ng)-16),jendp2
187 mask(i,j)=0.0_r8
188 END DO
189 END DO
190#elif defined RIVERPLUME2
191 DO j=jstrm2,jendp2
192 DO i=istrm2,iendp2
193 mask(i,j)=1.0_r8
194 END DO
195 END DO
196 DO i=istrm2,min(5,iendp2)
197 DO j=jstrm2,min(mm(ng)-11,jendp2)
198 mask(i,j)=0.0_r8
199 END DO
200 DO j=max(jstrm2,mm(ng)-9),jendp2
201 mask(i,j)=0.0_r8
202 END DO
203 END DO
204#elif defined SHOREFACE
205 DO j=jstrm2,jendp2
206 DO i=istrm2,iendp2
207 mask(i,j)=1.0_r8
208 END DO
209 END DO
210#else
211 ana_mask.h: no values provided for mask.
212#endif
213!
214 DO j=jstrt,jendt
215 DO i=istrt,iendt
216 rmask(i,j)=mask(i,j)
217 END DO
218 END DO
219!
220!-----------------------------------------------------------------------
221! Compute Land/Sea mask of U- and V-points.
222!-----------------------------------------------------------------------
223!
224 DO j=jstrt,jendt
225 DO i=istrp,iendt
226 umask(i,j)=mask(i-1,j)*mask(i,j)
227 END DO
228 END DO
229 DO j=jstrp,jendt
230 DO i=istrt,iendt
231 vmask(i,j)=mask(i,j-1)*mask(i,j)
232 END DO
233 END DO
234!
235!-----------------------------------------------------------------------
236! Compute Land/Sea mask of PSI-points.
237!-----------------------------------------------------------------------
238!
239 DO j=jstrp,jendt
240 DO i=istrp,iendt
241 pmask(i,j)=mask(i-1,j-1)*mask(i,j-1)* &
242 & mask(i-1,j )*mask(i,j )
243 END DO
244 END DO
245!
246!-----------------------------------------------------------------------
247! Report statitics.
248!-----------------------------------------------------------------------
249!
250 CALL stats_2dfld (ng, tile, inlm, p2dvar, stats(1), 0, &
251 & lbi, ubi, lbj, ubj, pmask)
252 IF (domain(ng)%NorthEast_Corner(tile)) THEN
253 WRITE (stdout,10) 'mask on PSI-points: mask_psi', &
254 & ng, stats(1)%min, stats(1)%max
255 END IF
256 CALL stats_2dfld (ng, tile, inlm, r2dvar, stats(2), 0, &
257 & lbi, ubi, lbj, ubj, rmask)
258 IF (domain(ng)%NorthEast_Corner(tile)) THEN
259 WRITE (stdout,10) 'mask on RHO-points: mask_rho', &
260 & ng, stats(2)%min, stats(2)%max
261 END IF
262 CALL stats_2dfld (ng, tile, inlm, u2dvar, stats(3), 0, &
263 & lbi, ubi, lbj, ubj, umask)
264 IF (domain(ng)%NorthEast_Corner(tile)) THEN
265 WRITE (stdout,10) 'mask on U-points: mask_u', &
266 & ng, stats(3)%min, stats(3)%max
267 END IF
268 CALL stats_2dfld (ng, tile, inlm, v2dvar, stats(4), 0, &
269 & lbi, ubi, lbj, ubj, vmask)
270 IF (domain(ng)%NorthEast_Corner(tile)) THEN
271 WRITE (stdout,10) 'mask on V-points: mask_v', &
272 & ng, stats(4)%min, stats(4)%max
273 END IF
274!
275!-----------------------------------------------------------------------
276! Exchange boundary data.
277!-----------------------------------------------------------------------
278!
279 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
280 CALL exchange_r2d_tile (ng, tile, &
281 & lbi, ubi, lbj, ubj, &
282 & rmask)
283 CALL exchange_p2d_tile (ng, tile, &
284 & lbi, ubi, lbj, ubj, &
285 & pmask)
286 CALL exchange_u2d_tile (ng, tile, &
287 & lbi, ubi, lbj, ubj, &
288 & umask)
289 CALL exchange_v2d_tile (ng, tile, &
290 & lbi, ubi, lbj, ubj, &
291 & vmask)
292 END IF
293
294#ifdef DISTRIBUTE
295 CALL mp_exchange2d (ng, tile, model, 4, &
296 & lbi, ubi, lbj, ubj, &
297 & nghostpoints, &
298 & ewperiodic(ng), nsperiodic(ng), &
299 & rmask, pmask, umask, vmask)
300#endif
301!
302 10 FORMAT (3x,' ANA_MASK - ',a,/,19x, &
303 & '(Grid = ',i2.2,', Min = ',1p,e15.8,0p, &
304 & ' Max = ',1p,e15.8,0p,')')
305!
306 RETURN
307 END SUBROUTINE ana_mask_tile
subroutine ana_mask(ng, tile, model)
Definition ana_mask.h:3
subroutine ana_mask_tile(ng, tile, model, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, pmask, rmask, umask, vmask)
Definition ana_mask.h:56
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_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
integer stdout
logical lanafile
character(len=256), dimension(39) ananame
integer, parameter inlm
Definition mod_param.F:662
integer nghostpoints
Definition mod_param.F:710
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
integer, dimension(:), allocatable lm
Definition mod_param.F:455
integer, parameter u2dvar
Definition mod_param.F:718
integer, parameter p2dvar
Definition mod_param.F:716
integer, dimension(:), allocatable mm
Definition mod_param.F:456
integer, parameter r2dvar
Definition mod_param.F:717
integer, parameter v2dvar
Definition mod_param.F:719
logical, dimension(:), allocatable ewperiodic
logical, dimension(:), allocatable nsperiodic
real(dp), parameter large
subroutine mp_exchange2d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, nghost, ew_periodic, ns_periodic, a, b, c, d)
subroutine, public stats_2dfld(ng, tile, model, gtype, s, extract_flag, lbi, ubi, lbj, ubj, f, fmask, debug)
Definition stats.F:47