MOM6
MOM_hor_index.F90
1 !> Defines the horizontal index type (hor_index_type) used for providing index ranges
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 use mom_domains, only : mom_domain_type, get_domain_extent, get_global_shape
7 use mom_error_handler, only : mom_error, mom_mesg, fatal
9 
10 implicit none ; private
11 
12 public :: hor_index_init, assignment(=)
13 public :: rotate_hor_index
14 
15 !> Container for horizontal index ranges for data, computational and global domains
16 type, public :: hor_index_type
17  integer :: isc !< The start i-index of cell centers within the computational domain
18  integer :: iec !< The end i-index of cell centers within the computational domain
19  integer :: jsc !< The start j-index of cell centers within the computational domain
20  integer :: jec !< The end j-index of cell centers within the computational domain
21 
22  integer :: isd !< The start i-index of cell centers within the data domain
23  integer :: ied !< The end i-index of cell centers within the data domain
24  integer :: jsd !< The start j-index of cell centers within the data domain
25  integer :: jed !< The end j-index of cell centers within the data domain
26 
27  integer :: isg !< The start i-index of cell centers within the global domain
28  integer :: ieg !< The end i-index of cell centers within the global domain
29  integer :: jsg !< The start j-index of cell centers within the global domain
30  integer :: jeg !< The end j-index of cell centers within the global domain
31 
32  integer :: iscb !< The start i-index of cell vertices within the computational domain
33  integer :: iecb !< The end i-index of cell vertices within the computational domain
34  integer :: jscb !< The start j-index of cell vertices within the computational domain
35  integer :: jecb !< The end j-index of cell vertices within the computational domain
36 
37  integer :: isdb !< The start i-index of cell vertices within the data domain
38  integer :: iedb !< The end i-index of cell vertices within the data domain
39  integer :: jsdb !< The start j-index of cell vertices within the data domain
40  integer :: jedb !< The end j-index of cell vertices within the data domain
41 
42  integer :: isgb !< The start i-index of cell vertices within the global domain
43  integer :: iegb !< The end i-index of cell vertices within the global domain
44  integer :: jsgb !< The start j-index of cell vertices within the global domain
45  integer :: jegb !< The end j-index of cell vertices within the global domain
46 
47  integer :: idg_offset !< The offset between the corresponding global and local i-indices.
48  integer :: jdg_offset !< The offset between the corresponding global and local j-indices.
49  logical :: symmetric !< True if symmetric memory is used.
50 
51  integer :: niglobal !< The global number of h-cells in the i-direction
52  integer :: njglobal !< The global number of h-cells in the j-direction
53 
54  integer :: turns !< Number of quarter-turn rotations from input to model
55 end type hor_index_type
56 
57 !> Copy the contents of one horizontal index type into another
58 interface assignment(=); module procedure HIT_assign ; end interface
59 
60 contains
61 
62 !> Sets various index values in a hor_index_type.
63 subroutine hor_index_init(Domain, HI, param_file, local_indexing, index_offset)
64  type(mom_domain_type), intent(in) :: domain !< The MOM domain from which to extract information.
65  type(hor_index_type), intent(inout) :: hi !< A horizontal index type to populate with data
66  type(param_file_type), intent(in) :: param_file !< Parameter file handle
67  logical, optional, intent(in) :: local_indexing !< If true, all tracer data domains start at 1
68  integer, optional, intent(in) :: index_offset !< A fixed additional offset to all indices
69 
70 ! This include declares and sets the variable "version".
71 #include "version_variable.h"
72 
73  ! get_domain_extent ensures that domains start at 1 for compatibility between
74  ! static and dynamically allocated arrays.
75  call get_domain_extent(domain, hi%isc, hi%iec, hi%jsc, hi%jec, &
76  hi%isd, hi%ied, hi%jsd, hi%jed, &
77  hi%isg, hi%ieg, hi%jsg, hi%jeg, &
78  hi%idg_offset, hi%jdg_offset, hi%symmetric, &
79  local_indexing=local_indexing)
80  call get_global_shape(domain, hi%niglobal, hi%njglobal)
81 
82  ! Read all relevant parameters and write them to the model log.
83  call log_version(param_file, "MOM_hor_index", version, &
84  "Sets the horizontal array index types.", all_default=.true.)
85 
86  hi%IscB = hi%isc ; hi%JscB = hi%jsc
87  hi%IsdB = hi%isd ; hi%JsdB = hi%jsd
88  hi%IsgB = hi%isg ; hi%JsgB = hi%jsg
89  if (hi%symmetric) then
90  hi%IscB = hi%isc-1 ; hi%JscB = hi%jsc-1
91  hi%IsdB = hi%isd-1 ; hi%JsdB = hi%jsd-1
92  hi%IsgB = hi%isg-1 ; hi%JsgB = hi%jsg-1
93  endif
94  hi%IecB = hi%iec ; hi%JecB = hi%jec
95  hi%IedB = hi%ied ; hi%JedB = hi%jed
96  hi%IegB = hi%ieg ; hi%JegB = hi%jeg
97 
98  hi%turns = 0
99 end subroutine hor_index_init
100 
101 !> HIT_assign copies one hor_index_type into another. It is accessed via an
102 !! assignment (=) operator.
103 subroutine hit_assign(HI1, HI2)
104  type(hor_index_type), intent(out) :: HI1 !< Horizontal index type to copy to
105  type(hor_index_type), intent(in) :: HI2 !< Horizontal index type to copy from
106  ! This subroutine copies all components of the horizontal array index type
107  ! variable on the RHS (HI2) to the variable on the LHS (HI1).
108 
109  hi1%isc = hi2%isc ; hi1%iec = hi2%iec ; hi1%jsc = hi2%jsc ; hi1%jec = hi2%jec
110  hi1%isd = hi2%isd ; hi1%ied = hi2%ied ; hi1%jsd = hi2%jsd ; hi1%jed = hi2%jed
111  hi1%isg = hi2%isg ; hi1%ieg = hi2%ieg ; hi1%jsg = hi2%jsg ; hi1%jeg = hi2%jeg
112 
113  hi1%IscB = hi2%IscB ; hi1%IecB = hi2%IecB ; hi1%JscB = hi2%JscB ; hi1%JecB = hi2%JecB
114  hi1%IsdB = hi2%IsdB ; hi1%IedB = hi2%IedB ; hi1%JsdB = hi2%JsdB ; hi1%JedB = hi2%JedB
115  hi1%IsgB = hi2%IsgB ; hi1%IegB = hi2%IegB ; hi1%JsgB = hi2%JsgB ; hi1%JegB = hi2%JegB
116 
117  hi1%niglobal = hi2%niglobal ; hi1%njglobal = hi2%njglobal
118  hi1%idg_offset = hi2%idg_offset ; hi1%jdg_offset = hi2%jdg_offset
119  hi1%symmetric = hi2%symmetric
120  hi1%turns = hi2%turns
121 end subroutine hit_assign
122 
123 !> Rotate the horizontal index ranges from the input to the output map.
124 subroutine rotate_hor_index(HI_in, turns, HI)
125  type(hor_index_type), intent(in) :: hi_in !< Unrotated horizontal indices
126  integer, intent(in) :: turns !< Number of quarter turns
127  type(hor_index_type), intent(inout) :: hi !< Rotated horizontal indices
128 
129  if (modulo(turns, 2) /= 0) then
130  hi%isc = hi_in%jsc
131  hi%iec = hi_in%jec
132  hi%jsc = hi_in%isc
133  hi%jec = hi_in%iec
134  hi%isd = hi_in%jsd
135  hi%ied = hi_in%jed
136  hi%jsd = hi_in%isd
137  hi%jed = hi_in%ied
138  hi%isg = hi_in%jsg
139  hi%ieg = hi_in%jeg
140  hi%jsg = hi_in%isg
141  hi%jeg = hi_in%ieg
142 
143  hi%IscB = hi_in%JscB
144  hi%IecB = hi_in%JecB
145  hi%JscB = hi_in%IscB
146  hi%JecB = hi_in%IecB
147  hi%IsdB = hi_in%JsdB
148  hi%IedB = hi_in%JedB
149  hi%JsdB = hi_in%IsdB
150  hi%JedB = hi_in%IedB
151  hi%IsgB = hi_in%JsgB
152  hi%IegB = hi_in%JegB
153  hi%JsgB = hi_in%IsgB
154  hi%JegB = hi_in%IegB
155 
156  hi%niglobal = hi_in%njglobal
157  hi%njglobal = hi_in%niglobal
158  hi%idg_offset = hi_in%jdg_offset
159  hi%jdg_offset = hi_in%idg_offset
160 
161  hi%symmetric = hi_in%symmetric
162  else
163  hi = hi_in
164  endif
165  hi%turns = hi_in%turns + turns
166 end subroutine rotate_hor_index
167 
168 !> \namespace mom_hor_index
169 !!
170 !! The hor_index_type provides the declarations and loop ranges for almost all data with horizontal extent.
171 !!
172 !! Declarations and loop ranges should always be coded with the symmetric memory model in mind.
173 !! The non-symmetric memory mode will then also work, albeit with a different (less efficient) communication pattern.
174 !!
175 !! Using the hor_index_type HI:
176 !! - declaration of h-point data is of the form `h(HI%%isd:HI%%ied,HI%%jsd:HI%%jed)`
177 !! - declaration of q-point data is of the form `q(HI%%IsdB:HI%%IedB,HI%%JsdB:HI%%JedB)`
178 !! - declaration of u-point data is of the form `u(HI%%IsdB:HI%%IedB,HI%%jsd:HI%%jed)`
179 !! - declaration of v-point data is of the form `v(HI%%isd:HI%%ied,HI%%JsdB:HI%%JedB)`.
180 !!
181 !! For more detail explanation of horizontal indexing see \ref Horizontal_Indexing.
182 
183 
184 end module mom_hor_index
A structure that can be parsed to read and document run-time parameters.
The MOM6 facility to parse input files for runtime parameters.
Defines the horizontal index type (hor_index_type) used for providing index ranges.
An overloaded interface to log the values of various types of parameters.
Container for horizontal index ranges for data, computational and global domains. ...
Describes the decomposed MOM domain and has routines for communications across PEs.
Definition: MOM_domains.F90:2
Routines for error handling and I/O management.
An overloaded interface to log version information about modules.
The MOM_domain_type contains information about the domain decompositoin.
An overloaded interface to read and log the values of various types of parameters.