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
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 !> Copy the contents of one horizontal index type into another
52 interface assignment(=); module procedure HIT_assign ; end interface
53 
54 contains
55 
56 !> Sets various index values in a hor_index_type.
57 subroutine hor_index_init(Domain, HI, param_file, local_indexing, index_offset)
58  type(mom_domain_type), intent(in) :: domain !< The MOM domain from which to extract information.
59  type(hor_index_type), intent(inout) :: hi !< A horizontal index type to populate with data
60  type(param_file_type), intent(in) :: param_file !< Parameter file handle
61  logical, optional, intent(in) :: local_indexing !< If true, all tracer data domains start at 1
62  integer, optional, intent(in) :: index_offset !< A fixed additional offset to all indices
63 
64 ! This include declares and sets the variable "version".
65 #include "version_variable.h"
66 
67  ! get_domain_extent ensures that domains start at 1 for compatibility between
68  ! static and dynamically allocated arrays.
69  call get_domain_extent(domain, hi%isc, hi%iec, hi%jsc, hi%jec, &
70  hi%isd, hi%ied, hi%jsd, hi%jed, &
71  hi%isg, hi%ieg, hi%jsg, hi%jeg, &
72  hi%idg_offset, hi%jdg_offset, hi%symmetric, &
73  local_indexing=local_indexing)
74 
75  ! Read all relevant parameters and write them to the model log.
76  call log_version(param_file, "MOM_hor_index", version, &
77  "Sets the horizontal array index types.")
78 
79  hi%IscB = hi%isc ; hi%JscB = hi%jsc
80  hi%IsdB = hi%isd ; hi%JsdB = hi%jsd
81  hi%IsgB = hi%isg ; hi%JsgB = hi%jsg
82  if (hi%symmetric) then
83  hi%IscB = hi%isc-1 ; hi%JscB = hi%jsc-1
84  hi%IsdB = hi%isd-1 ; hi%JsdB = hi%jsd-1
85  hi%IsgB = hi%isg-1 ; hi%JsgB = hi%jsg-1
86  endif
87  hi%IecB = hi%iec ; hi%JecB = hi%jec
88  hi%IedB = hi%ied ; hi%JedB = hi%jed
89  hi%IegB = hi%ieg ; hi%JegB = hi%jeg
90 
91 end subroutine hor_index_init
92 
93 !> HIT_assign copies one hor_index_type into another. It is accessed via an
94 !! assignment (=) operator.
95 subroutine hit_assign(HI1, HI2)
96  type(hor_index_type), intent(out) :: HI1 !< Horizontal index type to copy to
97  type(hor_index_type), intent(in) :: HI2 !< Horizontal index type to copy from
98  ! This subroutine copies all components of the horizontal array index type
99  ! variable on the RHS (HI2) to the variable on the LHS (HI1).
100 
101  hi1%isc = hi2%isc ; hi1%iec = hi2%iec ; hi1%jsc = hi2%jsc ; hi1%jec = hi2%jec
102  hi1%isd = hi2%isd ; hi1%ied = hi2%ied ; hi1%jsd = hi2%jsd ; hi1%jed = hi2%jed
103  hi1%isg = hi2%isg ; hi1%ieg = hi2%ieg ; hi1%jsg = hi2%jsg ; hi1%jeg = hi2%jeg
104 
105  hi1%IscB = hi2%IscB ; hi1%IecB = hi2%IecB ; hi1%JscB = hi2%JscB ; hi1%JecB = hi2%JecB
106  hi1%IsdB = hi2%IsdB ; hi1%IedB = hi2%IedB ; hi1%JsdB = hi2%JsdB ; hi1%JedB = hi2%JedB
107  hi1%IsgB = hi2%IsgB ; hi1%IegB = hi2%IegB ; hi1%JsgB = hi2%JsgB ; hi1%JegB = hi2%JegB
108 
109  hi1%idg_offset = hi2%idg_offset ; hi1%jdg_offset = hi2%jdg_offset
110  hi1%symmetric = hi2%symmetric
111 
112 end subroutine hit_assign
113 
114 !> \namespace mom_hor_index
115 !!
116 !! The hor_index_type provides the decalarations and loop ranges for almost all data with horizontal extent.
117 !!
118 !! Declarations and loop ranges should always be coded with the symmetric memory model in mind.
119 !! The non-symmetric memory mode will then also work, albeit with a different (less efficient) communication pattern.
120 !!
121 !! Using the hor_index_type HI:
122 !! - declaration of h-point data is of the form `h(HI%%isd:HI%%ied,HI%%jsd:HI%%jed)`
123 !! - declaration of q-point data is of the form `q(HI%%IsdB:HI%%IedB,HI%%JsdB:HI%%JedB)`
124 !! - declaration of u-point data is of the form `u(HI%%IsdB:HI%%IedB,HI%%jsd:HI%%jed)`
125 !! - declaration of v-point data is of the form `v(HI%%isd:HI%%ied,HI%%JsdB:HI%%JedB)`.
126 !!
127 !! For more detail explanation of horizontal indexing see \ref Horizontal_indexing.
128 
129 
130 end module mom_hor_index
mom_file_parser::log_version
An overloaded interface to log version information about modules.
Definition: MOM_file_parser.F90:109
mom_file_parser::param_file_type
A structure that can be parsed to read and document run-time parameters.
Definition: MOM_file_parser.F90:54
mom_file_parser::get_param
An overloaded interface to read and log the values of various types of parameters.
Definition: MOM_file_parser.F90:102
mom_hor_index
Defines the horizontal index type (hor_index_type) used for providing index ranges.
Definition: MOM_hor_index.F90:2
mom_domains
Describes the decomposed MOM domain and has routines for communications across PEs.
Definition: MOM_domains.F90:2
mom_file_parser
The MOM6 facility to parse input files for runtime parameters.
Definition: MOM_file_parser.F90:2
mom_hor_index::hor_index_type
Container for horizontal index ranges for data, computational and global domains.
Definition: MOM_hor_index.F90:15
mom_domains::mom_domain_type
The MOM_domain_type contains information about the domain decompositoin.
Definition: MOM_domains.F90:99
mom_file_parser::log_param
An overloaded interface to log the values of various types of parameters.
Definition: MOM_file_parser.F90:96
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2