MOM6
MOM_hor_index.F90
Go to the documentation of this file.
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 
7 use mom_error_handler, only : mom_error, mom_mesg, fatal
9 
10 implicit none ; private
11 
12 public :: hor_index_init, assignment(=)
13 
14 !> Container for horizontal index ranges for data, computational and global domains
15 type, public :: hor_index_type
16  integer :: isc !< The start i-index of cell centers within the computational domain
17  integer :: iec !< The end i-index of cell centers within the computational domain
18  integer :: jsc !< The start j-index of cell centers within the computational domain
19  integer :: jec !< The end j-index of cell centers within the computational domain
20 
21  integer :: isd !< The start i-index of cell centers within the data domain
22  integer :: ied !< The end i-index of cell centers within the data domain
23  integer :: jsd !< The start j-index of cell centers within the data domain
24  integer :: jed !< The end j-index of cell centers within the data domain
25 
26  integer :: isg !< The start i-index of cell centers within the global domain
27  integer :: ieg !< The end i-index of cell centers within the global domain
28  integer :: jsg !< The start j-index of cell centers within the global domain
29  integer :: jeg !< The end j-index of cell centers within the global domain
30 
31  integer :: iscb !< The start i-index of cell vertices within the computational domain
32  integer :: iecb !< The end i-index of cell vertices within the computational domain
33  integer :: jscb !< The start j-index of cell vertices within the computational domain
34  integer :: jecb !< The end j-index of cell vertices within the computational domain
35 
36  integer :: isdb !< The start i-index of cell vertices within the data domain
37  integer :: iedb !< The end i-index of cell vertices within the data domain
38  integer :: jsdb !< The start j-index of cell vertices within the data domain
39  integer :: jedb !< The end j-index of cell vertices within the data domain
40 
41  integer :: isgb !< The start i-index of cell vertices within the global domain
42  integer :: iegb !< The end i-index of cell vertices within the global domain
43  integer :: jsgb !< The start j-index of cell vertices within the global domain
44  integer :: jegb !< The end j-index of cell vertices within the global domain
45 
46  integer :: idg_offset !< The offset between the corresponding global and local i-indices.
47  integer :: jdg_offset !< The offset between the corresponding global and local j-indices.
48  logical :: symmetric !< True if symmetric memory is used.
49 end type hor_index_type
50 
51 interface assignment(=); module procedure hit_assign ; end interface
52 
53 contains
54 
55 !> Sets various index values in a hor_index_type.
56 subroutine hor_index_init(Domain, HI, param_file, local_indexing, index_offset)
57  type(mom_domain_type), intent(in) :: Domain !< The MOM domain from which to extract information.
58  type(hor_index_type), intent(inout) :: HI !< A horizontal index type to populate with data
59  type(param_file_type), intent(in) :: param_file !< Parameter file handle
60  logical, optional, intent(in) :: local_indexing !< If true, all tracer data domains start at 1
61  integer, optional, intent(in) :: index_offset !< A fixed additional offset to all indices
62 
63 ! This include declares and sets the variable "version".
64 #include "version_variable.h"
65 
66  ! get_domain_extent ensures that domains start at 1 for compatibility between
67  ! static and dynamically allocated arrays.
68  call get_domain_extent(domain, hi%isc, hi%iec, hi%jsc, hi%jec, &
69  hi%isd, hi%ied, hi%jsd, hi%jed, &
70  hi%isg, hi%ieg, hi%jsg, hi%jeg, &
71  hi%idg_offset, hi%jdg_offset, hi%symmetric, &
72  local_indexing=local_indexing)
73 
74  ! Read all relevant parameters and write them to the model log.
75  call log_version(param_file, "MOM_hor_index", version, &
76  "Sets the horizontal array index types.")
77 
78  hi%IscB = hi%isc ; hi%JscB = hi%jsc
79  hi%IsdB = hi%isd ; hi%JsdB = hi%jsd
80  hi%IsgB = hi%isg ; hi%JsgB = hi%jsg
81  if (hi%symmetric) then
82  hi%IscB = hi%isc-1 ; hi%JscB = hi%jsc-1
83  hi%IsdB = hi%isd-1 ; hi%JsdB = hi%jsd-1
84  hi%IsgB = hi%isg-1 ; hi%JsgB = hi%jsg-1
85  endif
86  hi%IecB = hi%iec ; hi%JecB = hi%jec
87  hi%IedB = hi%ied ; hi%JedB = hi%jed
88  hi%IegB = hi%ieg ; hi%JegB = hi%jeg
89 
90 end subroutine hor_index_init
91 
92 !> HIT_assign copies one hor_index_type into another. It is accessed via an
93 !! assignment (=) operator.
94 subroutine hit_assign(HI1, HI2)
95  type(hor_index_type), intent(out) :: HI1 !< Horizontal index type to copy to
96  type(hor_index_type), intent(in) :: HI2 !< Horizontal index type to copy from
97  ! This subroutine copies all components of the horizontal array index type
98  ! variable on the RHS (HI2) to the variable on the LHS (HI1).
99 
100  hi1%isc = hi2%isc ; hi1%iec = hi2%iec ; hi1%jsc = hi2%jsc ; hi1%jec = hi2%jec
101  hi1%isd = hi2%isd ; hi1%ied = hi2%ied ; hi1%jsd = hi2%jsd ; hi1%jed = hi2%jed
102  hi1%isg = hi2%isg ; hi1%ieg = hi2%ieg ; hi1%jsg = hi2%jsg ; hi1%jeg = hi2%jeg
103 
104  hi1%IscB = hi2%IscB ; hi1%IecB = hi2%IecB ; hi1%JscB = hi2%JscB ; hi1%JecB = hi2%JecB
105  hi1%IsdB = hi2%IsdB ; hi1%IedB = hi2%IedB ; hi1%JsdB = hi2%JsdB ; hi1%JedB = hi2%JedB
106  hi1%IsgB = hi2%IsgB ; hi1%IegB = hi2%IegB ; hi1%JsgB = hi2%JsgB ; hi1%JegB = hi2%JegB
107 
108  hi1%idg_offset = hi2%idg_offset ; hi1%jdg_offset = hi2%jdg_offset
109  hi1%symmetric = hi2%symmetric
110 
111 end subroutine hit_assign
112 
113 !> \namespace mom_hor_index
114 !!
115 !! The hor_index_type provides the decalarations and loop ranges for almost all data with horizontal extent.
116 !!
117 !! Declarations and loop ranges should always be coded with the symmetric memory model in mind.
118 !! The non-symmetric memory mode will then also work, albeit with a different (less efficient) communication pattern.
119 !!
120 !! Using the hor_index_type HI:
121 !! - declaration of h-point data is of the form `h(HI%%isd:HI%%ied,HI%%jsd:HI%%jed)`;
122 !! - declaration of q-point data is of the form `q(HI%%IsdB:HI%%IedB,HI%%JsdB:HI%%JedB)`;
123 !! - declaration of u-point data is of the form `u(HI%%IsdB:HI%%IedB,HI%%jsd:HI%%jed)`;
124 !! - declaration of v-point data is of the form `v(HI%%isd:HI%%ied,HI%%JsdB:HI%%JedB)`.
125 !!
126 !! For more detail explanation of horizontal indexing see \ref Horizontal_indexing.
127 
128 
129 end module mom_hor_index
subroutine, public hor_index_init(Domain, HI, param_file, local_indexing, index_offset)
Sets various index values in a hor_index_type.
Defines the horizontal index type (hor_index_type) used for providing index ranges.
Container for horizontal index ranges for data, computational and global domains. ...
subroutine hit_assign(HI1, HI2)
HIT_assign copies one hor_index_type into another. It is accessed via an assignment (=) operator...
subroutine, public mom_mesg(message, verb, all_print)
The MOM_domain_type contains information about the domain decompositoin.
subroutine, public get_domain_extent(Domain, isc, iec, jsc, jec, isd, ied, jsd, jed, isg, ieg, jsg, jeg, idg_offset, jdg_offset, symmetric, local_indexing, index_offset)
subroutine, public mom_error(level, message, all_print)