MOM_hor_index.F90

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