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

Functions/Subroutines

subroutine, public tile_indices (model, my_im, my_jm, my_lm, my_mm, my_bounds, my_domain, my_iobounds)
 
subroutine, public tile_obs_bounds (model, my_im, my_jm, my_lm, my_mm, my_domain)
 

Function/Subroutine Documentation

◆ tile_indices()

subroutine, public tile_indices_mod::tile_indices ( integer, intent(in) model,
integer, dimension(1:ngrids), intent(in) my_im,
integer, dimension(1:ngrids), intent(in) my_jm,
integer, dimension(1:ngrids), intent(in) my_lm,
integer, dimension(1:ngrids), intent(in) my_mm,
type (t_bounds), dimension(1:ngrids), intent(out) my_bounds,
type (t_domain), dimension(1:ngrids), intent(inout) my_domain,
type (t_iobounds), dimension(1:ngrids), intent(inout) my_iobounds )

Definition at line 62 of file tile_indices.F.

66!***********************************************************************
67!
68! Imported variable declarations.
69!
70 integer, intent(in) :: model
71 integer, intent(in) :: my_Im(1:Ngrids), my_Jm(1:Ngrids)
72 integer, intent(in) :: my_Lm(1:Ngrids), my_Mm(1:Ngrids)
73!
74 TYPE (T_BOUNDS), intent(out) :: my_BOUNDS(1:Ngrids)
75 TYPE (T_DOMAIN), intent(inout) :: my_DOMAIN(1:Ngrids)
76 TYPE (T_IOBOUNDS), intent(inout) :: my_IOBOUNDS(1:Ngrids)
77!
78! Local variable declarations.
79!
80 integer :: Itile, Jtile, Nghost
81 integer :: ng, tile
82!
83!-----------------------------------------------------------------------
84! Set lower and upper bounds indices per domain partition for all
85! nested grids.
86!-----------------------------------------------------------------------
87!
88! Set boundary edge I- or J-indices for each variable type.
89!
90 DO ng=1,ngrids
91 my_bounds(ng) % edge(iwest ,p2dvar) = 1
92 my_bounds(ng) % edge(iwest ,r2dvar) = 0
93 my_bounds(ng) % edge(iwest ,u2dvar) = 1
94 my_bounds(ng) % edge(iwest ,v2dvar) = 0
95
96 my_bounds(ng) % edge(ieast ,p2dvar) = my_lm(ng)+1
97 my_bounds(ng) % edge(ieast ,r2dvar) = my_lm(ng)+1
98 my_bounds(ng) % edge(ieast ,u2dvar) = my_lm(ng)+1
99 my_bounds(ng) % edge(ieast ,v2dvar) = my_lm(ng)+1
100
101 my_bounds(ng) % edge(isouth,p2dvar) = 1
102 my_bounds(ng) % edge(isouth,u2dvar) = 0
103 my_bounds(ng) % edge(isouth,r2dvar) = 0
104 my_bounds(ng) % edge(isouth,v2dvar) = 1
105
106 my_bounds(ng) % edge(inorth,p2dvar) = my_mm(ng)+1
107 my_bounds(ng) % edge(inorth,r2dvar) = my_mm(ng)+1
108 my_bounds(ng) % edge(inorth,u2dvar) = my_mm(ng)+1
109 my_bounds(ng) % edge(inorth,v2dvar) = my_mm(ng)+1
110 END DO
111!
112! Set logical switches needed when processing variables in tiles
113! adjacent to the domain boundary edges or corners. This needs to
114! be computed first since these switches are used in "get_tile".
115!
116 DO ng=1,ngrids
117 DO tile=-1,ntilei(ng)*ntilej(ng)-1
118 CALL get_domain_edges (ng, tile, &
119 & my_lm(ng), my_mm(ng), &
120 & my_domain(ng)% Eastern_Edge (tile), &
121 & my_domain(ng)% Western_Edge (tile), &
122 & my_domain(ng)% Northern_Edge (tile), &
123 & my_domain(ng)% Southern_Edge (tile), &
124 & my_domain(ng)% NorthEast_Corner(tile), &
125 & my_domain(ng)% NorthWest_Corner(tile), &
126 & my_domain(ng)% SouthEast_Corner(tile), &
127 & my_domain(ng)% SouthWest_Corner(tile), &
128 & my_domain(ng)% NorthEast_Test (tile), &
129 & my_domain(ng)% NorthWest_Test (tile), &
130 & my_domain(ng)% SouthEast_Test (tile), &
131 & my_domain(ng)% SouthWest_Test (tile))
132 END DO
133 END DO
134!
135! Set tile computational indices and arrays allocation bounds
136!
137 nghost=nghostpoints
138 DO ng=1,ngrids
139 my_bounds(ng) % LBij = 0
140 my_bounds(ng) % UBij = max(my_lm(ng)+1,my_mm(ng)+1)
141 DO tile=-1,ntilei(ng)*ntilej(ng)-1
142 my_bounds(ng) % tile(tile) = tile
143 CALL get_tile (ng, tile, &
144 & my_lm(ng), my_mm(ng), &
145 & itile, jtile, &
146 & my_bounds(ng)% Istr (tile), &
147 & my_bounds(ng)% Iend (tile), &
148 & my_bounds(ng)% Jstr (tile), &
149 & my_bounds(ng)% Jend (tile), &
150 & my_bounds(ng)% IstrM (tile), &
151 & my_bounds(ng)% IstrR (tile), &
152 & my_bounds(ng)% IstrU (tile), &
153 & my_bounds(ng)% IendR (tile), &
154 & my_bounds(ng)% JstrM (tile), &
155 & my_bounds(ng)% JstrR (tile), &
156 & my_bounds(ng)% JstrV (tile), &
157 & my_bounds(ng)% JendR (tile), &
158 & my_bounds(ng)% IstrB (tile), &
159 & my_bounds(ng)% IendB (tile), &
160 & my_bounds(ng)% IstrP (tile), &
161 & my_bounds(ng)% IendP (tile), &
162 & my_bounds(ng)% IstrT (tile), &
163 & my_bounds(ng)% IendT (tile), &
164 & my_bounds(ng)% JstrB (tile), &
165 & my_bounds(ng)% JendB (tile), &
166 & my_bounds(ng)% JstrP (tile), &
167 & my_bounds(ng)% JendP (tile), &
168 & my_bounds(ng)% JstrT (tile), &
169 & my_bounds(ng)% JendT (tile), &
170 & my_bounds(ng)% Istrm3 (tile), &
171 & my_bounds(ng)% Istrm2 (tile), &
172 & my_bounds(ng)% Istrm1 (tile), &
173 & my_bounds(ng)% IstrUm2(tile), &
174 & my_bounds(ng)% IstrUm1(tile), &
175 & my_bounds(ng)% Iendp1 (tile), &
176 & my_bounds(ng)% Iendp2 (tile), &
177 & my_bounds(ng)% Iendp2i(tile), &
178 & my_bounds(ng)% Iendp3 (tile), &
179 & my_bounds(ng)% Jstrm3 (tile), &
180 & my_bounds(ng)% Jstrm2 (tile), &
181 & my_bounds(ng)% Jstrm1 (tile), &
182 & my_bounds(ng)% JstrVm2(tile), &
183 & my_bounds(ng)% JstrVm1(tile), &
184 & my_bounds(ng)% Jendp1 (tile), &
185 & my_bounds(ng)% Jendp2 (tile), &
186 & my_bounds(ng)% Jendp2i(tile), &
187 & my_bounds(ng)% Jendp3 (tile))
188!
189 CALL get_bounds (ng, tile, 0, nghost, &
190 & my_im(ng), my_jm(ng), &
191 & my_lm(ng), my_mm(ng), &
192 & itile, jtile, &
193 & my_bounds(ng)% LBi(tile), &
194 & my_bounds(ng)% UBi(tile), &
195 & my_bounds(ng)% LBj(tile), &
196 & my_bounds(ng)% UBj(tile))
197 END DO
198 END DO
199!
200! Set I/O processing minimum (Imin, Jmax) and maximum (Imax, Jmax)
201! indices for non-overlapping (Nghost=0) and overlapping (Nghost>0)
202! tiles for each C-grid type variable.
203!
204 nghost=nghostpoints
205 DO ng=1,ngrids
206 DO tile=0,ntilei(ng)*ntilej(ng)-1
207 CALL get_bounds (ng, tile, p2dvar, 0 , &
208 & my_im(ng), my_jm(ng), &
209 & my_lm(ng), my_mm(ng), &
210 & itile, jtile, &
211 & my_bounds(ng)% Imin(1,0,tile), &
212 & my_bounds(ng)% Imax(1,0,tile), &
213 & my_bounds(ng)% Jmin(1,0,tile), &
214 & my_bounds(ng)% Jmax(1,0,tile))
215 CALL get_bounds (ng, tile, p2dvar, nghost, &
216 & my_im(ng), my_jm(ng), &
217 & my_lm(ng), my_mm(ng), &
218 & itile, jtile, &
219 & my_bounds(ng)% Imin(1,1,tile), &
220 & my_bounds(ng)% Imax(1,1,tile), &
221 & my_bounds(ng)% Jmin(1,1,tile), &
222 & my_bounds(ng)% Jmax(1,1,tile))
223!
224 CALL get_bounds (ng, tile, r2dvar, 0 , &
225 & my_im(ng), my_jm(ng), &
226 & my_lm(ng), my_mm(ng), &
227 & itile, jtile, &
228 & my_bounds(ng)% Imin(2,0,tile), &
229 & my_bounds(ng)% Imax(2,0,tile), &
230 & my_bounds(ng)% Jmin(2,0,tile), &
231 & my_bounds(ng)% Jmax(2,0,tile))
232 CALL get_bounds (ng, tile, r2dvar, nghost, &
233 & my_im(ng), my_jm(ng), &
234 & my_lm(ng), my_mm(ng), &
235 & itile, jtile, &
236 & my_bounds(ng)% Imin(2,1,tile), &
237 & my_bounds(ng)% Imax(2,1,tile), &
238 & my_bounds(ng)% Jmin(2,1,tile), &
239 & my_bounds(ng)% Jmax(2,1,tile))
240!
241 CALL get_bounds (ng, tile, u2dvar, 0 , &
242 & my_im(ng), my_jm(ng), &
243 & my_lm(ng), my_mm(ng), &
244 & itile, jtile, &
245 & my_bounds(ng)% Imin(3,0,tile), &
246 & my_bounds(ng)% Imax(3,0,tile), &
247 & my_bounds(ng)% Jmin(3,0,tile), &
248 & my_bounds(ng)% Jmax(3,0,tile))
249 CALL get_bounds (ng, tile, u2dvar, nghost, &
250 & my_im(ng), my_jm(ng), &
251 & my_lm(ng), my_mm(ng), &
252 & itile, jtile, &
253 & my_bounds(ng)% Imin(3,1,tile), &
254 & my_bounds(ng)% Imax(3,1,tile), &
255 & my_bounds(ng)% Jmin(3,1,tile), &
256 & my_bounds(ng)% Jmax(3,1,tile))
257!
258 CALL get_bounds (ng, tile, v2dvar, 0 , &
259 & my_im(ng), my_jm(ng), &
260 & my_lm(ng), my_mm(ng), &
261 & itile, jtile, &
262 & my_bounds(ng)% Imin(4,0,tile), &
263 & my_bounds(ng)% Imax(4,0,tile), &
264 & my_bounds(ng)% Jmin(4,0,tile), &
265 & my_bounds(ng)% Jmax(4,0,tile))
266 CALL get_bounds (ng, tile, v2dvar, nghost, &
267 & my_im(ng), my_jm(ng), &
268 & my_lm(ng), my_mm(ng), &
269 & itile, jtile, &
270 & my_bounds(ng)% Imin(4,1,tile), &
271 & my_bounds(ng)% Imax(4,1,tile), &
272 & my_bounds(ng)% Jmin(4,1,tile), &
273 & my_bounds(ng)% Jmax(4,1,tile))
274 END DO
275 END DO
276!
277! Set NetCDF IO bounds.
278!
279 DO ng=1,ngrids
280 CALL get_iobounds (ng, my_lm(ng), my_mm(ng), &
281 & my_bounds, my_iobounds)
282 END DO
283!
284 RETURN

References get_bounds_mod::get_bounds(), get_bounds_mod::get_domain_edges(), get_bounds_mod::get_iobounds(), get_bounds_mod::get_tile(), mod_scalars::ieast, mod_scalars::inorth, mod_scalars::isouth, mod_scalars::iwest, mod_param::nghostpoints, mod_param::ngrids, mod_param::ntilei, mod_param::ntilej, mod_param::p2dvar, mod_param::r2dvar, mod_param::u2dvar, and mod_param::v2dvar.

Referenced by inp_par_mod::inp_par().

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

◆ tile_obs_bounds()

subroutine, public tile_indices_mod::tile_obs_bounds ( integer, intent(in) model,
integer, dimension(1:ngrids), intent(in) my_im,
integer, dimension(1:ngrids), intent(in) my_jm,
integer, dimension(1:ngrids), intent(in) my_lm,
integer, dimension(1:ngrids), intent(in) my_mm,
type (t_domain), dimension(1:ngrids), intent(inout) my_domain )

Definition at line 288 of file tile_indices.F.

292!***********************************************************************
293!
294! Imported variable declarations.
295!
296 integer, intent(in) :: model
297 integer, intent(in) :: my_Im(1:Ngrids), my_Jm(1:Ngrids)
298 integer, intent(in) :: my_Lm(1:Ngrids), my_Mm(1:Ngrids)
299!
300 TYPE (T_DOMAIN), intent(inout) :: my_DOMAIN(1:Ngrids)
301!
302! Local variable declarations.
303!
304 integer :: Itile, Jtile, Uoff, Voff
305 integer :: ng, tile
306!
307 real(r8), parameter :: epsilon = 1.0e-8_r8
308!
309!-----------------------------------------------------------------------
310! Set minimum and maximum fractional coordinates for processing
311! observations. Either the full grid or only interior points will
312! be considered. The strategy here is to add a small value (epsilon)
313! to the eastern and northern boundary values of Xmax and Ymax so
314! observations at such boundaries locations are processed. This
315! is needed because the .lt. operator in the following conditional:
316!
317! IF (...
318! & ((Xmin.le.Xobs(iobs)).and.(Xobs(iobs).lt.Xmax)).and. &
319! & ((Ymin.le.Yobs(iobs)).and.(Yobs(iobs).lt.Ymax))) THEN
320!-----------------------------------------------------------------------
321!
322! Set RHO-points domain lower and upper bounds (integer).
323!
324 DO ng=1,ngrids
325#ifdef DISTRIBUTE
326 CALL get_bounds (ng, myrank, r2dvar, 0, &
327 & my_im(ng), my_jm(ng), &
328 & my_lm(ng), my_mm(ng), &
329 & itile, jtile, &
330 & rilb(ng), riub(ng), rjlb(ng), rjub(ng))
331# ifndef FULL_GRID
332 IF (itile.eq.0) THEN
333 rilb(ng)=rilb(ng)+1
334 END IF
335 IF (itile.eq.(ntilei(ng)-1)) THEN
336 riub(ng)=riub(ng)-1
337 END IF
338 IF (jtile.eq.0) THEN
339 rjlb(ng)=rjlb(ng)+1
340 END IF
341 IF (jtile.eq.(ntilej(ng)-1)) THEN
342 rjub(ng)=rjub(ng)-1
343 END IF
344# endif
345#else
346# ifdef FULL_GRID
347 rilb(ng)=0
348 riub(ng)=my_lm(ng)+1
349 rjlb(ng)=0
350 rjub(ng)=my_mm(ng)+1
351# else
352 rilb(ng)=1
353 riub(ng)=my_lm(ng)
354 rjlb(ng)=1
355 rjub(ng)=my_mm(ng)
356# endif
357#endif
358!
359! Minimum and maximum fractional coordinates for RHO-points.
360!
361 DO tile=0,ntilei(ng)*ntilej(ng)-1
362 CALL get_domain (ng, tile, r2dvar, 0, &
363 & my_im(ng), my_jm(ng), &
364 & my_lm(ng), my_mm(ng), &
365 & epsilon, &
366#ifdef FULL_GRID
367 & .true., &
368#else
369 & .false., &
370#endif
371 & my_domain(ng)% Xmin_rho(tile), &
372 & my_domain(ng)% Xmax_rho(tile), &
373 & my_domain(ng)% Ymin_rho(tile), &
374 & my_domain(ng)% Ymax_rho(tile))
375 END DO
376#ifdef DISTRIBUTE
377 rxmin(ng)=my_domain(ng)%Xmin_rho(myrank)
378 rxmax(ng)=my_domain(ng)%Xmax_rho(myrank)
379 rymin(ng)=my_domain(ng)%Ymin_rho(myrank)
380 rymax(ng)=my_domain(ng)%Ymax_rho(myrank)
381#else
382 rxmin(ng)=my_domain(ng)%Xmin_rho(0)
383 rxmax(ng)=my_domain(ng)%Xmax_rho(0)
384 rymin(ng)=my_domain(ng)%Ymin_rho(0)
385 rymax(ng)=my_domain(ng)%Ymax_rho(0)
386#endif
387 END DO
388!
389! Set U-points domain lower and upper bounds (integer).
390!
391 DO ng=1,ngrids
392 IF (ewperiodic(ng)) THEN
393 uoff=0
394 ELSE
395 uoff=1
396 END IF
397#ifdef DISTRIBUTE
398 CALL get_bounds (ng, myrank, u2dvar, 0, &
399 & my_im(ng), my_jm(ng), &
400 & my_lm(ng), my_mm(ng), &
401 & itile, jtile, &
402 & uilb(ng), uiub(ng), ujlb(ng), ujub(ng))
403# ifndef FULL_GRID
404 IF (itile.eq.0) THEN
405 uilb(ng)=uilb(ng)+uoff
406 END IF
407 IF (itile.eq.(ntilei(ng)-1)) THEN
408 uiub(ng)=uiub(ng)-1
409 END IF
410 IF (jtile.eq.0) THEN
411 ujlb(ng)=ujlb(ng)+1
412 END IF
413 IF (jtile.eq.(ntilej(ng)-1)) THEN
414 ujub(ng)=ujub(ng)-1
415 END IF
416# endif
417#else
418# ifdef FULL_GRID
419 uilb(ng)=1
420 uiub(ng)=my_lm(ng)+1
421 ujlb(ng)=0
422 ujub(ng)=my_mm(ng)+1
423# else
424 uilb(ng)=1+uoff
425 uiub(ng)=my_lm(ng)
426 ujlb(ng)=1
427 ujub(ng)=my_mm(ng)
428# endif
429#endif
430!
431! Minimum and maximum fractional coordinates for U-points.
432!
433 DO tile=0,ntilei(ng)*ntilej(ng)-1
434 CALL get_domain (ng, tile, u2dvar, 0, &
435 & my_im(ng), my_jm(ng), &
436 & my_lm(ng), my_mm(ng), &
437 & epsilon, &
438#ifdef FULL_GRID
439 & .true., &
440#else
441 & .false., &
442#endif
443 & my_domain(ng)% Xmin_u(tile), &
444 & my_domain(ng)% Xmax_u(tile), &
445 & my_domain(ng)% Ymin_u(tile), &
446 & my_domain(ng)% Ymax_u(tile))
447 END DO
448#ifdef DISTRIBUTE
449 uxmin(ng)=my_domain(ng)%Xmin_u(myrank)
450 uxmax(ng)=my_domain(ng)%Xmax_u(myrank)
451 uymin(ng)=my_domain(ng)%Ymin_u(myrank)
452 uymax(ng)=my_domain(ng)%Ymax_u(myrank)
453#else
454 uxmin(ng)=my_domain(ng)%Xmin_u(0)
455 uxmax(ng)=my_domain(ng)%Xmax_u(0)
456 uymin(ng)=my_domain(ng)%Ymin_u(0)
457 uymax(ng)=my_domain(ng)%Ymax_u(0)
458#endif
459 END DO
460!
461! Set V-points domain lower and upper bounds (integer).
462!
463 DO ng=1,ngrids
464 IF (nsperiodic(ng)) THEN
465 voff=0
466 ELSE
467 voff=1
468 END IF
469#ifdef DISTRIBUTE
470 CALL get_bounds (ng, myrank, v2dvar, 0, &
471 & my_im(ng), my_jm(ng), &
472 & my_lm(ng), my_mm(ng), &
473 & itile, jtile, &
474 & vilb(ng), viub(ng), vjlb(ng), vjub(ng))
475# ifndef FULL_GRID
476 IF (itile.eq.0) THEN
477 vilb(ng)=vilb(ng)+1
478 END IF
479 IF (itile.eq.(ntilei(ng)-1)) THEN
480 viub(ng)=viub(ng)-1
481 END IF
482 IF (jtile.eq.0) THEN
483 vjlb(ng)=vjlb(ng)+voff
484 END IF
485 IF (jtile.eq.(ntilej(ng)-1)) THEN
486 vjub(ng)=vjub(ng)-1
487 END IF
488# endif
489#else
490# ifdef FULL_GRID
491 vilb(ng)=0
492 viub(ng)=lm(ng)+1
493 vjlb(ng)=1
494 vjub(ng)=mm(ng)+1
495# else
496 vilb(ng)=1
497 viub(ng)=lm(ng)
498 vjlb(ng)=1+voff
499 vjub(ng)=mm(ng)
500# endif
501#endif
502!
503! Minimum and maximum fractional coordinates for V-points.
504!
505 DO tile=0,ntilei(ng)*ntilej(ng)-1
506 CALL get_domain (ng, tile, v2dvar, 0, &
507 & my_im(ng), my_jm(ng), &
508 & my_lm(ng), my_mm(ng), &
509 & epsilon, &
510#ifdef FULL_GRID
511 & .true., &
512#else
513 & .false., &
514#endif
515 & my_domain(ng)% Xmin_v(tile), &
516 & my_domain(ng)% Xmax_v(tile), &
517 & my_domain(ng)% Ymin_v(tile), &
518 & my_domain(ng)% Ymax_v(tile))
519 END DO
520#ifdef DISTRIBUTE
521 vxmin(ng)=my_domain(ng)%Xmin_v(myrank)
522 vxmax(ng)=my_domain(ng)%Xmax_v(myrank)
523 vymin(ng)=my_domain(ng)%Ymin_v(myrank)
524 vymax(ng)=my_domain(ng)%Ymax_v(myrank)
525#else
526 vxmin(ng)=my_domain(ng)%Xmin_v(0)
527 vxmax(ng)=my_domain(ng)%Xmax_v(0)
528 vymin(ng)=my_domain(ng)%Ymin_v(0)
529 vymax(ng)=my_domain(ng)%Ymax_v(0)
530#endif
531 END DO
532!
533 RETURN

References mod_scalars::ewperiodic, get_bounds_mod::get_bounds(), get_bounds_mod::get_domain(), mod_param::lm, mod_param::mm, mod_parallel::myrank, mod_param::ngrids, mod_scalars::nsperiodic, mod_param::ntilei, mod_param::ntilej, mod_param::r2dvar, mod_ncparam::rilb, mod_ncparam::riub, mod_ncparam::rjlb, mod_ncparam::rjub, mod_ncparam::rxmax, mod_ncparam::rxmin, mod_ncparam::rymax, mod_ncparam::rymin, mod_param::u2dvar, mod_ncparam::uilb, mod_ncparam::uiub, mod_ncparam::ujlb, mod_ncparam::ujub, mod_ncparam::uxmax, mod_ncparam::uxmin, mod_ncparam::uymax, mod_ncparam::uymin, mod_param::v2dvar, mod_ncparam::vilb, mod_ncparam::viub, mod_ncparam::vjlb, mod_ncparam::vjub, mod_ncparam::vxmax, mod_ncparam::vxmin, mod_ncparam::vymax, and mod_ncparam::vymin.

Referenced by inp_par_mod::inp_par().

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