mom_hor_index module reference¶
Defines the horizontal index type (hor_index_type()
) used for providing index ranges. ) used for providing index ranges.
Data Types¶
Container for horizontal index ranges for data, computational and global domains. |
Functions/Subroutines¶
Sets various index values in a |
|
HIT_assign copies one |
|
Rotate the horizontal index ranges from the input to the output map. |
Detailed Description¶
The hor_index_type()
provides the declarations and loop ranges for almost all data with horizontal extent.provides the declarations and loop ranges for almost all data with horizontal extent.
Declarations and loop ranges should always be coded with the symmetric memory model in mind. The non-symmetric memory mode will then also work, albeit with a different (less efficient) communication pattern.
Using the hor_index_type()
HI:HI:
* declaration of h-point data is of the form h(HI%isd:HI%ied,HI%jsd:HI%jed)
declaration of q-point data is of the form
q(HI%IsdB:HI%IedB,HI%JsdB:HI%JedB)
declaration of u-point data is of the form
u(HI%IsdB:HI%IedB,HI%jsd:HI%jed)
declaration of v-point data is of the form
v(HI%isd:HI%ied,HI%JsdB:HI%JedB)
.
For more detail explanation of horizontal indexing see Horizontal indexing and memory.
Type Documentation¶
-
type
mom_hor_index/
hor_index_type
¶ Container for horizontal index ranges for data, computational and global domains.
- Type fields:
%
isc
[integer] :: The start i-index of cell centers within the computational domain.%
iec
[integer] :: The end i-index of cell centers within the computational domain.%
jsc
[integer] :: The start j-index of cell centers within the computational domain.%
jec
[integer] :: The end j-index of cell centers within the computational domain.%
isd
[integer] :: The start i-index of cell centers within the data domain.%
ied
[integer] :: The end i-index of cell centers within the data domain.%
jsd
[integer] :: The start j-index of cell centers within the data domain.%
jed
[integer] :: The end j-index of cell centers within the data domain.%
isg
[integer] :: The start i-index of cell centers within the global domain.%
ieg
[integer] :: The end i-index of cell centers within the global domain.%
jsg
[integer] :: The start j-index of cell centers within the global domain.%
jeg
[integer] :: The end j-index of cell centers within the global domain.%
iscb
[integer] :: The start i-index of cell vertices within the computational domain.%
iecb
[integer] :: The end i-index of cell vertices within the computational domain.%
jscb
[integer] :: The start j-index of cell vertices within the computational domain.%
jecb
[integer] :: The end j-index of cell vertices within the computational domain.%
isdb
[integer] :: The start i-index of cell vertices within the data domain.%
iedb
[integer] :: The end i-index of cell vertices within the data domain.%
jsdb
[integer] :: The start j-index of cell vertices within the data domain.%
jedb
[integer] :: The end j-index of cell vertices within the data domain.%
isgb
[integer] :: The start i-index of cell vertices within the global domain.%
iegb
[integer] :: The end i-index of cell vertices within the global domain.%
jsgb
[integer] :: The start j-index of cell vertices within the global domain.%
jegb
[integer] :: The end j-index of cell vertices within the global domain.%
idg_offset
[integer] :: The offset between the corresponding global and local i-indices.%
jdg_offset
[integer] :: The offset between the corresponding global and local j-indices.%
symmetric
[logical] :: True if symmetric memory is used.%
niglobal
[integer] :: The global number of h-cells in the i-direction.%
njglobal
[integer] :: The global number of h-cells in the j-direction.%
turns
[integer] :: Number of quarter-turn rotations from input to model.
Function/Subroutine Documentation¶
-
subroutine
mom_hor_index/
hor_index_init
(Domain, HI, param_file, local_indexing, index_offset)¶ Sets various index values in a
hor_index_type()
. .- Parameters:
domain :: [in] The MOM domain from which to extract information.
hi :: [inout] A horizontal index type to populate with data
param_file :: [in] Parameter file handle
local_indexing :: [in] If true, all tracer data domains start at 1
index_offset :: [in] A fixed additional offset to all indices
- Called from:
mom_oda_driver_mod::init_oda
mom_grid::mom_grid_init
mom_io_file::open_file_nc
-
subroutine
mom_hor_index/
hit_assign
(HI1, HI2)¶ HIT_assign copies one
hor_index_type()
into another. It is accessed via an assignment (=) operator. into another. It is accessed via an assignment (=) operator.- Parameters:
hi1 :: [out] Horizontal index type to copy to
hi2 :: [in] Horizontal index type to copy from
-
subroutine
mom_hor_index/
rotate_hor_index
(HI_in, turns, HI)¶ Rotate the horizontal index ranges from the input to the output map.
- Parameters:
hi_in :: [in] Unrotated horizontal indices
turns :: [in] Number of quarter turns
hi :: [inout] Rotated horizontal indices
- Called from:
mom_checksums::bchksum::chksum_b_2d
mom_checksums::bchksum::chksum_b_3d
mom_checksums::hchksum::chksum_h_2d
mom_checksums::hchksum::chksum_h_3d
mom_checksums::bchksum_pair::chksum_pair_b_2d
mom_checksums::bchksum_pair::chksum_pair_b_3d
mom_checksums::hchksum_pair::chksum_pair_h_2d
mom_checksums::hchksum_pair::chksum_pair_h_3d
mom_checksums::uchksum::chksum_u_2d
mom_checksums::uchksum::chksum_u_3d
mom_checksums::uvchksum::chksum_uv_2d
mom_checksums::uvchksum::chksum_uv_3d
mom_checksums::vchksum::chksum_v_2d
mom_checksums::vchksum::chksum_v_3d
mom::initialize_mom