FMS  2025.04
Flexible Modeling System
fms_diag_bbox.F90
1 !***********************************************************************
2 !* Apache License 2.0
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* Licensed under the Apache License, Version 2.0 (the "License");
7 !* you may not use this file except in compliance with the License.
8 !* You may obtain a copy of the License at
9 !*
10 !* http://www.apache.org/licenses/LICENSE-2.0
11 !*
12 !* FMS is distributed in the hope that it will be useful, but WITHOUT
13 !* WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied;
14 !* without even the implied warranty of MERCHANTABILITY or FITNESS FOR A
15 !* PARTICULAR PURPOSE. See the License for the specific language
16 !* governing permissions and limitations under the License.
17 !***********************************************************************
18 
19 !> @defgroup fms_diag_bbox_mod fms_diag_bbox_mod
20 !> @ingroup diag_manager
21 !> @brief fms_diag_bbox_mod defines classes encapsulating bounding boxes
22 !! and interval bounds.
23 !!
24 !> @author Miguel Zuniga
25 !!
26 !> @file
27 !> @brief File for @ref fms_diag_bbox_mod
28 !> @addtogroup fms_diag_bbox_mod
29 !> @{
30 MODULE fms_diag_bbox_mod
31 
32  USE fms_mod, ONLY: error_mesg, fatal, fms_error_handler, string
33 
34  implicit none
35 
36 !> @brief Data structure holding a 3D bounding box. It is commonlyused to
37 !! represent the interval bounds or limits of a 3D sub-array such as the
38 !! array index bounds of the spatial component a diag_manager field output
39 !! buffer array.
40  TYPE, public :: fmsdiagibounds_type
41  INTEGER :: imin !< Lower i bound.
42  INTEGER :: imax !< Upper i bound.
43  INTEGER :: jmin !< Lower j bound.
44  INTEGER :: jmax !< Upper j bound.
45  INTEGER :: kmin !< Lower k bound.
46  INTEGER :: kmax !< Upper k bound.
47  logical :: has_halos !< .True. if the buffer has halos
48  integer :: nhalo_i !< Number of halos in i
49  integer :: nhalo_j !< Number of halos in j
50  contains
51  procedure :: reset => reset_bounds
52  procedure :: reset_bounds_from_array_4d
53  procedure :: reset_bounds_from_array_5d
54  procedure :: update_bounds
55  procedure :: set_bounds
56  procedure :: rebase_input
57  procedure :: rebase_output
58  procedure :: get_imin
59  procedure :: get_imax
60  procedure :: get_jmin
61  procedure :: get_jmax
62  procedure :: get_kmin
63  procedure :: get_kmax
64  procedure :: update_index
65  END TYPE fmsdiagibounds_type
66 
67  !> @brief Data structure holding starting and ending indices in the I, J, and
68  !! K dimensions. It also has extra members related to halo sizes and updated indices
69  !! in I and J dimensions.
70  type, public :: fmsdiagboundshalos_type
71  private
72  type(fmsdiagibounds_type) :: bounds3d !< Holds starting and ending indices of
73  !! the I, J, and K dimensions
74  integer :: hi !< Halo size in the I dimension
75  integer :: hj !< Halo size in the J dimension
76  integer :: fis !< Updated starting index in the I dimension
77  integer :: fie !< Updated ending index in the I dimension
78  integer :: fjs !< Updated starting index in the J dimension
79  integer :: fje !< Updated ending index in the J dimension
80  contains
81  procedure :: get_hi
82  procedure :: get_hj
83  procedure :: get_fis
84  procedure :: get_fie
85  procedure :: get_fjs
86  procedure :: get_fje
88 
90 
91  integer, parameter :: xdimension = 1 !< Parameter defining the x dimension
92  integer, parameter :: ydimension = 2 !< Parameter defining the y dimension
93  integer, parameter :: zdimension = 3 !< Parameter defininf the z dimension
94 
95 CONTAINS
96 
97 !> @brief The PEs grid points are divided further into "blocks". This function determines if a block
98 ! has data for a given subregion and dimension
99 !! @return .true. if the a subergion is inside a block
100 logical pure function determine_if_block_is_in_region(subregion_start, subregion_end, bounds, dim)
101  integer, intent(in) :: subregion_start !< Begining of the subregion
102  integer, intent(in) :: subregion_end !< Ending of the subregion
103  type(fmsdiagibounds_type), intent(in) :: bounds !< Starting and ending of the subregion
104  integer, intent(in) :: dim !< Dimension to check
105 
106  integer :: block_start !< Begining index of the block
107  integer :: block_end !< Ending index of the block
108 
110  select case (dim)
111  case (xdimension)
112  block_start = bounds%imin
113  block_end = bounds%imax
114  case (ydimension)
115  block_start = bounds%jmin
116  block_end = bounds%jmax
117  case (zdimension)
118  block_start = bounds%kmin
119  block_end = bounds%kmax
120  end select
121 
122  if (block_start < subregion_start .and. block_end < subregion_start) then
124  return
125  endif
126 
127  if (block_start > subregion_end .and. block_end > subregion_end) then
129  return
130  endif
131 
134 
135  !> @brief Gets imin of fmsDiagIbounds_type
136  !! @return copy of integer member imin
137  pure integer function get_imin (this) result(rslt)
138  class(fmsdiagibounds_type), intent(in) :: this !< The !< ibounds instance
139  rslt = this%imin
140  end function get_imin
141 
142  !> @brief Gets imax of fmsDiagIbounds_type
143  !! @return copy of integer member imax
144  pure integer function get_imax (this) result(rslt)
145  class(fmsdiagibounds_type), intent(in) :: this !< The !< ibounds instance
146  rslt = this%imax
147  end function get_imax
148 
149  !> @brief Gets jmin of fmsDiagIbounds_type
150  !! @return copy of integer member jmin
151  pure integer function get_jmin (this) result(rslt)
152  class(fmsdiagibounds_type), intent(in) :: this !< The !< ibounds instance
153  rslt = this%jmin
154  end function get_jmin
155 
156  !> @brief Gets jmax of fmsDiagIbounds_type
157  !! @return copy of integer member jmax
158  pure integer function get_jmax (this) result(rslt)
159  class(fmsdiagibounds_type), intent(in) :: this !< The !< ibounds instance
160  rslt = this%jmax
161  end function get_jmax
162 
163 
164  !> @brief Gets kmin of fmsDiagIbounds_type
165  !! @return copy of integer member kmin
166  pure integer function get_kmin (this) result(rslt)
167  class(fmsdiagibounds_type), intent(in) :: this !< The !< ibounds instance
168  rslt = this%kmin
169  end function get_kmin
170 
171  !> @brief Gets kmax of fmsDiagIbounds_type
172  !! @return copy of integer member kmax
173  pure integer function get_kmax (this) result(rslt)
174  class(fmsdiagibounds_type), intent(in) :: this !< The !< ibounds instance
175  rslt = this%kmax
176  end function get_kmax
177 
178  !> @brief Updates the starting and ending index of a given dimension
179  subroutine update_index(this, starting_index, ending_index, dim, ignore_halos)
180  class(fmsdiagibounds_type), intent(inout) :: this !< The bounding box to update
181  integer, intent(in) :: starting_index !< Starting index to update to
182  integer, intent(in) :: ending_index !< Ending index to update to
183  integer, intent(in) :: dim !< Dimension to update
184  logical, intent(in) :: ignore_halos !< If .true. halos will be ignored
185  !! i.e output buffers can ignore halos as
186  !! they do not get updates. The indices of the
187  !! Input buffers need to add the number of halos
188  !! so math is done only on the compute domain
189 
190  integer :: nhalox !< Number of halos in x
191  integer :: nhaloy !< Number of halos in y
192 
193  if (ignore_halos) then
194  nhalox = 0
195  nhaloy = 0
196  else
197  nhalox= this%nhalo_I
198  nhaloy= this%nhalo_J
199  endif
200  select case(dim)
201  case (xdimension)
202  this%imin = starting_index + nhalox
203  this%imax = ending_index + nhalox
204  case (ydimension)
205  this%jmin = starting_index + nhaloy
206  this%jmax = ending_index + nhaloy
207  case (zdimension)
208  this%kmin = starting_index
209  this%kmax = ending_index
210  end select
211  end subroutine
212 
213  !> @brief Gets the halo size of fmsDiagBoundsHalos_type in the I dimension
214  !! @return copy of integer member hi
215  pure integer function get_hi (this) result(rslt)
216  class(fmsdiagboundshalos_type), intent(in) :: this !< Calling object
217  rslt = this%hi
218  end function get_hi
219 
220  !> @brief Gets the halo size of fmsDiagBoundsHalos_type in the J dimension
221  !! @return copy of integer member hj
222  pure integer function get_hj (this) result(rslt)
223  class(fmsdiagboundshalos_type), intent(in) :: this !< Calling object
224  rslt = this%hj
225  end function get_hj
226 
227  !> @brief Gets the updated index `fis' of fmsDiagBoundsHalos_type in the I dimension
228  !! @return copy of integer member `fis'
229  pure integer function get_fis (this) result(rslt)
230  class(fmsdiagboundshalos_type), intent(in) :: this !< Calling object
231  rslt = this%fis
232  end function get_fis
233 
234  !> @brief Gets the updated index `fie' of fmsDiagBoundsHalos_type in the I dimension
235  !! @return copy of integer member `fie'
236  pure integer function get_fie (this) result(rslt)
237  class(fmsdiagboundshalos_type), intent(in) :: this !< Calling object
238  rslt = this%fie
239  end function get_fie
240 
241  !> @brief Gets the updated index `fjs' of fmsDiagBoundsHalos_type in the I dimension
242  !! @return copy of integer member `fjs'
243  pure integer function get_fjs (this) result(rslt)
244  class(fmsdiagboundshalos_type), intent(in) :: this !< Calling object
245  rslt = this%fjs
246  end function get_fjs
247 
248  !> @brief Gets the updated index `fje' of fmsDiagBoundsHalos_type in the I dimension
249  !! @return copy of integer member `fje'
250  pure integer function get_fje (this) result(rslt)
251  class(fmsdiagboundshalos_type), intent(in) :: this !< Calling object
252  rslt = this%fje
253  end function get_fje
254 
255  !> @brief Reset the instance bounding lower and upper bounds to lower_val and upper_val, respectively.
256  SUBROUTINE reset_bounds (this, lower_val, upper_val)
257  class(fmsdiagibounds_type), target, intent(inout) :: this !< ibounds instance
258  integer, intent(in) :: lower_val !< value for the lower bounds in each dimension
259  integer, intent(in) :: upper_val !< value for the upper bounds in each dimension
260  this%imin = lower_val
261  this%jmin = lower_val
262  this%kmin = lower_val
263  this%imax = upper_val
264  this%jmax = upper_val
265  this%kmax = upper_val
266  END SUBROUTINE reset_bounds
267 
268  !> @brief Update the first three (normally x, y, and z) min and
269  !! max boundaries (array indices) of the instance bounding box
270  !! the six specified bounds values.
271  SUBROUTINE update_bounds(this, lower_i, upper_i, lower_j, upper_j, lower_k, upper_k)
272  class(fmsdiagibounds_type), intent(inout) :: this !<The bounding box of the output field buffer inindex space.
273  INTEGER, INTENT(in) :: lower_i !< Lower i bound.
274  INTEGER, INTENT(in) :: upper_i !< Upper i bound.
275  INTEGER, INTENT(in) :: lower_j !< Lower j bound.
276  INTEGER, INTENT(in) :: upper_j !< Upper j bound.
277  INTEGER, INTENT(in) :: lower_k !< Lower k bound.
278  INTEGER, INTENT(in) :: upper_k !< Upper k bound.
279  this%imin = min(this%imin, lower_i)
280  this%imax = max(this%imax, upper_i)
281  this%jmin = min(this%jmin, lower_j)
282  this%jmax = max(this%jmax, upper_j)
283  this%kmin = min(this%kmin, lower_k)
284  this%kmax = max(this%kmax, upper_k)
285  END SUBROUTINE update_bounds
286 
287  !> @brief Sets the bounds of a bounding region
288  !! @return empty string if sucessful or error message if unsucessful
289  function set_bounds(this, field_data, lower_i, upper_i, lower_j, upper_j, lower_k, upper_k, has_halos) &
290  result(error_msg)
291  class(fmsdiagibounds_type), intent(inout) :: this !< The bounding box of the field
292  class(*), intent(in) :: field_data(:,:,:,:) !< Field data
293  INTEGER, INTENT(in) :: lower_i !< Lower i bound.
294  INTEGER, INTENT(in) :: upper_i !< Upper i bound.
295  INTEGER, INTENT(in) :: lower_j !< Lower j bound.
296  INTEGER, INTENT(in) :: upper_j !< Upper j bound.
297  INTEGER, INTENT(in) :: lower_k !< Lower k bound.
298  INTEGER, INTENT(in) :: upper_k !< Upper k bound.
299  LOGICAL, INTENT(in) :: has_halos !< .true. if the field has halos
300 
301  character(len=150) :: error_msg !< Error message to output
302 
303  integer :: nhalos_2 !< 2 times the number of halo points
304  integer :: nhalox !< Number of halos in x
305  integer :: nhaloy !< Number of halos in y
306 
307  error_msg = ""
308  this%kmin = lower_k
309  this%kmax = upper_k
310  this%has_halos = has_halos
311  this%nhalo_I = 0
312  this%nhalo_J = 0
313  if (has_halos) then
314  !upper_i-lower_i+1 is the size of the compute domain
315  !ubound(field_data,1) is the size of the data domain
316  nhalos_2 = ubound(field_data,1)-(upper_i-lower_i+1)
317  if (mod(nhalos_2, 2) .ne. 0) then
318  error_msg = "There are non-symmetric halos in the first dimension"
319  return
320  endif
321  nhalox = nhalos_2/2
322  this%nhalo_I = nhalox
323 
324  nhalos_2 = ubound(field_data,2)-(upper_j-lower_j + 1)
325  if (mod(nhalos_2, 2) .ne. 0) then
326  error_msg = "There are non-symmetric halos in the second dimension"
327  return
328  endif
329  nhaloy = nhalos_2/2
330  this%nhalo_J = nhaloy
331 
332  this%imin = 1 + nhalox
333  this%imax = ubound(field_data,1) - nhalox
334  this%jmin = 1 + nhaloy
335  this%jmax = ubound(field_data,2) - nhaloy
336  else
337  this%imin = lower_i
338  this%imax = upper_i
339  this%jmin = lower_j
340  this%jmax = upper_j
341  endif
342 
343  end function set_bounds
344  !> @brief Reset the instance bounding box with the bounds determined from the
345  !! first three dimensions of the 5D "array" argument
346  SUBROUTINE reset_bounds_from_array_4d(this, array)
347  class(fmsdiagibounds_type), INTENT(inout) :: this !< The instance of the bounding box.
348  class(*), INTENT( in), DIMENSION(:,:,:,:) :: array !< The 4D input array.
349  this%imin = lbound(array,1)
350  this%imax = ubound(array,1)
351  this%jmin = lbound(array,2)
352  this%jmax = ubound(array,2)
353  this%kmin = lbound(array,3)
354  this%kmax = ubound(array,3)
355 
356  this%has_halos = .false.
357  this%nhalo_I = 0
358  this%nhalo_J = 0
359  END SUBROUTINE reset_bounds_from_array_4d
360 
361  !> @brief Reset the instance bounding box with the bounds determined from the
362  !! first three dimensions of the 5D "array" argument
363  SUBROUTINE reset_bounds_from_array_5d(this, array)
364  class(fmsdiagibounds_type), INTENT(inout) :: this !< The instance of the bounding box.
365  CLASS(*), INTENT( in), DIMENSION(:,:,:,:,:) :: array !< The 5D input array.
366  this%imin = lbound(array,1)
367  this%imax = ubound(array,1)
368  this%jmin = lbound(array,2)
369  this%jmax = ubound(array,2)
370  this%kmin = lbound(array,3)
371  this%kmax = ubound(array,3)
372  END SUBROUTINE reset_bounds_from_array_5d
373 
374  !> @brief Updates indices based on presence/absence of input indices is, js, ks, ie, je, and ke.
375  ! Computes halo sizes in the I and J dimensions.
376  ! This routine is intended to be used in diag manager.
377  !> @return .false. if there is no error else .true.
378  function recondition_indices(indices, field, is_in, js_in, ks_in, &
379  ie_in, je_in, ke_in, err_msg) result(ierr)
380  type(fmsdiagboundshalos_type), intent(inout) :: indices !< Stores indices in order:
381  !! (/is, js, ks, ie, je, ke, hi, fis, fie, hj, fjs, fje/)
382  class(*), intent(in) :: field(:,:,:,:) !< Dummy variable; only the sizes of the first 3 dimensions are used
383  integer, intent(in), optional :: is_in, js_in, ks_in, ie_in, je_in, ke_in !< User input indices
384  character(len=*), intent(out), optional :: err_msg !< Error message to pass back to caller
385  logical :: ierr !< Error flag
386 
387  integer :: is, js, ks, ie, je, ke !< Local indices to update
388  integer :: hi !< halo size in the I dimension
389  integer :: hj !< halo size in the J dimension
390  integer :: twohi, twohj !< Temporary storages
391  integer :: fis, fie, fjs, fje !< ! Updated starting and ending indices in the I and J dimensions
392  integer :: n1, n2, n3 !< Sizes of the first 3 dimenstions indicies of the data
393 
394  ierr = .false.
395  if (present(err_msg)) err_msg = ''
396 
397  ! If is, js, or ks not present default them to 1
398  is = 1
399  js = 1
400  ks = 1
401 
402  IF ( PRESENT(is_in) ) is = is_in
403  IF ( PRESENT(js_in) ) js = js_in
404  IF ( PRESENT(ks_in) ) ks = ks_in
405 
406  n1 = SIZE(field, 1)
407  n2 = SIZE(field, 2)
408  n3 = SIZE(field, 3)
409 
410  ie = is + n1 - 1
411  je = js + n2 - 1
412  ke = ks + n3 - 1
413 
414  IF ( PRESENT(ie_in) ) ie = ie_in
415  IF ( PRESENT(je_in) ) je = je_in
416  IF ( PRESENT(ke_in) ) ke = ke_in
417 
418  twohi = n1 - (ie - is + 1)
419  IF ( mod(twohi, 2) /= 0 ) THEN
420  ierr = fms_error_handler('diag_util_mod:recondition_indices', &
421  'non-symmetric halos in first dimension', err_msg)
422  IF (ierr) RETURN
423  END IF
424 
425  twohj = n2 - (je - js + 1)
426  IF ( mod(twohj, 2) /= 0 ) THEN
427  ierr = fms_error_handler('diag_util_mod:recondition_indices', &
428  'non-symmetric halos in second dimension', err_msg)
429  IF (ierr) RETURN
430  END IF
431 
432  hi = twohi/2
433  hj = twohj/2
434 
435  ! The next line is necessary to ensure that is, ie, js, ie are relative to field(1:,1:)
436  ! But this works only when there is no windowing.
437  IF ( PRESENT(ie_in) .AND. PRESENT(je_in) ) THEN
438  is = 1 + hi
439  ie = n1 - hi
440  js = 1 + hj
441  je = n2 - hj
442  END IF
443 
444  ! Used for field, mask and rmask bounds
445  fis = 1 + hi
446  fie = n1 - hi
447  fjs = 1 + hj
448  fje = n2 - hj
449 
450  ! Update indices
451  indices%bounds3D%imin = is
452  indices%bounds3D%imax = ie
453  indices%bounds3D%jmin = js
454  indices%bounds3D%jmax = je
455  indices%bounds3D%kmin = ks
456  indices%bounds3D%kmax = ke
457  indices%hi = hi
458  indices%hj = hj
459  indices%fis = fis
460  indices%fie = fie
461  indices%fjs = fjs
462  indices%fje = fje
463  end function recondition_indices
464 
465  !> @brief Rebase the ouput bounds for a given dimension based on the starting and ending indices of
466  !! a subregion. This is for when blocking is used.
467  subroutine rebase_output(bounds_out, starting, ending, dim)
468  class(fmsdiagibounds_type), INTENT(inout) :: bounds_out !< Bounds to rebase
469  integer, intent(in) :: starting !< Starting index of the dimension
470  integer, intent(in) :: ending !< Ending index of the dimension
471  integer, intent(in) :: dim !< Dimension to update
472 
473  !> The starting index is going to be either "starting" if only a section of the
474  !! block is in the subregion or bounds_out%[]min if the whole section of the block is in the
475  !! subregion. The -starting+1 s needed so that indices start as 1 since the output buffer has
476  !! indices 1:size of a subregion
477 
478  !> The ending index is going to be either bounds_out%[]max if the whole section of the block
479  !! is in the subregion or bounds_out%[]min + size of the subregion if only a section of the
480  !! block is in the susbregion
481  select case (dim)
482  case (xdimension)
483  bounds_out%imin = max(starting, bounds_out%imin)-starting+1
484  bounds_out%imax = min(bounds_out%imax, bounds_out%imin + ending-starting)
485  case (ydimension)
486  bounds_out%jmin = max(starting, bounds_out%jmin)-starting+1
487  bounds_out%jmax = min(bounds_out%jmax, bounds_out%jmin + ending-starting)
488  case (zdimension)
489  bounds_out%kmin =max(starting, bounds_out%kmin)-starting+1
490  bounds_out%kmax = min(bounds_out%kmax, bounds_out%kmin + ending-starting)
491  end select
492  end subroutine
493 
494  !> @brief Rebase the input bounds for a given dimension based on the starting and ending indices
495  !! of a subregion. This is for when blocking is used
496  subroutine rebase_input(bounds_in, bounds, starting, ending, dim)
497  class(fmsdiagibounds_type), INTENT(inout) :: bounds_in !< Bounds to rebase
498  class(fmsdiagibounds_type), INTENT(in) :: bounds !< Original indices (i.e is_in, ie_in,
499  !! passed into diag_manager)
500  integer, intent(in) :: starting !< Starting index of the dimension
501  integer, intent(in) :: ending !< Ending index of the dimension
502  integer, intent(in) :: dim !< Dimension to update
503 
504  !> The starting index is going to be either "starting" if only a section of the
505  !! block is in the subregion or starting-bounds%imin+1 if the whole section of the block is in the
506  !! subregion.
507 
508  !> The ending index is going to be either bounds_out%[]max if the whole section of the block
509  !! is in the subregion or bounds%[]min + size of the subregion if only a section of the
510  !! block is in the susbregion
511  select case (dim)
512  case (xdimension)
513  bounds_in%imin = min(abs(starting-bounds%imin+1), starting)
514  bounds_in%imax = min(bounds_in%imax, (bounds_in%imin + ending-starting))
515  case (ydimension)
516  bounds_in%jmin = min(abs(starting-bounds%jmin+1), starting)
517  bounds_in%jmax = min(bounds_in%jmax, (bounds_in%jmin + ending-starting))
518  case (zdimension)
519  bounds_in%kmin = min(abs(starting-bounds%kmin+1), starting)
520  bounds_in%kmax = min(bounds_in%kmax, (bounds_in%kmin + ending-starting))
521  end select
522  end subroutine
523 
524  END MODULE fms_diag_bbox_mod
525  !> @}
526  ! close documentation grouping
subroutine rebase_input(bounds_in, bounds, starting, ending, dim)
Rebase the input bounds for a given dimension based on the starting and ending indices of a subregion...
logical function, public recondition_indices(indices, field, is_in, js_in, ks_in, ie_in, je_in, ke_in, err_msg)
Updates indices based on presence/absence of input indices is, js, ks, ie, je, and ke.
subroutine update_bounds(this, lower_i, upper_i, lower_j, upper_j, lower_k, upper_k)
Update the first three (normally x, y, and z) min and max boundaries (array indices) of the instance ...
integer, parameter ydimension
Parameter defining the y dimension.
pure integer function get_imin(this)
Gets imin of fmsDiagIbounds_type.
subroutine reset_bounds_from_array_4d(this, array)
Reset the instance bounding box with the bounds determined from the first three dimensions of the 5D ...
pure integer function get_fjs(this)
Gets the updated index ‘fjs’ of fmsDiagBoundsHalos_type in the I dimension.
integer, parameter xdimension
Parameter defining the x dimension.
subroutine rebase_output(bounds_out, starting, ending, dim)
Rebase the ouput bounds for a given dimension based on the starting and ending indices of a subregion...
pure integer function get_hi(this)
Gets the halo size of fmsDiagBoundsHalos_type in the I dimension.
pure integer function get_imax(this)
Gets imax of fmsDiagIbounds_type.
pure integer function get_fje(this)
Gets the updated index ‘fje’ of fmsDiagBoundsHalos_type in the I dimension.
pure integer function get_kmax(this)
Gets kmax of fmsDiagIbounds_type.
pure integer function get_fie(this)
Gets the updated index ‘fie’ of fmsDiagBoundsHalos_type in the I dimension.
logical pure function, public determine_if_block_is_in_region(subregion_start, subregion_end, bounds, dim)
The PEs grid points are divided further into "blocks". This function determines if a block.
pure integer function get_hj(this)
Gets the halo size of fmsDiagBoundsHalos_type in the J dimension.
subroutine reset_bounds(this, lower_val, upper_val)
Reset the instance bounding lower and upper bounds to lower_val and upper_val, respectively.
subroutine update_index(this, starting_index, ending_index, dim, ignore_halos)
Updates the starting and ending index of a given dimension.
subroutine reset_bounds_from_array_5d(this, array)
Reset the instance bounding box with the bounds determined from the first three dimensions of the 5D ...
pure integer function get_fis(this)
Gets the updated index ‘fis’ of fmsDiagBoundsHalos_type in the I dimension.
pure integer function get_jmin(this)
Gets jmin of fmsDiagIbounds_type.
pure integer function get_jmax(this)
Gets jmax of fmsDiagIbounds_type.
pure integer function get_kmin(this)
Gets kmin of fmsDiagIbounds_type.
integer, parameter zdimension
Parameter defininf the z dimension.
character(len=150) function set_bounds(this, field_data, lower_i, upper_i, lower_j, upper_j, lower_k, upper_k, has_halos)
Sets the bounds of a bounding region.
Data structure holding starting and ending indices in the I, J, and K dimensions. It also has extra m...
logical function, public fms_error_handler(routine, message, err_msg)
Facilitates the control of fatal error conditions.
Definition: fms.F90:468
subroutine, public error_mesg(routine, message, level)
Print notes, warnings and error messages; terminates program for warning and error messages....
Definition: fms.F90:441
Data structure holding a 3D bounding box. It is commonlyused to represent the interval bounds or limi...