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

Functions/Subroutines

subroutine, public get_bounds (ng, tile, gtype, nghost, my_im, my_jm, my_lm, my_mm, itile, jtile, lbi, ubi, lbj, ubj)
 
subroutine, public get_domain (ng, tile, gtype, nghost, my_im, my_jm, my_lm, my_mm, epsilon, lfullgrid, xmin, xmax, ymin, ymax)
 
subroutine, public get_domain_edges (ng, tile, my_lm, my_mm, eastern_edge, western_edge, northern_edge, southern_edge, northeast_corner, northwest_corner, southeast_corner, southwest_corner, northeast_test, northwest_test, southeast_test, southwest_test)
 
subroutine, public get_iobounds (ng, my_lm, my_mm, my_bounds, my_iobounds)
 
subroutine, public get_tile (ng, tile, my_lm, my_mm, itile, jtile, istr, iend, jstr, jend, istrm, istrr, istru, iendr, jstrm, jstrr, jstrv, jendr, istrb, iendb, istrp, iendp, istrt, iendt, jstrb, jendb, jstrp, jendp, jstrt, jendt, istrm3, istrm2, istrm1, istrum2, istrum1, iendp1, iendp2, iendp2i, iendp3, jstrm3, jstrm2, jstrm1, jstrvm2, jstrvm1, jendp1, jendp2, jendp2i, jendp3)
 
subroutine, public tile_bounds_1d (ng, tile, imax, istr, iend)
 
subroutine, public tile_bounds_2d (ng, tile, imax, jmax, itile, jtile, istr, iend, jstr, jend)
 
subroutine, private var_bounds (ng, tile, my_lm, my_mm, my_istr, my_iend, my_jstr, my_jend, istr, iend, jstr, jend, istrm, istrr, istru, iendr, jstrm, jstrr, jstrv, jendr, istrb, iendb, istrp, iendp, istrt, iendt, jstrb, jendb, jstrp, jendp, jstrt, jendt, istrm3, istrm2, istrm1, istrum2, istrum1, iendp1, iendp2, iendp2i, iendp3, jstrm3, jstrm2, jstrm1, jstrvm2, jstrvm1, jendp1, jendp2, jendp2i, jendp3)
 

Function/Subroutine Documentation

◆ get_bounds()

subroutine, public get_bounds_mod::get_bounds ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) gtype,
integer, intent(in) nghost,
integer, intent(in) my_im,
integer, intent(in) my_jm,
integer, intent(in) my_lm,
integer, intent(in) my_mm,
integer, intent(out) itile,
integer, intent(out) jtile,
integer, intent(out) lbi,
integer, intent(out) ubi,
integer, intent(out) lbj,
integer, intent(out) ubj )

Definition at line 32 of file get_bounds.F.

37!
38!=======================================================================
39! !
40! This routine compute grid bounds in the I- and J-directions. !
41! !
42! On Input: !
43! !
44! ng Nested grid number. !
45! tile Domain partition. !
46! gtype C-grid type. If zero, compute array allocation bounds.!
47! Otherwise, compute bounds for IO processing. !
48! Nghost Number of ghost-points in the halo region: !
49! Nghost = 0, compute non-overlapping bounds. !
50! Nghost > 0, compute overlapping bounds. !
51! my_Im Number of global grid points in the I-direction !
52! for each nested grid, [1:Ngrids]. !
53! my_Jm Number of global grid points in the J-direction !
54! for each nested grid, [1:Ngrids]. !
55! my_Lm Number of computational points in the I-direction !
56! for each nested grid, [1:Ngrids]. !
57! my_Mm Number of computational points in the J-direction. !
58! for each nested grid, [1:Ngrids]. !
59! !
60! On Output: !
61! !
62! Itile I-tile coordinate (a value from 0 to NtileI(ng)). !
63! Jtile J-tile coordinate (a value from 0 to NtileJ(ng)). !
64! LBi I-dimension Lower bound. !
65! UBi I-dimension Upper bound. !
66! LBj J-dimension Lower bound. !
67! UBj J-dimension Upper bound. !
68! !
69!=======================================================================
70!
71! Imported variable declarations.
72!
73 integer, intent(in) :: ng, tile, gtype, Nghost
74 integer, intent(in) :: my_Im, my_Jm
75 integer, intent(in) :: my_Lm, my_Mm
76 integer, intent(out) :: Itile, Jtile, LBi, UBi, LBj, UBj
77!
78! Local variable declarations.
79!
80 integer :: Imin, Imax, Jmin, Jmax
81#ifdef DISTRIBUTE
82 integer :: Iend, Istr, Jend, Jstr
83 integer :: IstrM, IstrR, IstrU, IendR
84 integer :: JstrM, JstrR, JstrV, JendR
85 integer :: IstrB, IendB, IstrP, IendP, IstrT, IendT
86 integer :: JstrB, JendB, JstrP, JendP, JstrT, JendT
87 integer :: Istrm3, Istrm2, Istrm1, IstrUm2, IstrUm1
88 integer :: Iendp1, Iendp2, Iendp2i, Iendp3
89 integer :: Jstrm3, Jstrm2, Jstrm1, JstrVm2, JstrVm1
90 integer :: Jendp1, Jendp2, Jendp2i, Jendp3
91 integer :: MyType
92!
93!-----------------------------------------------------------------------
94! Set array bounds in the I- and J-direction for distributed-memory
95! configurations.
96!-----------------------------------------------------------------------
97!
98! Set first and last grid-points according to staggered C-grid
99! classification. If gtype = 0, it returns the values needed for
100! array allocation. Otherwise, it returns the values needed for IO
101! processing.
102!
103 mytype=abs(gtype)
104 IF (mytype.eq.0) THEN
105# ifdef NESTING
106 IF (refinedgrid(ng).and.(refinescale(ng).gt.0) ) THEN
107 imin=-nghostpoints
108 imax=my_lm+nghostpoints
109 jmin=-nghostpoints
110 jmax=my_mm+nghostpoints
111 ELSE
112 IF (compositegrid(iwest,ng).or.ewperiodic(ng)) THEN
113 imin=-nghostpoints
114 ELSE
115 imin=0
116 END IF
117 IF (compositegrid(ieast,ng).or.ewperiodic(ng)) THEN
118 imax=my_lm+nghostpoints
119 ELSE
120 imax=my_im+1
121 END IF
122 IF (compositegrid(isouth,ng).or.nsperiodic(ng)) THEN
123 jmin=-nghostpoints
124 ELSE
125 jmin=0
126 END IF
127 IF (compositegrid(inorth,ng).or.nsperiodic(ng)) THEN
128 jmax=my_mm+nghostpoints
129 ELSE
130 jmax=my_jm+1
131 END IF
132 END IF
133# else
134 IF (ewperiodic(ng)) THEN
135 IF (nsperiodic(ng)) THEN
136 imin=-nghostpoints
137 imax=my_im+nghostpoints
138 jmin=-nghostpoints
139 jmax=my_jm+nghostpoints
140 ELSE
141 imin=-nghostpoints
142 imax=my_im+nghostpoints
143 jmin=0
144 jmax=my_jm+1
145 END IF
146 ELSE
147 IF (nsperiodic(ng)) THEN
148 imin=0
149 imax=my_im+1
150 jmin=-nghostpoints
151 jmax=my_jm+nghostpoints
152 ELSE
153 imin=0
154 imax=my_im+1
155 jmin=0
156 jmax=my_jm+1
157 END IF
158 END IF
159# endif
160 ELSE
161 IF ((mytype.eq.p2dvar).or.(mytype.eq.u2dvar).or. &
162 & (mytype.eq.p3dvar).or.(mytype.eq.u3dvar)) THEN
163 imin=1
164 ELSE
165 imin=0
166 END IF
167 imax=my_lm+1
168 IF ((mytype.eq.p2dvar).or.(mytype.eq.v2dvar).or. &
169 & (mytype.eq.p3dvar).or.(mytype.eq.v3dvar)) THEN
170 jmin=1
171 ELSE
172 jmin=0
173 END IF
174 jmax=my_mm+1
175 END IF
176!
177! Set physical, overlapping (Nghost>0) or non-overlapping (Nghost=0)
178! grid bounds according to tile rank.
179!
180 CALL get_tile (ng, tile, &
181 & my_lm, my_mm, &
182 & itile, jtile, &
183 & istr, iend, jstr, jend, &
184 & istrm, istrr, istru, iendr, &
185 & jstrm, jstrr, jstrv, jendr, &
186 & istrb, iendb, istrp, iendp, istrt, iendt, &
187 & jstrb, jendb, jstrp, jendp, jstrt, jendt, &
188 & istrm3, istrm2, istrm1, istrum2, istrum1, &
189 & iendp1, iendp2, iendp2i, iendp3, &
190 & jstrm3, jstrm2, jstrm1, jstrvm2, jstrvm1, &
191 & jendp1, jendp2, jendp2i, jendp3)
192!
193 IF ((itile.eq.-1).or.(itile.eq.0)) THEN
194 lbi=imin
195 ELSE
196 lbi=istr-nghost
197 END IF
198 IF ((itile.eq.-1).or.(itile.eq.(ntilei(ng)-1))) THEN
199 ubi=imax
200 ELSE
201 ubi=iend+nghost
202 END IF
203 IF ((jtile.eq.-1).or.(jtile.eq.0)) THEN
204 lbj=jmin
205 ELSE
206 lbj=jstr-nghost
207 END IF
208 IF ((jtile.eq.-1).or.(jtile.eq.(ntilej(ng)-1))) THEN
209 ubj=jmax
210 ELSE
211 ubj=jend+nghost
212 END IF
213#else
214!
215!-----------------------------------------------------------------------
216! Set array allocation bounds in the I- and J-direction for serial and
217! shared-memory configurations.
218!-----------------------------------------------------------------------
219!
220 IF (tile.eq.-1) THEN
221 itile=-1
222 jtile=-1
223 ELSE
224 jtile=tile/ntilei(ng)
225 itile=tile-jtile*ntilei(ng)
226 END IF
227!
228# ifdef NESTING
229 IF (refinedgrid(ng).and.(refinescale(ng).gt.0) ) THEN
230 lbi=-nghostpoints
231 ubi=my_lm+nghostpoints
232 lbj=-nghostpoints
233 ubj=my_mm+nghostpoints
234 ELSE
235 IF (compositegrid(iwest,ng).or.ewperiodic(ng)) THEN
236 lbi=-nghostpoints
237 ELSE
238 lbi=0
239 END IF
240 IF (compositegrid(ieast,ng).or.ewperiodic(ng)) THEN
241 ubi=my_lm+nghostpoints
242 ELSE
243 ubi=my_im+1
244 END IF
245 IF (compositegrid(isouth,ng).or.nsperiodic(ng)) THEN
246 lbj=-nghostpoints
247 ELSE
248 lbj=0
249 END IF
250 IF (compositegrid(inorth,ng).or.nsperiodic(ng)) THEN
251 ubj=my_mm+nghostpoints
252 ELSE
253 ubj=my_jm+1
254 END IF
255 END IF
256# else
257
258 IF (ewperiodic(ng)) THEN
259 IF (nsperiodic(ng)) THEN
260 lbi=-nghostpoints
261 ubi=my_im+nghostpoints
262 lbj=-nghostpoints
263 ubj=my_jm+nghostpoints
264 ELSE
265 lbi=-nghostpoints
266 ubi=my_im+nghostpoints
267 lbj=0
268 ubj=my_jm+1
269 END IF
270 ELSE
271 IF (nsperiodic(ng)) THEN
272 lbi=0
273 ubi=my_im+1
274 lbj=-nghostpoints
275 ubj=my_jm+nghostpoints
276 ELSE
277 lbi=0
278 ubi=my_im+1
279 lbj=0
280 ubj=my_jm+1
281 END IF
282 END IF
283# endif
284#endif
285!
286 RETURN

References mod_scalars::compositegrid, mod_scalars::ewperiodic, get_tile(), mod_scalars::ieast, mod_scalars::inorth, mod_scalars::isouth, mod_scalars::iwest, mod_param::nghostpoints, mod_scalars::nsperiodic, mod_param::ntilei, mod_param::ntilej, mod_param::p2dvar, mod_param::p3dvar, mod_scalars::refinedgrid, mod_scalars::refinescale, mod_param::u2dvar, mod_param::u3dvar, mod_param::v2dvar, and mod_param::v3dvar.

Referenced by get_domain(), tile_indices_mod::tile_indices(), and tile_indices_mod::tile_obs_bounds().

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

◆ get_domain()

subroutine, public get_bounds_mod::get_domain ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) gtype,
integer, intent(in) nghost,
integer, intent(in) my_im,
integer, intent(in) my_jm,
integer, intent(in) my_lm,
integer, intent(in) my_mm,
real(r8), intent(in) epsilon,
logical, intent(in) lfullgrid,
real(r8), intent(out) xmin,
real(r8), intent(out) xmax,
real(r8), intent(out) ymin,
real(r8), intent(out) ymax )

Definition at line 289 of file get_bounds.F.

294!
295!=======================================================================
296! !
297! This routine computes tile minimum and maximum fractional grid !
298! coordinates. !
299! !
300! On Input: !
301! !
302! ng Nested grid number. !
303! tile Domain partition. !
304! Nghost Number of ghost-points in the halo region: !
305! Nghost = 0, compute non-overlapping coordinates. !
306! Nghost > 0, compute overlapping bounds. !
307! gtype C-grid type !
308! my_Im Number of global grid points in the I-direction. !
309! my_Jm Number of global grid points in the J-direction. !
310! my_Lm Number of computational points in the I-direction. !
311! my_Mm Number of computational points in the J-direction. !
312! epsilon Small value to add to Xmax and Ymax when the tile !
313! is lying on the eastern and northern boundaries !
314! of the grid. This is usefull when processing !
315! observations. !
316! Lfullgrid Switch to include interior and boundaries points !
317! (TRUE) or just interior points (FALSE). !
318! !
319! On Output: !
320! !
321! Xmin Minimum tile fractional X-coordinate. !
322! Xmax Maximum tile fractional X-coordinate. !
323! Ymin Minimum tile fractional Y-coordinate. !
324! Ymax Maximum tile fractional Y-coordinate. !
325! !
326!=======================================================================
327!
328! Imported variable declarations.
329!
330 logical, intent(in) :: Lfullgrid
331!
332 integer, intent(in) :: ng, tile, gtype, Nghost
333 integer, intent(in) :: my_Im, my_Jm, my_Lm, my_Mm
334!
335 real(r8), intent(in) :: epsilon
336 real(r8), intent(out) :: Xmin, Xmax, Ymin, Ymax
337!
338! Local variable declarations.
339!
340 integer :: Imin, Imax, Jmin, Jmax
341 integer :: Itile, Jtile
342!
343!-----------------------------------------------------------------------
344! Computes tile minimum and maximum fractional-grid coordinates.
345!-----------------------------------------------------------------------
346!
347 CALL get_bounds (ng, tile, gtype, nghost, &
348 & my_im, my_jm, &
349 & my_lm, my_mm, &
350 & itile, jtile, &
351 & imin, imax, jmin, jmax)
352!
353! Include interior and boundary points.
354!
355 IF (lfullgrid) THEN
356 IF ((itile.eq.0).and. &
357 & ((gtype.eq.r2dvar).or.(gtype.eq.r3dvar).or. &
358 & (gtype.eq.v2dvar).or.(gtype.eq.v3dvar))) THEN
359 xmin=real(imin,r8)
360 ELSE
361 xmin=real(imin,r8)-0.5_r8
362 END IF
363 IF (itile.eq.(ntilei(ng)-1)) THEN
364 IF ((gtype.eq.u2dvar).or.(gtype.eq.u3dvar)) THEN
365 xmax=real(imax,r8)-0.5_r8
366 ELSE
367 xmax=real(imax,r8)
368 END IF
369 ELSE
370 xmax=real(imax,r8)+0.5_r8
371 END IF
372 IF ((jtile.eq.0).and. &
373 & ((gtype.eq.r2dvar).or.(gtype.eq.r3dvar).or. &
374 & (gtype.eq.u2dvar).or.(gtype.eq.u3dvar))) THEN
375 ymin=real(jmin,r8)
376 ELSE
377 ymin=real(jmin,r8)-0.5_r8
378 END IF
379 IF (jtile.eq.(ntilej(ng)-1)) THEN
380 IF ((gtype.eq.v2dvar).or.(gtype.eq.v3dvar)) THEN
381 ymax=real(jmax,r8)-0.5_r8
382 ELSE
383 ymax=real(jmax,r8)
384 END IF
385 ELSE
386 ymax=real(jmax,r8)+0.5_r8
387 END IF
388!
389! Include only interior points.
390!
391 ELSE
392 IF (itile.eq.0) THEN
393 IF ((gtype.eq.u2dvar).or.(gtype.eq.u3dvar)) THEN
394 xmin=real(imin,r8)
395 ELSE
396 xmin=real(imin,r8)+0.5_r8
397 END IF
398 ELSE
399 xmin=real(imin,r8)-0.5_r8
400 END IF
401 IF (itile.eq.(ntilei(ng)-1)) THEN
402 IF ((gtype.eq.u2dvar).or.(gtype.eq.u3dvar)) THEN
403 xmax=real(imax,r8)-1.0_r8
404 ELSE
405 xmax=real(imax,r8)-0.5_r8
406 END IF
407 ELSE
408 xmax=real(imax,r8)+0.5_r8
409 END IF
410 IF (jtile.eq.0) THEN
411 IF ((gtype.eq.v2dvar).or.(gtype.eq.v3dvar)) THEN
412 ymin=real(jmin,r8)
413 ELSE
414 ymin=real(jmin,r8)+0.5
415 END IF
416 ELSE
417 ymin=real(jmin,r8)-0.5_r8
418 END IF
419 IF (jtile.eq.(ntilej(ng)-1)) THEN
420 IF ((gtype.eq.v2dvar).or.(gtype.eq.v3dvar)) THEN
421 ymax=real(jmax,r8)-1.0_r8
422 ELSE
423 ymax=real(jmax,r8)-0.5_r8
424 END IF
425 ELSE
426 ymax=real(jmax,r8)+0.5_r8
427 END IF
428 END IF
429!
430! If tile lie at the grid eastern or northen boundary, add provided
431! offset value to allow processing at those boundaries.
432!
433 IF (itile.eq.(ntilei(ng)-1)) THEN
434 xmax=xmax+epsilon
435 END IF
436 IF (jtile.eq.(ntilej(ng)-1)) THEN
437 ymax=ymax+epsilon
438 END IF
439!
440 RETURN

References get_bounds(), mod_param::ntilei, mod_param::ntilej, mod_param::r2dvar, mod_param::r3dvar, mod_param::u2dvar, mod_param::u3dvar, mod_param::v2dvar, and mod_param::v3dvar.

Referenced by tile_indices_mod::tile_obs_bounds().

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

◆ get_domain_edges()

subroutine, public get_bounds_mod::get_domain_edges ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) my_lm,
integer, intent(in) my_mm,
logical, intent(out) eastern_edge,
logical, intent(out) western_edge,
logical, intent(out) northern_edge,
logical, intent(out) southern_edge,
logical, intent(out) northeast_corner,
logical, intent(out) northwest_corner,
logical, intent(out) southeast_corner,
logical, intent(out) southwest_corner,
logical, intent(out) northeast_test,
logical, intent(out) northwest_test,
logical, intent(out) southeast_test,
logical, intent(out) southwest_test )

Definition at line 443 of file get_bounds.F.

457!
458!=======================================================================
459! !
460! This routine sets the logical switches (T/F) needed for processing !
461! model variables in tiles adjacent to the domain boundary edges. It !
462! facilitates complicated nesting configurations. !
463! !
464! On Input: !
465! !
466! ng Nested grid number. !
467! tile Domain partition. !
468! my_Lm Number of computational points in the I-direction. !
469! my_Mm Number of computational points in the J-direction. !
470! !
471! On Output: !
472! !
473! Eastern_Edge tile next to the domain eastern boundary !
474! Western_Edge tile next to the domain western boundary !
475! Northern_Edge tile next to the domain northern boundary !
476! Southern_Edge tile next to the domain southern boundary !
477! !
478! NorthEast_Corner tile next to the domain northeastern corner !
479! NorthWest_Corner tile next to the domain northwestern corner !
480! SouthEast_Corner tile next to the domain southeastern corner !
481! SouthWest_Corner tile next to the domain southwestern corner !
482! !
483! NorthEast_Test test for tiles in the northeastern corner !
484! NorthWest_Test test for tiles in the northwestern corner !
485! SouthEast_Test test for tiles in the southeastern corner !
486! SouthWest_Test test for tiles in the southwestern corner !
487! !
488!=======================================================================
489!
490! Imported variable declarations.
491!
492 integer, intent(in) :: ng, tile
493 integer, intent(in) :: my_Lm, My_Mm
494!
495 logical, intent(out) :: Eastern_Edge
496 logical, intent(out) :: Western_Edge
497 logical, intent(out) :: Northern_Edge
498 logical, intent(out) :: Southern_Edge
499 logical, intent(out) :: NorthEast_Corner
500 logical, intent(out) :: NorthWest_Corner
501 logical, intent(out) :: SouthEast_Corner
502 logical, intent(out) :: SouthWest_Corner
503 logical, intent(out) :: NorthEast_Test
504 logical, intent(out) :: NorthWest_Test
505 logical, intent(out) :: SouthEast_Test
506 logical, intent(out) :: SouthWest_Test
507!
508! Local variable declarations.
509!
510 integer :: Istr, Iend, Jstr, Jend
511 integer :: Itile, Jtile
512!
513!-----------------------------------------------------------------------
514! Compute Itile and Jtile.
515!-----------------------------------------------------------------------
516!
517 IF (tile.eq.-1) THEN
518 itile=-1
519 jtile=-1
520 istr=1
521 iend=my_lm
522 jstr=1
523 jend=my_mm
524 ELSE
525 CALL tile_bounds_2d (ng, tile, &
526 & my_lm, my_mm, &
527 & itile, jtile, &
528 & istr, iend, jstr, jend)
529 END IF
530!
531!-----------------------------------------------------------------------
532! Set the logical switches (T/F) needed for processing model variables
533! in tiles adjacent to the domain boundary edges.
534!-----------------------------------------------------------------------
535!
536! HGA: Need to add the logic for composed grids.
537
538 IF (tile.eq.-1) THEN
539!
540! Set switches for the full grid (tile=-1) to TRUE, since it contains
541! all the boundary edges and corners. This is a special case use for
542! other purposes and need only in routine "var_bounds".
543!
544 western_edge=.true.
545 eastern_edge=.true.
546 southern_edge=.true.
547 northern_edge=.true.
548!
549 southwest_test=.true.
550 southeast_test=.true.
551 northwest_test=.true.
552 northeast_test=.true.
553!
554 southwest_corner=.true.
555 southeast_corner=.true.
556 northwest_corner=.true.
557 northeast_corner=.true.
558 ELSE
559!
560! Is the tile adjacent to the western or eastern domain
561! boundary edge?
562!
563 IF (itile.eq.0) THEN
564 western_edge=.true. ! (Istr.eq.1)
565 ELSE
566 western_edge=.false. ! elsewhere
567 END IF
568!
569 IF (itile.eq.(ntilei(ng)-1)) THEN
570 eastern_edge=.true. ! (Iend.eq.my_Lm))
571 ELSE
572 eastern_edge=.false. ! elsewhere
573 END IF
574!
575! Is the tile adjacent to the southern or northern domain
576! boundary edge?
577!
578 IF (jtile.eq.0) THEN
579 southern_edge=.true. ! (Jstr.eq.1)
580 ELSE
581 southern_edge=.false. ! elsewhere
582 END IF
583!
584 IF (jtile.eq.(ntilej(ng)-1)) THEN
585 northern_edge=.true. ! (Jend.eq.my_Mm)
586 ELSE
587 northern_edge=.false. ! elsewhere
588 END IF
589!
590! Is the tile adjacent to the southwestern domain corner?
591!
592 IF ((itile.eq.0).and. &
593 & (jtile.eq.0)) THEN
594 southwest_corner=.true. ! (Istr.eq.1).and.
595 southwest_test =.true. ! (Jstr.eq.1)
596 ELSE
597 southwest_corner=.false. ! elsewhere
598#ifdef DISTRIBUTE
599 southwest_test =.true.
600#else
601 southwest_test =.false.
602#endif
603 END IF
604!
605! Is the tile adjacent to the southeastern domain corner?
606!
607 IF ((itile.eq.(ntilei(ng)-1)).and. &
608 & (jtile.eq.0)) THEN
609 southeast_corner=.true. ! (Iend.eq.my_Lm).and.
610 southeast_test =.true. ! (Jstr.eq.1)
611 ELSE
612 southeast_corner=.false. ! elsewhere
613#ifdef DISTRIBUTE
614 southeast_test =.true.
615#else
616 southeast_test =.false.
617#endif
618 END IF
619!
620! Is the tile adjacent to the northwestern domain corner?
621!
622 IF ((itile.eq.0).and. &
623 & (jtile.eq.(ntilej(ng)-1))) THEN
624 northwest_corner=.true. ! (Istr.eq.1).and.
625 northwest_test =.true. ! (Jend.eq.my_Mm)
626 ELSE
627 northwest_corner=.false. ! elsewhere
628#ifdef DISTRIBUTE
629 northwest_test =.true.
630#else
631 northwest_test =.false.
632#endif
633 END IF
634!
635! Is the tile adjacent to the northeastern domain corner?
636!
637 IF ((itile.eq.(ntilei(ng)-1)).and. &
638 & (jtile.eq.(ntilej(ng)-1))) THEN
639 northeast_corner=.true. ! (Iend.eq.my_Lm).and.
640 northeast_test =.true. ! (Jend.eq.my_Mm)
641 ELSE
642 northeast_corner=.false. ! elsewhere
643#ifdef DISTRIBUTE
644 northeast_test =.true.
645#else
646 northeast_test =.false.
647#endif
648 END IF
649 END IF
650!
651 RETURN

References mod_param::ntilei, mod_param::ntilej, and tile_bounds_2d().

Referenced by tile_indices_mod::tile_indices().

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

◆ get_iobounds()

subroutine, public get_bounds_mod::get_iobounds ( integer, intent(in) ng,
integer, intent(in) my_lm,
integer, intent(in) my_mm,
type (t_bounds), dimension(1:ngrids), intent(in) my_bounds,
type (t_iobounds), dimension(1:ngrids), intent(out) my_iobounds )

Definition at line 654 of file get_bounds.F.

656!
657!=======================================================================
658! !
659! This routine computes the horizontal lower bound, upper bound, and !
660! grid size for IO (NetCDF) variables. Nested grids require special !
661! attention due to their connetivity. !
662! !
663! On Input: !
664! !
665! ng Nested grid number. !
666! my_Lm Number of computational points in the I-direction. !
667! my_Mm Number of computational points in the J-direction. !
668! my_BOUNDS Lower and upper bounds indices structure per domain !
669! partition for all grids. !
670! !
671! On Output: !
672! !
673! my_IOBOUNDS I/O lower and upper bounds indices structure per !
674! domain partition for all grids. !
675! !
676! ILB_psi I-direction lower bound (PSI) !
677! IUB_psi I-direction upper bound (PSI) !
678! JLB_psi J-direction lower bound (PSI) !
679! JUB_psi J-direction upper bound (PSI) !
680! !
681! ILB_rho I-direction lower bound (RHO) !
682! IUB_rho I-direction upper bound (RHO) !
683! JLB_rho J-direction lower bound (RHO) !
684! JUB_rho J-direction upper bound (RHO) !
685! !
686! ILB_u I-direction lower bound (U) !
687! IUB_u I-direction upper bound (U) !
688! JLB_u J-direction lower bound (U) !
689! JUB_u J-direction upper bound (U) !
690! !
691! ILB_v I-direction lower bound (V) !
692! IUB_v I-direction upper bound (V) !
693! JLB_v J-direction lower bound (V) !
694! JUB_v J-direction upper bound (V) !
695! !
696! xi_psi Number of I-direction points (PSI) !
697! xi_rho Number of I-direction points (RHO) !
698! xi_u Number of I-direction points (U) !
699! xi_v Number of I-direction points (V) !
700! !
701! eta_psi Number of J-direction points (PSI) !
702! eta_rho Number of J-direction points (RHO) !
703! eta_u Number of I-direction points (U) !
704! eta_v Number of I-direction points (V) !
705! !
706!=======================================================================
707!
708! Imported variable declarations.
709!
710 integer, intent(in) :: ng
711 integer, intent(in) :: my_Lm, my_Mm
712!
713 TYPE (T_BOUNDS), intent(in) :: my_BOUNDS(1:Ngrids)
714 TYPE (T_IOBOUNDS), intent(out) :: my_IOBOUNDS(1:Ngrids)
715!
716!-----------------------------------------------------------------------
717! Set IO lower/upper bounds and grid size for each C-grid type
718! variable.
719!-----------------------------------------------------------------------
720!
721! Recall that in non-nested applications the horizontal range,
722! including interior and boundary points, for all variable types
723! are:
724!
725! PSI-type [xi_psi, eta_psi] = [1:my_Lm+1, 1:my_Mm+1]
726! RHO-type [xi_rho, eta_rho] = [0:my_Lm+1, 0:my_Mm+1]
727! U-type [xi_u, eta_u ] = [1:my_Lm+1, 0:my_Mm+1]
728! V-type [xi_v, eta_v ] = [0:my_Lm+1, 1:my_Mm+1]
729!
730 my_iobounds(ng) % ILB_psi = 1
731 my_iobounds(ng) % IUB_psi = my_lm+1
732 my_iobounds(ng) % JLB_psi = 1
733 my_iobounds(ng) % JUB_psi = my_mm+1
734!
735 my_iobounds(ng) % ILB_rho = 0
736 my_iobounds(ng) % IUB_rho = my_lm+1
737 my_iobounds(ng) % JLB_rho = 0
738 my_iobounds(ng) % JUB_rho = my_mm+1
739!
740 my_iobounds(ng) % ILB_u = 1
741 my_iobounds(ng) % IUB_u = my_lm+1
742 my_iobounds(ng) % JLB_u = 0
743 my_iobounds(ng) % JUB_u = my_mm+1
744!
745 my_iobounds(ng) % ILB_v = 0
746 my_iobounds(ng) % IUB_v = my_lm+1
747 my_iobounds(ng) % JLB_v = 1
748 my_iobounds(ng) % JUB_v = my_mm+1
749!
750! Set IO NetCDF files horizontal dimension size. Recall that NetCDF
751! does not support arrays with zero index as an array element.
752!
753 my_iobounds(ng) % IorJ = my_bounds(ng) % UBij - &
754 & my_bounds(ng) % LBij + 1
755!
756 my_iobounds(ng) % xi_psi = my_iobounds(ng) % IUB_psi - &
757 & my_iobounds(ng) % ILB_psi + 1
758 my_iobounds(ng) % xi_rho = my_iobounds(ng) % IUB_rho - &
759 & my_iobounds(ng) % ILB_rho + 1
760 my_iobounds(ng) % xi_u = my_iobounds(ng) % IUB_u - &
761 & my_iobounds(ng) % ILB_u + 1
762 my_iobounds(ng) % xi_v = my_iobounds(ng) % IUB_v - &
763 & my_iobounds(ng) % ILB_v + 1
764!
765 my_iobounds(ng) % eta_psi = my_iobounds(ng) % JUB_psi - &
766 & my_iobounds(ng) % JLB_psi + 1
767 my_iobounds(ng) % eta_rho = my_iobounds(ng) % JUB_rho - &
768 & my_iobounds(ng) % JLB_rho + 1
769 my_iobounds(ng) % eta_u = my_iobounds(ng) % JUB_u - &
770 & my_iobounds(ng) % JLB_u + 1
771 my_iobounds(ng) % eta_v = my_iobounds(ng) % JUB_v - &
772 & my_iobounds(ng) % JLB_v + 1
773!
774 RETURN

Referenced by tile_indices_mod::tile_indices().

Here is the caller graph for this function:

◆ get_tile()

subroutine, public get_bounds_mod::get_tile ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) my_lm,
integer, intent(in) my_mm,
integer, intent(out) itile,
integer, intent(out) jtile,
integer, intent(out) istr,
integer, intent(out) iend,
integer, intent(out) jstr,
integer, intent(out) jend,
integer, intent(out) istrm,
integer, intent(out) istrr,
integer, intent(out) istru,
integer, intent(out) iendr,
integer, intent(out) jstrm,
integer, intent(out) jstrr,
integer, intent(out) jstrv,
integer, intent(out) jendr,
integer, intent(out) istrb,
integer, intent(out) iendb,
integer, intent(out) istrp,
integer, intent(out) iendp,
integer, intent(out) istrt,
integer, intent(out) iendt,
integer, intent(out) jstrb,
integer, intent(out) jendb,
integer, intent(out) jstrp,
integer, intent(out) jendp,
integer, intent(out) jstrt,
integer, intent(out) jendt,
integer, intent(out) istrm3,
integer, intent(out) istrm2,
integer, intent(out) istrm1,
integer, intent(out) istrum2,
integer, intent(out) istrum1,
integer, intent(out) iendp1,
integer, intent(out) iendp2,
integer, intent(out) iendp2i,
integer, intent(out) iendp3,
integer, intent(out) jstrm3,
integer, intent(out) jstrm2,
integer, intent(out) jstrm1,
integer, intent(out) jstrvm2,
integer, intent(out) jstrvm1,
integer, intent(out) jendp1,
integer, intent(out) jendp2,
integer, intent(out) jendp2i,
integer, intent(out) jendp3 )

Definition at line 777 of file get_bounds.F.

789!
790!=======================================================================
791! !
792! This routine computes the starting and ending horizontal indices !
793! for each sub-domain partition or tile. !
794! !
795! On Input: !
796! !
797! ng Nested grid number (integer). !
798! tile Sub-domain partition. !
799! my_Lm Number of computational points in the I-direction. !
800! my_Mm Number of computational points in the J-direction. !
801! !
802! On Output: !
803! !
804! Itile I-tile coordinate (a value from 0 to NtileI(ng)). !
805! Jtile J-tile coordinate (a value from 0 to NtileJ(ng)). !
806! !
807! Istr Starting tile index in the I-direction. !
808! Iend Ending tile index in the I-direction. !
809! Jstr Starting tile index in the J-direction. !
810! Jend Ending tile index in the J-direction. !
811! !
812! IstrR Starting tile index in the I-direction (RHO-points) !
813! IstrU Starting tile index in the I-direction (U-points) !
814! IendR Ending tile index in the I-direction (RHO_points) !
815! !
816! JstrR Starting tile index in the J-direction (RHO-points) !
817! JstrV Starting tile index in the J-direction (V-points) !
818! JendR Ending tile index in the J-direction (RHO_points) !
819! !
820! IstrB Starting nest tile in the I-direction (RHO-, V-obc) !
821! IendB Ending nest tile in the I-direction (RHO-, V-obc) !
822! IstrM Starting nest tile in the I-direction (PSI-, U-obc) !
823! IstrP Starting nest tile in the I-direction (PSI, U-points) !
824! IendP Ending nest tile in the I-direction (PSI) !
825! IstrT Starting nest tile in the I-direction (RHO-points) !
826! IendT Ending nest tile in the I-direction (RHO_points) !
827! !
828! JstrB Starting nest tile in the J-direction (RHO-, U-obc) !
829! JendB Ending nest tile in the J-direction (RHO-, U-obc) !
830! JstrM Starting nest tile in the J-direction (PSI-, V-obc) !
831! JstrP Starting nest tile in the J-direction (PSI, V-points) !
832! JendP Ending nest tile in the J-direction (PSI) !
833! JstrT Starting nest tile in the J-direction (RHO-points) !
834! JendT Ending nest tile in the J-direction (RHO-points) !
835! !
836! Istrm3 Starting private I-halo computations, Istr-3 !
837! Istrm2 Starting private I-halo computations, Istr-2 !
838! Istrm1 Starting private I-halo computations, Istr-1 !
839! IstrUm2 Starting private I-halo computations, IstrU-2 !
840! IstrUm1 Starting private I-halo computations, IstrU-1 !
841! Iendp1 Ending private I-halo computations, Iend+1 !
842! Iendp2 Ending private I-halo computations, Iend+2 !
843! Iendp2i Ending private I-halo computations, Iend+2 interior !
844! Iendp3 Ending private I-halo computations, Iend+3 !
845! !
846! Jstrm3 Starting private J-halo computations, Jstr-3 !
847! Jstrm2 Starting private J-halo computations, Jstr-2 !
848! Jstrm1 Starting private J-halo computations, Jstr-1 !
849! JstrVm2 Starting private J-halo computations, JstrV-2 !
850! JstrVm1 Starting private J-halo computations, JstrV-1 !
851! Jendp1 Ending private J-halo computations, Jend+1 !
852! Jendp2 Ending private J-halo computations, Jend+2 !
853! Jendp2i Ending private J-halo computations, Jend+2 interior !
854! Jendp3 Ending private J-halo computations, Jend+3 !
855! !
856!======================================================================!
857!
858! Imported variable declarations.
859!
860 integer, intent(in) :: ng, tile
861 integer, intent(in) :: my_Lm, my_Mm
862 integer, intent(out) :: Itile, Jtile
863 integer, intent(out) :: Iend, Istr, Jend, Jstr
864 integer, intent(out) :: IstrM, IstrR, IstrU, IendR
865 integer, intent(out) :: JstrM, JstrR, JstrV, JendR
866 integer, intent(out) :: IstrB, IendB, IstrP, IendP, IstrT, IendT
867 integer, intent(out) :: JstrB, JendB, JstrP, JendP, JstrT, JendT
868 integer, intent(out) :: Istrm3, Istrm2, Istrm1, IstrUm2, IstrUm1
869 integer, intent(out) :: Iendp1, Iendp2, Iendp2i, Iendp3
870 integer, intent(out) :: Jstrm3, Jstrm2, Jstrm1, JstrVm2, JstrVm1
871 integer, intent(out) :: Jendp1, Jendp2, Jendp2i, Jendp3
872!
873! Local variable declarations.
874!
875 integer :: my_Istr, my_Iend, my_Jstr, my_Jend
876!
877!-----------------------------------------------------------------------
878! Set physical non-overlapping grid bounds according to tile rank.
879!-----------------------------------------------------------------------
880!
881! Non-tiled grid bounds. This is used in serial or shared-memory
882! modes to compute values in the full grid outside of parallel
883! regions.
884!
885 IF (tile.eq.-1) THEN
886 itile=-1
887 jtile=-1
888 my_istr=1
889 my_iend=my_lm
890 my_jstr=1
891 my_jend=my_mm
892!
893! Tiled grids bounds.
894!
895 ELSE
896 CALL tile_bounds_2d (ng, tile, &
897 & my_lm, my_mm, &
898 & itile, jtile, &
899 & my_istr, my_iend, my_jstr, my_jend)
900 END IF
901!
902! Compute C-staggered variables bounds from tile bounds.
903!
904 CALL var_bounds (ng, tile, &
905 & my_lm, my_mm, &
906 & my_istr, my_iend, my_jstr, my_jend, &
907 & istr, iend, jstr, jend, &
908 & istrm, istrr, istru, iendr, &
909 & jstrm, jstrr, jstrv, jendr, &
910 & istrb, iendb, istrp, iendp, istrt, iendt, &
911 & jstrb, jendb, jstrp, jendp, jstrt, jendt, &
912 & istrm3, istrm2, istrm1, istrum2, istrum1, &
913 & iendp1, iendp2, iendp2i, iendp3, &
914 & jstrm3, jstrm2, jstrm1, jstrvm2, jstrvm1, &
915 & jendp1, jendp2, jendp2i, jendp3)
916!
917 RETURN

References tile_bounds_2d(), and var_bounds().

Referenced by get_bounds(), and tile_indices_mod::tile_indices().

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

◆ tile_bounds_1d()

subroutine, public get_bounds_mod::tile_bounds_1d ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) imax,
integer, intent(out) istr,
integer, intent(out) iend )

Definition at line 920 of file get_bounds.F.

921!
922!=======================================================================
923! !
924! This routine computes the starting and ending indices for the 1D !
925! decomposition between all available threads or partitions. !
926! !
927! 1 _____________________ Imax !
928! !
929! On Input: !
930! !
931! ng Nested grid number !
932! tile Thread or partition !
933! Imax Global number of points !
934! !
935! On Output: !
936! !
937! Istr Starting partition index !
938! Iend Ending partition index !
939! !
940!======================================================================!
941!
942! Imported variable declarations.
943!
944 integer, intent(in) :: ng, tile, Imax
945 integer, intent(out) :: Iend, Istr
946!
947! Local variable declarations.
948!
949 integer :: ChunkSize, Margin, Nnodes
950!
951!-----------------------------------------------------------------------
952! Compute 1D decomposition starting and ending indices.
953!-----------------------------------------------------------------------
954!
955 nnodes=ntilei(ng)*ntilej(ng)
956 chunksize=(imax+nnodes-1)/nnodes
957 margin=(nnodes*chunksize-imax)/2
958!
959 IF (imax.ge.nnodes) THEN
960 istr=1+tile*chunksize-margin
961 iend=istr+chunksize-1
962 istr=max(istr,1)
963 iend=min(iend,imax)
964 ELSE
965 istr=1
966 iend=imax
967 END IF
968!
969 RETURN

References mod_param::ntilei, and mod_param::ntilej.

Referenced by nf_fread2d_mod::nf_fread2d::nf90_fread2d(), nf_fread3d_mod::nf_fread3d::nf90_fread3d(), nf_fread4d_mod::nf_fread4d::nf90_fread4d(), nf_fwrite2d_mod::nf_fwrite2d::nf90_fwrite2d(), nf_fwrite3d_mod::nf_fwrite3d::nf90_fwrite3d(), and nf_fwrite4d_mod::nf_fwrite4d::nf90_fwrite4d().

Here is the caller graph for this function:

◆ tile_bounds_2d()

subroutine, public get_bounds_mod::tile_bounds_2d ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) imax,
integer, intent(in) jmax,
integer, intent(out) itile,
integer, intent(out) jtile,
integer, intent(out) istr,
integer, intent(out) iend,
integer, intent(out) jstr,
integer, intent(out) jend )

Definition at line 972 of file get_bounds.F.

974!
975!=======================================================================
976! !
977! This routine computes the starting and ending horizontal indices !
978! for each sub-domain partition or tile for a grid bounded between !
979! (1,1) and (Imax,Jmax): !
980! !
981! _________ (Imax,Jmax) !
982! | | !
983! | | !
984! |_________| !
985! (1,1) !
986! !
987! On Input: !
988! !
989! ng Nested grid number !
990! tile Sub-domain partition !
991! Imax Global number of points in the I-direction !
992! Jmax Global number of points in the J-direction !
993! !
994! On Output: !
995! !
996! Itile I-tile coordinate (a value from 0 to NtileI(ng)) !
997! Jtile J-tile coordinate (a value from 0 to NtileJ(ng)) !
998! Istr Starting tile index in the I-direction !
999! Iend Ending tile index in the I-direction !
1000! Jstr Starting tile index in the J-direction !
1001! Jend Ending tile index in the J-direction !
1002! !
1003!======================================================================!
1004!
1005! Imported variable declarations.
1006!
1007 integer, intent(in) :: ng, tile, Imax, Jmax
1008 integer, intent(out) :: Itile, Jtile
1009 integer, intent(out) :: Iend, Istr, Jend, Jstr
1010!
1011! Local variable declarations.
1012!
1013 integer :: ChunkSizeI, ChunkSizeJ, MarginI, MarginJ
1014!
1015!-----------------------------------------------------------------------
1016! Compute tile decomposition for a horizontal grid bounded between
1017! (1,1) and (Imax,Jmax).
1018!-----------------------------------------------------------------------
1019!
1020 chunksizei=(imax+ntilei(ng)-1)/ntilei(ng)
1021 chunksizej=(jmax+ntilej(ng)-1)/ntilej(ng)
1022 margini=(ntilei(ng)*chunksizei-imax)/2
1023 marginj=(ntilej(ng)*chunksizej-jmax)/2
1024 jtile=tile/ntilei(ng)
1025 itile=tile-jtile*ntilei(ng)
1026!
1027! Tile bounds in the I-direction.
1028!
1029 istr=1+itile*chunksizei-margini
1030 iend=istr+chunksizei-1
1031 istr=max(istr,1)
1032 iend=min(iend,imax)
1033!
1034! Tile bounds in the J-direction.
1035!
1036 jstr=1+jtile*chunksizej-marginj
1037 jend=jstr+chunksizej-1
1038 jstr=max(jstr,1)
1039 jend=min(jend,jmax)
1040!
1041 RETURN

References mod_param::ntilei, and mod_param::ntilej.

Referenced by get_domain_edges(), get_tile(), and nf_fread2d_mod::nf_fread2d::nf90_fread2d().

Here is the caller graph for this function:

◆ var_bounds()

subroutine, private get_bounds_mod::var_bounds ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) my_lm,
integer, intent(in) my_mm,
integer, intent(in) my_istr,
integer, intent(in) my_iend,
integer, intent(in) my_jstr,
integer, intent(in) my_jend,
integer, intent(out) istr,
integer, intent(out) iend,
integer, intent(out) jstr,
integer, intent(out) jend,
integer, intent(out) istrm,
integer, intent(out) istrr,
integer, intent(out) istru,
integer, intent(out) iendr,
integer, intent(out) jstrm,
integer, intent(out) jstrr,
integer, intent(out) jstrv,
integer, intent(out) jendr,
integer, intent(out) istrb,
integer, intent(out) iendb,
integer, intent(out) istrp,
integer, intent(out) iendp,
integer, intent(out) istrt,
integer, intent(out) iendt,
integer, intent(out) jstrb,
integer, intent(out) jendb,
integer, intent(out) jstrp,
integer, intent(out) jendp,
integer, intent(out) jstrt,
integer, intent(out) jendt,
integer, intent(out) istrm3,
integer, intent(out) istrm2,
integer, intent(out) istrm1,
integer, intent(out) istrum2,
integer, intent(out) istrum1,
integer, intent(out) iendp1,
integer, intent(out) iendp2,
integer, intent(out) iendp2i,
integer, intent(out) iendp3,
integer, intent(out) jstrm3,
integer, intent(out) jstrm2,
integer, intent(out) jstrm1,
integer, intent(out) jstrvm2,
integer, intent(out) jstrvm1,
integer, intent(out) jendp1,
integer, intent(out) jendp2,
integer, intent(out) jendp2i,
integer, intent(out) jendp3 )
private

Definition at line 1044 of file get_bounds.F.

1056!
1057!=======================================================================
1058! !
1059! This routine computes the computational grid starting and ending !
1060! horizontal indices for each C-staggered variable in terms of the !
1061! physical grid sub-domain partition or tile. !
1062! !
1063! On Input: !
1064! !
1065! ng Nested grid number !
1066! tile Domain partition !
1067! my_Lm Number of computational points in the I-direction. !
1068! my_Mm Number of computational points in the J-direction. !
1069! my_Istr Physical grid starting tile index in the I-direction !
1070! my_Iend Physical grid ending tile index in the I-direction !
1071! my_Jstr Physical grid starting tile index in the J-direction !
1072! my_Jend Physical grid ending tile index in the J-direction !
1073! !
1074! On Output: !
1075! !
1076! Istr Starting tile index in the I-direction !
1077! Iend Ending tile index in the I-direction !
1078! Jstr Starting tile index in the J-direction !
1079! Jend Ending tile index in the J-direction !
1080! !
1081! IstrR Starting tile index in the I-direction (RHO-points) !
1082! IstrU Starting tile index in the I-direction (U-points) !
1083! IendR Ending tile index in the I-direction (RHO_points) !
1084! !
1085! JstrR Starting tile index in the J-direction (RHO-points) !
1086! JstrV Starting tile index in the J-direction (V-points) !
1087! JendR Ending tile index in the J-direction (RHO_points) !
1088! !
1089! IstrB Starting nest tile in the I-direction (RHO-, V-obc) !
1090! IendB Ending nest tile in the I-direction (RHO-, V-obc) !
1091! IstrM Starting nest tile in the I-direction (PSI-, U-obc) !
1092! IstrP Starting nest tile in the I-direction (PSI, U-points) !
1093! IendP Ending nest tile in the I-direction (PSI) !
1094! IstrT Starting nest tile in the I-direction (RHO-points) !
1095! IendT Ending nest tile in the I-direction (RHO_points) !
1096! !
1097! JstrB Starting nest tile in the J-direction (RHO-, U-obc) !
1098! JendB Ending nest tile in the J-direction (RHO-, U-obc) !
1099! JstrM Starting nest tile in the J-direction (PSI-, V-obc) !
1100! JstrP Starting nest tile in the J-direction (PSI, V-points) !
1101! JendP Ending nest tile in the J-direction (PSI) !
1102! JstrT Starting nest tile in the J-direction (RHO-points) !
1103! JendT Ending nest tile in the J-direction (RHO-points) !
1104! !
1105! Istrm3 Starting private I-halo computations, Istr-3 !
1106! Istrm2 Starting private I-halo computations, Istr-2 !
1107! Istrm1 Starting private I-halo computations, Istr-1 !
1108! IstrUm2 Starting private I-halo computations, IstrU-2 !
1109! IstrUm1 Starting private I-halo computations, IstrU-1 !
1110! Iendp1 Ending private I-halo computations, Iend+1 !
1111! Iendp2 Ending private I-halo computations, Iend+2 !
1112! Iendp2i Ending private I-halo computations, Iend+2 interior !
1113! Iendp3 Ending private I-halo computations, Iend+3 !
1114! !
1115! Jstrm3 Starting private J-halo computations, Jstr-3 !
1116! Jstrm2 Starting private J-halo computations, Jstr-2 !
1117! Jstrm1 Starting private J-halo computations, Jstr-1 !
1118! JstrVm2 Starting private J-halo computations, JstrV-2 !
1119! JstrVm1 Starting private J-halo computations, JstrV-1 !
1120! Jendp1 Ending private J-halo computations, Jend+1 !
1121! Jendp2 Ending private J-halo computations, Jend+2 !
1122! Jendp2i Ending private J-halo computations, Jend+2 interior !
1123! Jendp3 Ending private J-halo computations, Jend+3 !
1124! !
1125!======================================================================!
1126!
1127! Imported variable declarations.
1128!
1129 integer, intent(in) :: ng, tile
1130 integer, intent(in) :: my_Lm, my_Mm
1131 integer, intent(in) :: my_Istr, my_Iend, my_Jstr, my_Jend
1132!
1133 integer, intent(out) :: Iend, Istr, Jend, Jstr
1134 integer, intent(out) :: IstrM, IstrR, IstrU, IendR
1135 integer, intent(out) :: JstrM, JstrR, JstrV, JendR
1136 integer, intent(out) :: IstrB, IendB, IstrP, IendP, IstrT, IendT
1137 integer, intent(out) :: JstrB, JendB, JstrP, JendP, JstrT, JendT
1138 integer, intent(out) :: Istrm3, Istrm2, Istrm1, IstrUm2, IstrUm1
1139 integer, intent(out) :: Iendp1, Iendp2, Iendp2i, Iendp3
1140 integer, intent(out) :: Jstrm3, Jstrm2, Jstrm1, JstrVm2, JstrVm1
1141 integer, intent(out) :: Jendp1, Jendp2, Jendp2i, Jendp3
1142!
1143!=======================================================================
1144! Compute lower and upper bounds over a particular domain partition or
1145! tile for RHO-, U-, and V-variables.
1146!=======================================================================
1147!
1148! ROMS uses at staggered stencil:
1149!
1150! -------v(i,j+1,k)------- ------W(i,j,k)-------
1151! | | | |
1152! u(i,j,k) r(i,j,k) u(i+1,j,k) | r(i,j,k) |
1153! | | | |
1154! --------v(i,j,k)-------- -----W(i,j,k-1)------
1155!
1156! horizontal stencil vertical stencil
1157! C-grid
1158!
1159!
1160! M r..u..r..u..r..u..r..u..r..u..r..u..r..u..r..u..r..u..r..u..r
1161! : :
1162! M v p++v++p++v++p++v++p++v++p++v++p++v++p++v++p++v++p++v++p v
1163! : + | | | | | | | | + :
1164! Mm r u r u r u r u r u r u r u r u r u r u r
1165! : + | | | | | | | | + :
1166! Mm v p--v--p--v--p--v--p--v--p--v--p--v--p--v--p--v--p--v--p v
1167! : + | | | | | | | | + :
1168! r u r u r u r u r u r u r u r u r u r u r
1169! : + | | | | | | | | + :
1170! v p--v--p--v--p--v--p--v--p--v--p--v--p--v--p--v--p--v--p v
1171! : + | | | | | | | | + :
1172! r u r u r u r u r u r u r u r u r u r u r
1173! : + | | | | | | | | + :
1174! v p--v--p--v--p--v--p--v--p--v--p--v--p--v--p--v--p--v--p v
1175! : + | | | | | | | | + :
1176! 2 r u r u r u r u r u r u r u r u r u r u r
1177! : + | | | | | | | | + :
1178! 2 v p--v--p--v--p--v--p--v--p--v--p--v--p--v--p--v--p--v--p v
1179! : + | | | | | | | | + :
1180! 1 r u r u r u r u r u r u r u r u r u r u r
1181! : + | | | | | | | | + :
1182! 1 v p++v++p++v++p++v++p++v++p++v++p++v++p++v++p++v++p++v++p v
1183! : :
1184! 0 r..u..r..u..r..u..r..u..r..u..r..u..r..u..r..u..r..u..r..u..r
1185! 1 2 Lm L
1186! 0 1 2 Lm L
1187!
1188! interior Boundary Conditions
1189! computations W E S N
1190!
1191! RH0-type variables: [1:Lm, 1:Mm] [0,:] [L,:] [:,0] [:,M]
1192! PSI-type variables: [2:Lm, 2:Mm] [1,:] [L,:] [:,1] [:,M]
1193! U-type variables: [2:Lm, 1:Mm] [1,:] [L,:] [:,0] [:,M]
1194! V-type variables: [1:Lm, 2:Mm] [0,:] [L,:] [:,1] [:,M]
1195!
1196! Compute derived bounds for the loop indices over a subdomain tile.
1197! The extended bounds (labelled by suffix R) are designed to cover
1198! also the outer grid points (outlined above with :), if the subdomain
1199! tile is adjacent to the physical boundary (outlined above with +).
1200! Notice that IstrR, IendR, JstrR, JendR tile bounds computed here
1201! DO NOT COVER ghost points (outlined below with *) associated with
1202! periodic boundaries (if any) or the computational margins of MPI
1203! subdomains.
1204!
1205! Left/Top Tile Right/Top Tile
1206!
1207! JendR r..u..r..u..r..u..r..u * * * * u..r..u..r..u..r..u..r
1208! : Istr Iend Istr Iend :
1209! v p++v++p++v++p++v++p * * Jend * * p++v++p++v++p++v++p v
1210! : + | | | | | | + :
1211! r u r u r u r u * * * * u r u r u r u r
1212! : + | | | | | | + :
1213! v p--v--p--v--p--v--p * * * * p--v--p--v--p--v--p v
1214! : + | | | | | | + :
1215! r u r u r u r u * * * * u r u r u r u r
1216! : + | | | | | | + :
1217! v p--v--p--v--p--v--p * * Jstr * * p--v--p--v--p--v--p v
1218!
1219! * * * * * * * * * * * * * * * * * * * *
1220!
1221! * * * * * * * * * * * * * * * * * * * *
1222!
1223! IstrR IstrU IendR
1224! IstrM
1225!
1226! IstrB=Istr
1227! * * * * * * * * * * * IstrM=Istr
1228! Ghost Points IstrP=Istr
1229! * * * * * * * * * * * IstrR=Istr
1230! IstrT=Istr
1231! * * p--v--p--v--p--v--p * * Jend IstrU=Istr
1232! | | | | IendB=Iend
1233! Interior * * u r u r u r u * * IendP=Iend
1234! Tile | | | | IendR=Iend
1235! * * p--v--p--v--p--v--p * * IendT=Iend
1236! | | | | JstrB=Jstr
1237! * * u r u r u r u * * JstrM=Jstr
1238! | | | | JstrP=Jstr
1239! * * p--v--p--v--p--v--p * * Jstr JstrR=Jstr
1240! JstrT=Jstr
1241! * * * * * * * * * * * JstrV=Jstr
1242! JendB=Jend
1243! * * * * * * * * * * * JendP=Jend
1244! JendR=Jend
1245! Istr Iend JendT=Jend
1246!
1247!
1248!
1249! * * * * * * * * * * * * * * * * * * * *
1250!
1251! * * * * * * * * * * * * * * * * * * * *
1252! Istr Iend
1253! v p--v--p--v--p--v--p * * Jend * * p--v--p--v--p--v--p v
1254! : + | | | | | | + :
1255! r u r u r u r u * * * * u r u r u r u r
1256! : + | | | | | | + :
1257! JstrV v p--v--p--v--p--v--p * * * * p--v--p--v--p--v--p v JstrM
1258! : + | | | | | | + :
1259! r u r u r u r u * * * * u r u r u r u r
1260! : + | | | | | | + :
1261! v p++v++p++v++p++v++p * * Jstr * * p++v++p++v++p++v++p v
1262! : :
1263! r..u..r..u..r..u..r..u * * * * u..r..u..r..u..r..u..r
1264!
1265! IstrR IstrU IendR
1266! IstrM
1267!
1268! Left/Bottom Tile Right/Bottom Tile
1269!
1270!
1271! It also computes loop-bounds for U- and V-type variables which
1272! belong to the interior of the computational domain. These are
1273! labelled by suffixes U,V and they step one grid point inward from
1274! the side of the subdomain adjacent to the physical boundary.
1275! Conversely, for an internal subdomain which does not include a
1276! segment of the physical boundary, all bounds with suffixes R,U,V
1277! are set to the same values of corresponding non-suffixed bounds.
1278!
1279! In nested grids there are additional indices to process the
1280! overlap regions and ranges when calling lateral boundary
1281! conditions routines (see diagrams below):
1282!
1283! IstrT starting overlap I-direction (RHO-points)
1284! IendT ending overlap I-direction (RHO-points)
1285! JstrT starting overlap J-direction (RHO-points)
1286! JendT ending overlap J-direction (RHO-points)
1287!
1288! IstrP starting overlap I-direction (PSI-, U-points)
1289! IendP ending overlap I-direction (PSI-points)
1290! JstrP starting overlap J-direction (PSI-, V-points)
1291! JendP ending overlap J-direction (PSI-points)
1292!
1293! IstrB starting boundary I-direction (RHO-, V-points), IstrT+1
1294! IendB ending boundary I-direction (RHO-, V-points), IendT-1
1295! JstrB starting boundary J-direction (RHO-, U-points), JstrT+1
1296! JendB ending boundary J-direction (RHO-, U-points), JendT-1
1297!
1298! IstrM starting boundary I-direction (PSI-, U-points), IstrP+1
1299! JstrM starting boundary J-direction (PSI-, V-points), JstrP+1
1300!
1301! If not nesting, these indices are set to:
1302!
1303! IstrT = IstrR IstrP = Istr IstrB = Istr IstrM = IstrU
1304! IendT = IendR IendP = Iend IendB = Iend
1305! JstrT = JstrR JstrP = Jstr JstrB = Jstr JstrM = JstrV
1306! JendT = JendR JendP = Jend JendB = Jend
1307!
1308! The following diagram shows the lower left corner in nested grids:
1309!
1310! +
1311! r u r u r u r u r u r u r u r u r u r
1312! : : : + | | | | |
1313! v..p..v..p..v..p..v p--v--p--v--p--v--p--v--p--v--p--v
1314! : : : + | | | | |
1315! 2 r u r u r u r u r u r u r u r u r u r
1316! : : : + | | | | |
1317! 2 v..p..v..p..v..p..v p--v--p--v--p--v--p--v--p--v--p--v JstrV
1318! : : : + | | | | |
1319! 1 r u r u r u r u r u r u r u r u r u r Jstr
1320! : : : + | | | | |
1321! 1 v..p..v..p..v..p..v p++v++p++v++p++v++p++v++p++v++p++v++
1322! : : :
1323! 0 r u r u r u r u r u r u r u r u r u r JstrR
1324! : : : : : : : : :
1325! 0 v..p..v..p..v..p..v..p..v..p..v..p..v..p..v..p..v..p..v
1326! : : : : : : : : :
1327! -1 r u r u r u r u r u r u r u r u r u r
1328! : : : : : : : : :
1329! -1 v..p..v..p..v..p..v..p..v..p..v..p..v..p..v..p..v..p..v JstrM
1330! : : : : : : : : :
1331! -2 r u r u r u r u r u r u r u r u r u r JstrB
1332! : : : : : : : : :
1333! -2 v..p..v..p..v..p..v..p..v..p..v..p..v..p..v..p..v..p..v JstrP
1334! : : : : : : : : :
1335! -3 r u r u r u r u r u r u r u r u r u r JstrT
1336!
1337! -2 -1 0 1 2
1338! -3 -2 -1 0 1 2
1339!
1340! IstrP IstrM IstrU
1341! IstrB Istr
1342! IstrT IstrR
1343!
1344!
1345! The following diagram shows the upper right in nested grids:
1346!
1347!
1348! IendR IendT
1349! Iend IendB
1350! IendP
1351!
1352! Lm L L+1 L+2
1353! Lm L L+1 L+2
1354!
1355! M+2 r u r u r u r u r u r u r u r u r JendT
1356! : : : : : : : :
1357! M+1 v p..v..p..v..p..v..p..v..p..v..p..v..p..v..p..v JendP
1358! : : : : : : : :
1359! M+1 r u r u r u r u r u r u r u r u r JendB
1360! : : : : : : : :
1361! M+1 v p..v..p..v..p..v..p..v..p..v..p..v..p..v..p..v
1362! : : : : : : : :
1363! M r u r u r u r u r u r u r u r u r JendR
1364! : :
1365! M ++v++p++v++p++v++p++v++p++v++p++v++p v..p..v..p..v
1366! | | | | | + : :
1367! Mm r u r u r u r u r u r u r u r u r Jend
1368! | | | | | + : :
1369! Mm v--p--v--p--v--p--v--p--v--p--v--p v..p..v..p..v
1370! | | | | | + : :
1371! r u r u r u r u r u r u r u r u r
1372! | | | | | + : :
1373! v--p--v--p--v--p--v--p--v--p--v--p v..p..v..p..v
1374! | | | | | + : :
1375! r u r u r u r u r u r u r u r u r
1376! +
1377!
1378!-----------------------------------------------------------------------
1379! Starting I-tile indices.
1380!-----------------------------------------------------------------------
1381!
1382 IF (domain(ng)%Western_Edge(tile)) THEN
1383 IF (ewperiodic(ng)) THEN
1384 istr =my_istr
1385 istrp=my_istr
1386 istrr=my_istr
1387 istrt=istrr
1388 istru=my_istr
1389 istrb=my_istr
1390 istrm=istru
1391#ifdef NESTING
1392 ELSE IF (compositegrid(iwest,ng)) THEN
1393 istr =my_istr
1394 istrp=-nghostpoints+1
1395 istrr=my_istr
1396 istrt=-nghostpoints
1397 istru=my_istr
1398 istrb=istrt ! boundary conditions maybe avoided
1399 istrm=istrp ! boundary conditions maybe avoided
1400 ELSE IF (refinedgrid(ng).and.(refinescale(ng).gt.0)) THEN
1401 istr =my_istr
1402 istrr=my_istr-1
1403 istru=my_istr+1
1404 istrp=-nghostpoints+1
1405 istrt=-nghostpoints
1406 istrb=istrt ! boundary conditions are avoided
1407 istrm=istrp ! boundary conditions are avoided
1408#endif
1409 ELSE
1410 istr =my_istr
1411 istrp=my_istr
1412 istrr=my_istr-1
1413 istrt=istrr
1414 istru=my_istr+1
1415 istrb=istrt+1
1416 istrm=istrp+1
1417 END IF
1418 ELSE
1419 istr =my_istr
1420 istrp=my_istr
1421 istrr=my_istr
1422 istrt=istrr
1423 istru=my_istr
1424 istrb=my_istr
1425 istrm=istru
1426 END IF
1427!
1428! Special case, Istrm3: used when MAX(0,Istr-3) is needed.
1429!
1430 IF (domain(ng)%Western_Edge(tile)) THEN
1431 IF (ewperiodic(ng)) THEN
1432 istrm3=my_istr-3
1433#ifdef NESTING
1434 ELSE IF (compositegrid(iwest,ng)) THEN
1435 istrm3=my_istr-3
1436 ELSE IF (refinedgrid(ng).and.(refinescale(ng).gt.0)) THEN
1437 istrm3=my_istr-3
1438#endif
1439 ELSE
1440 istrm3=max(0,my_istr-3)
1441 END IF
1442 ELSE
1443 istrm3=my_istr-3
1444 END IF
1445!
1446! Special case, Istrm2: used when MAX(0,Istr-2) is needed.
1447!
1448 IF (domain(ng)%Western_Edge(tile)) THEN
1449 IF (ewperiodic(ng)) THEN
1450 istrm2=my_istr-2
1451#ifdef NESTING
1452 ELSE IF (compositegrid(iwest,ng)) THEN
1453 istrm2=my_istr-2
1454 ELSE IF (refinedgrid(ng).and.(refinescale(ng).gt.0)) THEN
1455 istrm2=my_istr-2
1456#endif
1457 ELSE
1458 istrm2=max(0,my_istr-2)
1459 END IF
1460 ELSE
1461 istrm2=my_istr-2
1462 END IF
1463!
1464! Special case, IstrUm2: used when MAX(1,IstrU-2) is needed.
1465!
1466 IF (domain(ng)%Western_Edge(tile)) THEN
1467 IF (ewperiodic(ng)) THEN
1468 istrum2=istru-2
1469#ifdef NESTING
1470 ELSE IF (compositegrid(iwest,ng)) THEN
1471 istrum2=istru-2
1472 ELSE IF (refinedgrid(ng).and.(refinescale(ng).gt.0)) THEN
1473 istrum2=istru-2
1474#endif
1475 ELSE
1476 istrum2=max(1,istru-2)
1477 END IF
1478 ELSE
1479 istrum2=istru-2
1480 END IF
1481!
1482! Special case, Istrm1: used when MAX(1,Istr-1) is needed.
1483!
1484 IF (domain(ng)%Western_Edge(tile)) THEN
1485 IF (ewperiodic(ng)) THEN
1486 istrm1=my_istr-1
1487#ifdef NESTING
1488 ELSE IF (compositegrid(iwest,ng)) THEN
1489 istrm1=my_istr-1
1490 ELSE IF (refinedgrid(ng).and.(refinescale(ng).gt.0)) THEN
1491 istrm1=my_istr-1
1492#endif
1493 ELSE
1494 istrm1=max(1,my_istr-1)
1495 END IF
1496 ELSE
1497 istrm1=my_istr-1
1498 END IF
1499!
1500! Special case, IstrUm1: used when MAX(2,IstrU-1) is needed.
1501!
1502 IF (domain(ng)%Western_Edge(tile)) THEN
1503 IF (ewperiodic(ng)) THEN
1504 istrum1=istru-1
1505#ifdef NESTING
1506 ELSE IF (compositegrid(iwest,ng)) THEN
1507 istrum1=istru-1
1508 ELSE IF (refinedgrid(ng).and.(refinescale(ng).gt.0)) THEN
1509 istrum1=istru-1
1510#endif
1511 ELSE
1512 istrum1=max(2,istru-1)
1513 END IF
1514 ELSE
1515 istrum1=istru-1
1516 END IF
1517!
1518!-----------------------------------------------------------------------
1519! Ending I-tile indices.
1520!-----------------------------------------------------------------------
1521!
1522 IF (domain(ng)%Eastern_Edge(tile)) THEN
1523 IF (ewperiodic(ng)) THEN
1524 iend =my_iend
1525 iendr=my_iend
1526 iendp=iendr ! check this one, same as IendR?
1527 iendt=iendr
1528 iendb=my_iend
1529#ifdef NESTING
1530 ELSE IF (compositegrid(ieast,ng)) THEN
1531 iend =my_iend
1532 iendr=my_iend
1533 iendp=my_iend+nghostpoints ! check this one, same as IendT?
1534 iendt=my_iend+nghostpoints
1535 iendb=iendt ! boundary conditions maybe avoided
1536 ELSE IF (refinedgrid(ng).and.(refinescale(ng).gt.0)) THEN
1537 iend =my_iend
1538 iendr=my_iend+1
1539 iendp=my_iend+nghostpoints ! check this one, same as IendT
1540 iendt=my_iend+nghostpoints
1541 iendb=iendt ! boundary conditions are avoided
1542#endif
1543 ELSE
1544 iend =my_iend
1545 iendr=my_iend+1
1546 iendp=iendr
1547 iendt=iendr
1548 iendb=iendt-1
1549 END IF
1550 ELSE
1551 iend =my_iend
1552 iendr=my_iend
1553 iendp=iendr
1554 iendt=iendr
1555 iendb=my_iend
1556 END IF
1557!
1558! Special case, Iendp1: used when MIN(Iend+1,my_Lm) is needed.
1559!
1560 IF (domain(ng)%Eastern_Edge(tile)) THEN
1561 IF (ewperiodic(ng)) THEN
1562 iendp1=my_iend+1
1563#ifdef NESTING
1564 ELSE IF (compositegrid(ieast,ng)) THEN
1565 iendp1=my_iend+1
1566 ELSE IF (refinedgrid(ng).and.(refinescale(ng).gt.0)) THEN
1567 iendp1=my_iend+1
1568#endif
1569 ELSE
1570 iendp1=min(my_iend+1,my_lm)
1571 END IF
1572 ELSE
1573 iendp1=my_iend+1
1574 END IF
1575!
1576! Special case, Iendp2i: used when MIN(Iend+2,my_Lm) is needed.
1577!
1578 IF (domain(ng)%Eastern_Edge(tile)) THEN
1579 IF (ewperiodic(ng)) THEN
1580 iendp2i=my_iend+2
1581#ifdef NESTING
1582 ELSE IF (compositegrid(ieast,ng)) THEN
1583 iendp2i=my_iend+2
1584 ELSE IF (refinedgrid(ng).and.(refinescale(ng).gt.0)) THEN
1585 iendp2i=my_iend+2
1586#endif
1587 ELSE
1588 iendp2i=min(my_iend+2,my_lm)
1589 END IF
1590 ELSE
1591 iendp2i=my_iend+2
1592 END IF
1593!
1594! Special case, Iendp2: used when MIN(Iend+2,my_Lm+1) is needed.
1595!
1596 IF (domain(ng)%Eastern_Edge(tile)) THEN
1597 IF (ewperiodic(ng)) THEN
1598 iendp2=my_iend+2
1599#ifdef NESTING
1600 ELSE IF (compositegrid(ieast,ng)) THEN
1601 iendp2=my_iend+2
1602 ELSE IF (refinedgrid(ng).and.(refinescale(ng).gt.0)) THEN
1603 iendp2=my_iend+2
1604#endif
1605 ELSE
1606 iendp2=min(my_iend+2,my_lm+1)
1607 END IF
1608 ELSE
1609 iendp2=my_iend+2
1610 END IF
1611!
1612! Special case, Iendp3: used when MIN(Iend+3,my_Lm+1) is needed.
1613!
1614 IF (domain(ng)%Eastern_Edge(tile)) THEN
1615 IF (ewperiodic(ng)) THEN
1616 iendp3=my_iend+3
1617#ifdef NESTING
1618 ELSE IF (compositegrid(ieast,ng)) THEN
1619 iendp3=my_iend+3
1620 ELSE IF (refinedgrid(ng).and.(refinescale(ng).gt.0)) THEN
1621 iendp3=my_iend+3
1622#endif
1623 ELSE
1624 iendp3=min(my_iend+3,my_lm+1)
1625 END IF
1626 ELSE
1627 iendp3=my_iend+3
1628 END IF
1629!
1630!-----------------------------------------------------------------------
1631! Starting J-tile indices.
1632!-----------------------------------------------------------------------
1633!
1634 IF (domain(ng)%Southern_Edge(tile)) THEN
1635 IF (nsperiodic(ng)) THEN
1636 jstr =my_jstr
1637 jstrp=my_jstr
1638 jstrr=my_jstr
1639 jstrt=jstrr
1640 jstrv=my_jstr
1641 jstrb=my_jstr
1642 jstrm=jstrv
1643#ifdef NESTING
1644 ELSE IF (compositegrid(isouth,ng)) THEN
1645 jstr =my_jstr
1646 jstrp=-nghostpoints+1
1647 jstrr=my_jstr
1648 jstrt=-nghostpoints
1649 jstrv=my_jstr
1650 jstrb=jstrt ! boundary conditions maybe avoided
1651 jstrm=jstrp ! boundary conditions maybe avoided
1652 ELSE IF (refinedgrid(ng).and.(refinescale(ng).gt.0)) THEN
1653 jstr =my_jstr
1654 jstrr=my_jstr-1
1655 jstrv=my_jstr+1
1656 jstrp=-nghostpoints+1
1657 jstrt=-nghostpoints
1658 jstrb=jstrt ! boundary conditions are avoided
1659 jstrm=jstrp ! boundary conditions are avoided
1660#endif
1661 ELSE
1662 jstr =my_jstr
1663 jstrp=my_jstr
1664 jstrr=my_jstr-1
1665 jstrt=jstrr
1666 jstrv=my_jstr+1
1667 jstrb=jstrt+1
1668 jstrm=jstrp+1
1669 END IF
1670 ELSE
1671 jstr =my_jstr
1672 jstrp=my_jstr
1673 jstrr=my_jstr
1674 jstrt=jstrr
1675 jstrv=my_jstr
1676 jstrb=my_jstr
1677 jstrm=jstrv
1678 END IF
1679!
1680! Special case, Jstrm3: used when MAX(0,Jstr-3) is needed.
1681!
1682 IF (domain(ng)%Southern_Edge(tile)) THEN
1683 IF (nsperiodic(ng)) THEN
1684 jstrm3=my_jstr-3
1685#ifdef NESTING
1686 ELSE IF (compositegrid(isouth,ng)) THEN
1687 jstrm3=my_jstr-3
1688 ELSE IF (refinedgrid(ng).and.(refinescale(ng).gt.0)) THEN
1689 jstrm3=my_jstr-3
1690#endif
1691 ELSE
1692 jstrm3=max(0,my_jstr-3)
1693 END IF
1694 ELSE
1695 jstrm3=my_jstr-3
1696 END IF
1697!
1698! Special case, Jstrm2: used when MAX(0,Jstr-2) is needed.
1699!
1700 IF (domain(ng)%Southern_Edge(tile)) THEN
1701 IF (nsperiodic(ng)) THEN
1702 jstrm2=my_jstr-2
1703#ifdef NESTING
1704 ELSE IF (compositegrid(isouth,ng)) THEN
1705 jstrm2=my_jstr-2
1706 ELSE IF (refinedgrid(ng).and.(refinescale(ng).gt.0)) THEN
1707 jstrm2=my_jstr-2
1708#endif
1709 ELSE
1710 jstrm2=max(0,my_jstr-2)
1711 END IF
1712 ELSE
1713 jstrm2=my_jstr-2
1714 END IF
1715!
1716! Special case, JstrVm2: used when MAX(1,JstrV-2) is needed.
1717!
1718 IF (domain(ng)%Southern_Edge(tile)) THEN
1719 IF (nsperiodic(ng)) THEN
1720 jstrvm2=jstrv-2
1721#ifdef NESTING
1722 ELSE IF (compositegrid(isouth,ng)) THEN
1723 jstrvm2=jstrv-2
1724 ELSE IF (refinedgrid(ng).and.(refinescale(ng).gt.0)) THEN
1725 jstrvm2=jstrv-2
1726#endif
1727 ELSE
1728 jstrvm2=max(1,jstrv-2)
1729 END IF
1730 ELSE
1731 jstrvm2=jstrv-2
1732 END IF
1733!
1734! Special case, Jstrm1: used when MAX(1,Jstr-1) is needed.
1735!
1736 IF (domain(ng)%Southern_Edge(tile)) THEN
1737 IF (nsperiodic(ng)) THEN
1738 jstrm1=my_jstr-1
1739#ifdef NESTING
1740 ELSE IF (compositegrid(isouth,ng)) THEN
1741 jstrm1=my_jstr-1
1742 ELSE IF (refinedgrid(ng).and.(refinescale(ng).gt.0)) THEN
1743 jstrm1=my_jstr-1
1744#endif
1745 ELSE
1746 jstrm1=max(1,my_jstr-1)
1747 END IF
1748 ELSE
1749 jstrm1=my_jstr-1
1750 END IF
1751!
1752! Special case, JstrVm1: used when MAX(2,JstrV-1) is needed.
1753!
1754 IF (domain(ng)%Southern_Edge(tile)) THEN
1755 IF (nsperiodic(ng)) THEN
1756 jstrvm1=jstrv-1
1757#ifdef NESTING
1758 ELSE IF (compositegrid(isouth,ng)) THEN
1759 jstrvm1=jstrv-1
1760 ELSE IF (refinedgrid(ng).and.(refinescale(ng).gt.0)) THEN
1761 jstrvm1=jstrv-1
1762#endif
1763 ELSE
1764 jstrvm1=max(2,jstrv-1)
1765 END IF
1766 ELSE
1767 jstrvm1=jstrv-1
1768 END IF
1769!
1770!-----------------------------------------------------------------------
1771! Ending J-tile indices.
1772!-----------------------------------------------------------------------
1773!
1774 IF (domain(ng)%Northern_Edge(tile)) THEN
1775 IF (nsperiodic(ng)) THEN
1776 jend =my_jend
1777 jendr=my_jend
1778 jendp=jendr ! check this one, same as JendR?
1779 jendt=jendr
1780 jendb=my_jend
1781#ifdef NESTING
1782 ELSE IF (compositegrid(inorth,ng)) THEN
1783 jend =my_jend
1784 jendr=my_jend
1785 jendp=my_jend+nghostpoints ! check this one, same as JendT?
1786 jendt=my_jend+nghostpoints
1787 jendb=jendt ! boundary conditions maybe avoided
1788 ELSE IF (refinedgrid(ng).and.(refinescale(ng).gt.0)) THEN
1789 jend =my_jend
1790 jendr=my_jend+1
1791 jendp=my_jend+nghostpoints ! check this one, same as JendT
1792 jendt=my_jend+nghostpoints
1793 jendb=jendt ! boundary conditions are avoided
1794#endif
1795 ELSE
1796 jend =my_jend
1797 jendr=my_jend+1
1798 jendp=jendr
1799 jendt=jendr
1800 jendb=jendt-1
1801 END IF
1802 ELSE
1803 jend =my_jend
1804 jendr=my_jend
1805 jendp=jendr
1806 jendp=jendr
1807 jendt=jendr
1808 jendb=my_jend
1809 END IF
1810!
1811! Special case, Jendp1: used when MIN(Jend+1,my_Mm) is needed.
1812!
1813 IF (domain(ng)%Northern_Edge(tile)) THEN
1814 IF (nsperiodic(ng)) THEN
1815 jendp1=my_jend+1
1816#ifdef NESTING
1817 ELSE IF (compositegrid(inorth,ng)) THEN
1818 jendp1=my_jend+1
1819 ELSE IF (refinedgrid(ng).and.(refinescale(ng).gt.0)) THEN
1820 jendp1=my_jend+1
1821#endif
1822 ELSE
1823 jendp1=min(my_jend+1,my_mm)
1824 END IF
1825 ELSE
1826 jendp1=my_jend+1
1827 END IF
1828!
1829! Special case, Jendp2i: used when MIN(Jend+2,my_Mm) is needed.
1830!
1831 IF (domain(ng)%Northern_Edge(tile)) THEN
1832 IF (nsperiodic(ng)) THEN
1833 jendp2i=my_jend+2
1834#ifdef NESTING
1835 ELSE IF (compositegrid(inorth,ng)) THEN
1836 jendp2i=my_jend+2
1837 ELSE IF (refinedgrid(ng).and.(refinescale(ng).gt.0)) THEN
1838 jendp2i=my_jend+2
1839#endif
1840 ELSE
1841 jendp2i=min(my_jend+2,my_mm)
1842 END IF
1843 ELSE
1844 jendp2i=my_jend+2
1845 END IF
1846!
1847! Special case, Jendp2: used when MIN(Jend+2,my_Mm+1) is needed.
1848!
1849 IF (domain(ng)%Northern_Edge(tile)) THEN
1850 IF (nsperiodic(ng)) THEN
1851 jendp2=my_jend+2
1852#ifdef NESTING
1853 ELSE IF (compositegrid(inorth,ng)) THEN
1854 jendp2=my_jend+2
1855 ELSE IF (refinedgrid(ng).and.(refinescale(ng).gt.0)) THEN
1856 jendp2=my_jend+2
1857#endif
1858 ELSE
1859 jendp2=min(my_jend+2,my_mm+1)
1860 END IF
1861 ELSE
1862 jendp2=my_jend+2
1863 END IF
1864!
1865! Special case, Jendp3: used when MIN(Jend+3,my_Mm+1) is needed.
1866!
1867 IF (domain(ng)%Northern_Edge(tile)) THEN
1868 IF (nsperiodic(ng)) THEN
1869 jendp3=my_jend+3
1870#ifdef NESTING
1871 ELSE IF (compositegrid(inorth,ng)) THEN
1872 jendp3=my_jend+3
1873 ELSE IF (refinedgrid(ng).and.(refinescale(ng).gt.0)) THEN
1874 jendp3=my_jend+3
1875#endif
1876 ELSE
1877 jendp3=min(my_jend+3,my_mm+1)
1878 END IF
1879 ELSE
1880 jendp3=my_jend+3
1881 END IF
1882!
1883 RETURN

References mod_scalars::compositegrid, mod_param::domain, mod_scalars::ewperiodic, mod_scalars::ieast, mod_scalars::inorth, mod_scalars::isouth, mod_scalars::iwest, mod_param::nghostpoints, mod_scalars::nsperiodic, mod_scalars::refinedgrid, and mod_scalars::refinescale.

Referenced by get_tile().

Here is the caller graph for this function: