MOM6
MOM_grid.F90
1 !> Provides the ocean grid type
2 module mom_grid
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 use mom_hor_index, only : hor_index_type, hor_index_init
7 use mom_domains, only : mom_domain_type, get_domain_extent, compute_block_extent
8 use mom_domains, only : get_global_shape, get_domain_extent_dsamp2
9 use mom_error_handler, only : mom_error, mom_mesg, fatal
12 
13 implicit none ; private
14 
15 #include <MOM_memory.h>
16 
17 public mom_grid_init, mom_grid_end, set_derived_metrics, set_first_direction
18 public ispointincell, hor_index_type, get_global_grid_size, rescale_grid_bathymetry
19 
20 ! A note on unit descriptions in comments: MOM6 uses units that can be rescaled for dimensional
21 ! consistency testing. These are noted in comments with units like Z, H, L, and T, along with
22 ! their mks counterparts with notation like "a velocity [Z T-1 ~> m s-1]". If the units
23 ! vary with the Boussinesq approximation, the Boussinesq variant is given first.
24 
25 !> Ocean grid type. See mom_grid for details.
26 type, public :: ocean_grid_type
27  type(mom_domain_type), pointer :: domain => null() !< Ocean model domain
28  type(mom_domain_type), pointer :: domain_aux => null() !< A non-symmetric auxiliary domain type.
29  type(hor_index_type) :: hi !< Horizontal index ranges
30  type(hor_index_type) :: hid2 !< Horizontal index ranges for level-2-downsampling
31 
32  integer :: isc !< The start i-index of cell centers within the computational domain
33  integer :: iec !< The end i-index of cell centers within the computational domain
34  integer :: jsc !< The start j-index of cell centers within the computational domain
35  integer :: jec !< The end j-index of cell centers within the computational domain
36 
37  integer :: isd !< The start i-index of cell centers within the data domain
38  integer :: ied !< The end i-index of cell centers within the data domain
39  integer :: jsd !< The start j-index of cell centers within the data domain
40  integer :: jed !< The end j-index of cell centers within the data domain
41 
42  integer :: isg !< The start i-index of cell centers within the global domain
43  integer :: ieg !< The end i-index of cell centers within the global domain
44  integer :: jsg !< The start j-index of cell centers within the global domain
45  integer :: jeg !< The end j-index of cell centers within the global domain
46 
47  integer :: iscb !< The start i-index of cell vertices within the computational domain
48  integer :: iecb !< The end i-index of cell vertices within the computational domain
49  integer :: jscb !< The start j-index of cell vertices within the computational domain
50  integer :: jecb !< The end j-index of cell vertices within the computational domain
51 
52  integer :: isdb !< The start i-index of cell vertices within the data domain
53  integer :: iedb !< The end i-index of cell vertices within the data domain
54  integer :: jsdb !< The start j-index of cell vertices within the data domain
55  integer :: jedb !< The end j-index of cell vertices within the data domain
56 
57  integer :: isgb !< The start i-index of cell vertices within the global domain
58  integer :: iegb !< The end i-index of cell vertices within the global domain
59  integer :: jsgb !< The start j-index of cell vertices within the global domain
60  integer :: jegb !< The end j-index of cell vertices within the global domain
61 
62  integer :: isd_global !< The value of isd in the global index space (decompoistion invariant).
63  integer :: jsd_global !< The value of isd in the global index space (decompoistion invariant).
64  integer :: idg_offset !< The offset between the corresponding global and local i-indices.
65  integer :: jdg_offset !< The offset between the corresponding global and local j-indices.
66  integer :: ke !< The number of layers in the vertical.
67  logical :: symmetric !< True if symmetric memory is used.
68  logical :: nonblocking_updates !< If true, non-blocking halo updates are
69  !! allowed. The default is .false. (for now).
70  integer :: first_direction !< An integer that indicates which direction is
71  !! to be updated first in directionally split
72  !! parts of the calculation. This can be altered
73  !! during the course of the run via calls to
74  !! set_first_direction.
75 
76  real allocable_, dimension(NIMEM_,NJMEM_) :: &
77  mask2dt, & !< 0 for land points and 1 for ocean points on the h-grid [nondim].
78  geolatt, & !< The geographic latitude at q points in degrees of latitude or m.
79  geolont, & !< The geographic longitude at q points in degrees of longitude or m.
80  dxt, & !< dxT is delta x at h points [L ~> m].
81  idxt, & !< 1/dxT [L-1 ~> m-1].
82  dyt, & !< dyT is delta y at h points [L ~> m].
83  idyt, & !< IdyT is 1/dyT [L-1 ~> m-1].
84  areat, & !< The area of an h-cell [L2 ~> m2].
85  iareat, & !< 1/areaT [L-2 ~> m-2].
86  sin_rot, & !< The sine of the angular rotation between the local model grid's northward
87  !! and the true northward directions.
88  cos_rot !< The cosine of the angular rotation between the local model grid's northward
89  !! and the true northward directions.
90 
91  real allocable_, dimension(NIMEMB_PTR_,NJMEM_) :: &
92  mask2dcu, & !< 0 for boundary points and 1 for ocean points on the u grid [nondim].
93  geolatcu, & !< The geographic latitude at u points in degrees of latitude or m.
94  geoloncu, & !< The geographic longitude at u points in degrees of longitude or m.
95  dxcu, & !< dxCu is delta x at u points [L ~> m].
96  idxcu, & !< 1/dxCu [L-1 ~> m-1].
97  dycu, & !< dyCu is delta y at u points [L ~> m].
98  idycu, & !< 1/dyCu [L-1 ~> m-1].
99  dy_cu, & !< The unblocked lengths of the u-faces of the h-cell [L ~> m].
100  iareacu, & !< The masked inverse areas of u-grid cells [L-2 ~> m-2].
101  areacu !< The areas of the u-grid cells [L2 ~> m2].
102 
103  real allocable_, dimension(NIMEM_,NJMEMB_PTR_) :: &
104  mask2dcv, & !< 0 for boundary points and 1 for ocean points on the v grid [nondim].
105  geolatcv, & !< The geographic latitude at v points in degrees of latitude or m.
106  geoloncv, & !< The geographic longitude at v points in degrees of longitude or m.
107  dxcv, & !< dxCv is delta x at v points [L ~> m].
108  idxcv, & !< 1/dxCv [L-1 ~> m-1].
109  dycv, & !< dyCv is delta y at v points [L ~> m].
110  idycv, & !< 1/dyCv [L-1 ~> m-1].
111  dx_cv, & !< The unblocked lengths of the v-faces of the h-cell [L ~> m].
112  iareacv, & !< The masked inverse areas of v-grid cells [L-2 ~> m-2].
113  areacv !< The areas of the v-grid cells [L2 ~> m2].
114 
115  real allocable_, dimension(NIMEMB_PTR_,NJMEMB_PTR_) :: &
116  mask2dbu, & !< 0 for boundary points and 1 for ocean points on the q grid [nondim].
117  geolatbu, & !< The geographic latitude at q points in degrees of latitude or m.
118  geolonbu, & !< The geographic longitude at q points in degrees of longitude or m.
119  dxbu, & !< dxBu is delta x at q points [L ~> m].
120  idxbu, & !< 1/dxBu [L-1 ~> m-1].
121  dybu, & !< dyBu is delta y at q points [L ~> m].
122  idybu, & !< 1/dyBu [L-1 ~> m-1].
123  areabu, & !< areaBu is the area of a q-cell [L2 ~> m2]
124  iareabu !< IareaBu = 1/areaBu [L-2 ~> m-2].
125 
126  real, pointer, dimension(:) :: &
127  gridlatt => null(), & !< The latitude of T points for the purpose of labeling the output axes.
128  !! On many grids this is the same as geoLatT.
129  gridlatb => null() !< The latitude of B points for the purpose of labeling the output axes.
130  !! On many grids this is the same as geoLatBu.
131  real, pointer, dimension(:) :: &
132  gridlont => null(), & !< The longitude of T points for the purpose of labeling the output axes.
133  !! On many grids this is the same as geoLonT.
134  gridlonb => null() !< The longitude of B points for the purpose of labeling the output axes.
135  !! On many grids this is the same as geoLonBu.
136  character(len=40) :: &
137  x_axis_units, & !< The units that are used in labeling the x coordinate axes.
138  y_axis_units !< The units that are used in labeling the y coordinate axes.
139 
140  real allocable_, dimension(NIMEM_,NJMEM_) :: &
141  bathyt !< Ocean bottom depth at tracer points, in depth units [Z ~> m].
142  real :: z_ref !< A reference value for all geometric height fields, such as bathyT [Z ~> m].
143 
144  logical :: bathymetry_at_vel !< If true, there are separate values for the
145  !! basin depths at velocity points. Otherwise the effects of
146  !! of topography are entirely determined from thickness points.
147  real allocable_, dimension(NIMEMB_PTR_,NJMEM_) :: &
148  dblock_u, & !< Topographic depths at u-points at which the flow is blocked [Z ~> m].
149  dopen_u !< Topographic depths at u-points at which the flow is open at width dy_Cu [Z ~> m].
150  real allocable_, dimension(NIMEM_,NJMEMB_PTR_) :: &
151  dblock_v, & !< Topographic depths at v-points at which the flow is blocked [Z ~> m].
152  dopen_v !< Topographic depths at v-points at which the flow is open at width dx_Cv [Z ~> m].
153  real allocable_, dimension(NIMEMB_PTR_,NJMEMB_PTR_) :: &
154  coriolisbu !< The Coriolis parameter at corner points [T-1 ~> s-1].
155  real allocable_, dimension(NIMEM_,NJMEM_) :: &
156  df_dx, & !< Derivative d/dx f (Coriolis parameter) at h-points [T-1 L-1 ~> s-1 m-1].
157  df_dy !< Derivative d/dy f (Coriolis parameter) at h-points [T-1 L-1 ~> s-1 m-1].
158 
159  ! These variables are global sums that are useful for 1-d diagnostics and should not be rescaled.
160  real :: areat_global !< Global sum of h-cell area [m2]
161  real :: iareat_global !< Global sum of inverse h-cell area (1/areaT_global) [m-2].
162 
163  type(unit_scale_type), pointer :: us => null() !< A dimensional unit scaling type
164 
165 
166  ! These variables are for block structures.
167  integer :: nblocks !< The number of sub-PE blocks on this PE
168  type(hor_index_type), pointer :: block(:) => null() !< Index ranges for each block
169 
170  ! These parameters are run-time parameters that are used during some
171  ! initialization routines (but not all)
172  real :: south_lat !< The latitude (or y-coordinate) of the first v-line
173  real :: west_lon !< The longitude (or x-coordinate) of the first u-line
174  real :: len_lat = 0. !< The latitudinal (or y-coord) extent of physical domain
175  real :: len_lon = 0. !< The longitudinal (or x-coord) extent of physical domain
176  real :: rad_earth = 6.378e6 !< The radius of the planet [m].
177  real :: max_depth !< The maximum depth of the ocean in depth units [Z ~> m].
178 end type ocean_grid_type
179 
180 contains
181 
182 !> MOM_grid_init initializes the ocean grid array sizes and grid memory.
183 subroutine mom_grid_init(G, param_file, US, HI, global_indexing, bathymetry_at_vel)
184  type(ocean_grid_type), intent(inout) :: g !< The horizontal grid type
185  type(param_file_type), intent(in) :: param_file !< Parameter file handle
186  type(unit_scale_type), optional, pointer :: us !< A dimensional unit scaling type
187  type(hor_index_type), &
188  optional, intent(in) :: hi !< A hor_index_type for array extents
189  logical, optional, intent(in) :: global_indexing !< If true use global index
190  !! values instead of having the data domain on each
191  !! processor start at 1.
192  logical, optional, intent(in) :: bathymetry_at_vel !< If true, there are
193  !! separate values for the ocean bottom depths at
194  !! velocity points. Otherwise the effects of topography
195  !! are entirely determined from thickness points.
196 
197  ! Local variables
198  real :: mean_sealev_scale
199  integer :: isd, ied, jsd, jed, nk
200  integer :: isdb, iedb, jsdb, jedb
201  integer :: ied_max, jed_max
202  integer :: niblock, njblock, nihalo, njhalo, nblocks, n, i, j
203  logical :: local_indexing ! If false use global index values instead of having
204  ! the data domain on each processor start at 1.
205  ! This include declares and sets the variable "version".
206 # include "version_variable.h"
207 
208  integer, allocatable, dimension(:) :: ibegin, iend, jbegin, jend
209  character(len=40) :: mod_nm = "MOM_grid" ! This module's name.
210 
211 
212  ! Read all relevant parameters and write them to the model log.
213  call get_param(param_file, mod_nm, "REFERENCE_HEIGHT", g%Z_ref, default=0.0, do_not_log=.true.)
214  call log_version(param_file, mod_nm, version, &
215  "Parameters providing information about the lateral grid.", &
216  log_to_all=.true., layout=.true., all_default=(g%Z_ref==0.0))
217 
218  call get_param(param_file, mod_nm, "NIBLOCK", niblock, "The number of blocks "// &
219  "in the x-direction on each processor (for openmp).", default=1, &
220  layoutparam=.true.)
221  call get_param(param_file, mod_nm, "NJBLOCK", njblock, "The number of blocks "// &
222  "in the y-direction on each processor (for openmp).", default=1, &
223  layoutparam=.true.)
224  if (present(us)) then ; if (associated(us)) g%US => us ; endif
225 
226  mean_sealev_scale = 1.0 ; if (associated(g%US)) mean_sealev_scale = g%US%m_to_Z
227  call get_param(param_file, mod_nm, "REFERENCE_HEIGHT", g%Z_ref, &
228  "A reference value for geometric height fields, such as bathyT.", &
229  units="m", default=0.0, scale=mean_sealev_scale)
230 
231  if (present(hi)) then
232  g%HI = hi
233 
234  g%isc = hi%isc ; g%iec = hi%iec ; g%jsc = hi%jsc ; g%jec = hi%jec
235  g%isd = hi%isd ; g%ied = hi%ied ; g%jsd = hi%jsd ; g%jed = hi%jed
236  g%isg = hi%isg ; g%ieg = hi%ieg ; g%jsg = hi%jsg ; g%jeg = hi%jeg
237 
238  g%IscB = hi%IscB ; g%IecB = hi%IecB ; g%JscB = hi%JscB ; g%JecB = hi%JecB
239  g%IsdB = hi%IsdB ; g%IedB = hi%IedB ; g%JsdB = hi%JsdB ; g%JedB = hi%JedB
240  g%IsgB = hi%IsgB ; g%IegB = hi%IegB ; g%JsgB = hi%JsgB ; g%JegB = hi%JegB
241 
242  g%idg_offset = hi%idg_offset ; g%jdg_offset = hi%jdg_offset
243  g%isd_global = g%isd + hi%idg_offset ; g%jsd_global = g%jsd + hi%jdg_offset
244  g%symmetric = hi%symmetric
245  else
246  local_indexing = .true.
247  if (present(global_indexing)) local_indexing = .not.global_indexing
248  call hor_index_init(g%Domain, g%HI, param_file, &
249  local_indexing=local_indexing)
250 
251  ! get_domain_extent ensures that domains start at 1 for compatibility between
252  ! static and dynamically allocated arrays, unless global_indexing is true.
253  call get_domain_extent(g%Domain, g%isc, g%iec, g%jsc, g%jec, &
254  g%isd, g%ied, g%jsd, g%jed, &
255  g%isg, g%ieg, g%jsg, g%jeg, &
256  g%idg_offset, g%jdg_offset, g%symmetric, &
257  local_indexing=local_indexing)
258  g%isd_global = g%isd+g%idg_offset ; g%jsd_global = g%jsd+g%jdg_offset
259  endif
260 
261  g%nonblocking_updates = g%Domain%nonblocking_updates
262 
263  ! Set array sizes for fields that are discretized at tracer cell boundaries.
264  g%IscB = g%isc ; g%JscB = g%jsc
265  g%IsdB = g%isd ; g%JsdB = g%jsd
266  g%IsgB = g%isg ; g%JsgB = g%jsg
267  if (g%symmetric) then
268  g%IscB = g%isc-1 ; g%JscB = g%jsc-1
269  g%IsdB = g%isd-1 ; g%JsdB = g%jsd-1
270  g%IsgB = g%isg-1 ; g%JsgB = g%jsg-1
271  endif
272  g%IecB = g%iec ; g%JecB = g%jec
273  g%IedB = g%ied ; g%JedB = g%jed
274  g%IegB = g%ieg ; g%JegB = g%jeg
275 
276  call mom_mesg(" MOM_grid.F90, MOM_grid_init: allocating metrics", 5)
277 
278  call allocate_metrics(g)
279 
280  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
281  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
282 
283  g%bathymetry_at_vel = .false.
284  if (present(bathymetry_at_vel)) g%bathymetry_at_vel = bathymetry_at_vel
285  if (g%bathymetry_at_vel) then
286  alloc_(g%Dblock_u(isdb:iedb, jsd:jed)) ; g%Dblock_u(:,:) = 0.0
287  alloc_(g%Dopen_u(isdb:iedb, jsd:jed)) ; g%Dopen_u(:,:) = 0.0
288  alloc_(g%Dblock_v(isd:ied, jsdb:jedb)) ; g%Dblock_v(:,:) = 0.0
289  alloc_(g%Dopen_v(isd:ied, jsdb:jedb)) ; g%Dopen_v(:,:) = 0.0
290  endif
291 
292 ! setup block indices.
293  nihalo = g%Domain%nihalo
294  njhalo = g%Domain%njhalo
295  nblocks = niblock * njblock
296  if (nblocks < 1) call mom_error(fatal, "MOM_grid_init: " // &
297  "nblocks(=NI_BLOCK*NJ_BLOCK) must be no less than 1")
298 
299  allocate(ibegin(niblock), iend(niblock), jbegin(njblock), jend(njblock))
300  call compute_block_extent(g%HI%isc,g%HI%iec,niblock,ibegin,iend)
301  call compute_block_extent(g%HI%jsc,g%HI%jec,njblock,jbegin,jend)
302  !-- make sure the last block is the largest.
303  do i = 1, niblock-1
304  if (iend(i)-ibegin(i) > iend(niblock)-ibegin(niblock) ) call mom_error(fatal, &
305  "MOM_grid_init: the last block size in x-direction is not the largest")
306  enddo
307  do j = 1, njblock-1
308  if (jend(j)-jbegin(j) > jend(njblock)-jbegin(njblock) ) call mom_error(fatal, &
309  "MOM_grid_init: the last block size in y-direction is not the largest")
310  enddo
311 
312  g%nblocks = nblocks
313  allocate(g%Block(nblocks))
314  ied_max = 1 ; jed_max = 1
315  do n = 1,nblocks
316  ! Copy all information from the array index type describing the local grid.
317  g%Block(n) = g%HI
318 
319  i = mod((n-1), niblock) + 1
320  j = (n-1)/niblock + 1
321  !--- isd and jsd are always 1 for each block to permit array reuse.
322  g%Block(n)%isd = 1 ; g%Block(n)%jsd = 1
323  g%Block(n)%isc = g%Block(n)%isd+nihalo
324  g%Block(n)%jsc = g%Block(n)%jsd+njhalo
325  g%Block(n)%iec = g%Block(n)%isc + iend(i) - ibegin(i)
326  g%Block(n)%jec = g%Block(n)%jsc + jend(j) - jbegin(j)
327  g%Block(n)%ied = g%Block(n)%iec + nihalo
328  g%Block(n)%jed = g%Block(n)%jec + njhalo
329  g%Block(n)%IscB = g%Block(n)%isc; g%Block(n)%IecB = g%Block(n)%iec
330  g%Block(n)%JscB = g%Block(n)%jsc; g%Block(n)%JecB = g%Block(n)%jec
331  ! For symmetric memory domains, the first block will have the extra point
332  ! at the lower boundary of its computational domain.
333  if (g%symmetric) then
334  if (i==1) g%Block(n)%IscB = g%Block(n)%IscB-1
335  if (j==1) g%Block(n)%JscB = g%Block(n)%JscB-1
336  endif
337  g%Block(n)%IsdB = g%Block(n)%isd; g%Block(n)%IedB = g%Block(n)%ied
338  g%Block(n)%JsdB = g%Block(n)%jsd; g%Block(n)%JedB = g%Block(n)%jed
339  !--- For symmetric memory domain, every block will have an extra point
340  !--- at the lower boundary of its data domain.
341  if (g%symmetric) then
342  g%Block(n)%IsdB = g%Block(n)%IsdB-1
343  g%Block(n)%JsdB = g%Block(n)%JsdB-1
344  endif
345  g%Block(n)%idg_offset = (ibegin(i) - g%Block(n)%isc) + g%HI%idg_offset
346  g%Block(n)%jdg_offset = (jbegin(j) - g%Block(n)%jsc) + g%HI%jdg_offset
347  ! Find the largest values of ied and jed so that all blocks will have the
348  ! same size in memory.
349  ied_max = max(ied_max, g%Block(n)%ied)
350  jed_max = max(jed_max, g%Block(n)%jed)
351  enddo
352 
353  ! Reset all of the data domain sizes to match the largest for array reuse,
354  ! recalling that all block have isd=jed=1 for array reuse.
355  do n = 1,nblocks
356  g%Block(n)%ied = ied_max ; g%Block(n)%IedB = ied_max
357  g%Block(n)%jed = jed_max ; g%Block(n)%JedB = jed_max
358  enddo
359 
360  !-- do some bounds error checking
361  if ( g%block(nblocks)%ied+g%block(nblocks)%idg_offset > g%HI%ied + g%HI%idg_offset ) &
362  call mom_error(fatal, "MOM_grid_init: G%ied_bk > G%ied")
363  if ( g%block(nblocks)%jed+g%block(nblocks)%jdg_offset > g%HI%jed + g%HI%jdg_offset ) &
364  call mom_error(fatal, "MOM_grid_init: G%jed_bk > G%jed")
365 
366  call get_domain_extent_dsamp2(g%Domain, g%HId2%isc, g%HId2%iec, g%HId2%jsc, g%HId2%jec,&
367  g%HId2%isd, g%HId2%ied, g%HId2%jsd, g%HId2%jed,&
368  g%HId2%isg, g%HId2%ieg, g%HId2%jsg, g%HId2%jeg)
369 
370  ! Set array sizes for fields that are discretized at tracer cell boundaries.
371  g%HId2%IscB = g%HId2%isc ; g%HId2%JscB = g%HId2%jsc
372  g%HId2%IsdB = g%HId2%isd ; g%HId2%JsdB = g%HId2%jsd
373  g%HId2%IsgB = g%HId2%isg ; g%HId2%JsgB = g%HId2%jsg
374  if (g%symmetric) then
375  g%HId2%IscB = g%HId2%isc-1 ; g%HId2%JscB = g%HId2%jsc-1
376  g%HId2%IsdB = g%HId2%isd-1 ; g%HId2%JsdB = g%HId2%jsd-1
377  g%HId2%IsgB = g%HId2%isg-1 ; g%HId2%JsgB = g%HId2%jsg-1
378  endif
379  g%HId2%IecB = g%HId2%iec ; g%HId2%JecB = g%HId2%jec
380  g%HId2%IedB = g%HId2%ied ; g%HId2%JedB = g%HId2%jed
381  g%HId2%IegB = g%HId2%ieg ; g%HId2%JegB = g%HId2%jeg
382 
383 end subroutine mom_grid_init
384 
385 !> rescale_grid_bathymetry permits a change in the internal units for the bathymetry on the grid,
386 !! both rescaling the depths and recording the new internal units.
387 subroutine rescale_grid_bathymetry(G, m_in_new_units)
388  type(ocean_grid_type), intent(inout) :: g !< The horizontal grid structure
389  real, intent(in) :: m_in_new_units !< The new internal representation of 1 m depth.
390 
391  ! Local variables
392  real :: rescale
393  integer :: i, j, isd, ied, jsd, jed, isdb, iedb, jsdb, jedb
394 
395  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
396  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
397 
398  if (m_in_new_units == 1.0) return
399  if (m_in_new_units < 0.0) &
400  call mom_error(fatal, "rescale_grid_bathymetry: Negative depth units are not permitted.")
401  if (m_in_new_units == 0.0) &
402  call mom_error(fatal, "rescale_grid_bathymetry: Zero depth units are not permitted.")
403 
404  rescale = 1.0 / m_in_new_units
405  do j=jsd,jed ; do i=isd,ied
406  g%bathyT(i,j) = rescale*g%bathyT(i,j)
407  enddo ; enddo
408  if (g%bathymetry_at_vel) then ; do j=jsd,jed ; do i=isdb,iedb
409  g%Dblock_u(i,j) = rescale*g%Dblock_u(i,j) ; g%Dopen_u(i,j) = rescale*g%Dopen_u(i,j)
410  enddo ; enddo ; endif
411  if (g%bathymetry_at_vel) then ; do j=jsdb,jedb ; do i=isd,ied
412  g%Dblock_v(i,j) = rescale*g%Dblock_v(i,j) ; g%Dopen_v(i,j) = rescale*g%Dopen_v(i,j)
413  enddo ; enddo ; endif
414  g%max_depth = rescale*g%max_depth
415 
416 end subroutine rescale_grid_bathymetry
417 
418 !> set_derived_metrics calculates metric terms that are derived from other metrics.
419 subroutine set_derived_metrics(G, US)
420  type(ocean_grid_type), intent(inout) :: g !< The horizontal grid structure
421  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
422 ! Various inverse grid spacings and derived areas are calculated within this
423 ! subroutine.
424  integer :: i, j, isd, ied, jsd, jed
425  integer :: isdb, iedb, jsdb, jedb
426 
427  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
428  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
429 
430  do j=jsd,jed ; do i=isd,ied
431  if (g%dxT(i,j) < 0.0) g%dxT(i,j) = 0.0
432  if (g%dyT(i,j) < 0.0) g%dyT(i,j) = 0.0
433  g%IdxT(i,j) = adcroft_reciprocal(g%dxT(i,j))
434  g%IdyT(i,j) = adcroft_reciprocal(g%dyT(i,j))
435  g%IareaT(i,j) = adcroft_reciprocal(g%areaT(i,j))
436  enddo ; enddo
437 
438  do j=jsd,jed ; do i=isdb,iedb
439  if (g%dxCu(i,j) < 0.0) g%dxCu(i,j) = 0.0
440  if (g%dyCu(i,j) < 0.0) g%dyCu(i,j) = 0.0
441  g%IdxCu(i,j) = adcroft_reciprocal(g%dxCu(i,j))
442  g%IdyCu(i,j) = adcroft_reciprocal(g%dyCu(i,j))
443  enddo ; enddo
444 
445  do j=jsdb,jedb ; do i=isd,ied
446  if (g%dxCv(i,j) < 0.0) g%dxCv(i,j) = 0.0
447  if (g%dyCv(i,j) < 0.0) g%dyCv(i,j) = 0.0
448  g%IdxCv(i,j) = adcroft_reciprocal(g%dxCv(i,j))
449  g%IdyCv(i,j) = adcroft_reciprocal(g%dyCv(i,j))
450  enddo ; enddo
451 
452  do j=jsdb,jedb ; do i=isdb,iedb
453  if (g%dxBu(i,j) < 0.0) g%dxBu(i,j) = 0.0
454  if (g%dyBu(i,j) < 0.0) g%dyBu(i,j) = 0.0
455 
456  g%IdxBu(i,j) = adcroft_reciprocal(g%dxBu(i,j))
457  g%IdyBu(i,j) = adcroft_reciprocal(g%dyBu(i,j))
458  ! areaBu has usually been set to a positive area elsewhere.
459  if (g%areaBu(i,j) <= 0.0) g%areaBu(i,j) = g%dxBu(i,j) * g%dyBu(i,j)
460  g%IareaBu(i,j) = adcroft_reciprocal(g%areaBu(i,j))
461  enddo ; enddo
462 end subroutine set_derived_metrics
463 
464 !> Adcroft_reciprocal(x) = 1/x for |x|>0 or 0 for x=0.
465 function adcroft_reciprocal(val) result(I_val)
466  real, intent(in) :: val !< The value being inverted.
467  real :: i_val !< The Adcroft reciprocal of val.
468 
469  i_val = 0.0 ; if (val /= 0.0) i_val = 1.0/val
470 end function adcroft_reciprocal
471 
472 !> Returns true if the coordinates (x,y) are within the h-cell (i,j)
473 logical function ispointincell(G, i, j, x, y)
474  type(ocean_grid_type), intent(in) :: g !< Grid type
475  integer, intent(in) :: i !< i index of cell to test
476  integer, intent(in) :: j !< j index of cell to test
477  real, intent(in) :: x !< x coordinate of point
478  real, intent(in) :: y !< y coordinate of point
479  ! Local variables
480  real :: xne, xnw, xse, xsw, yne, ynw, yse, ysw
481  real :: p0, p1, p2, p3, l0, l1, l2, l3
482  ispointincell = .false.
483  xne = g%geoLonBu(i ,j ) ; yne = g%geoLatBu(i ,j )
484  xnw = g%geoLonBu(i-1,j ) ; ynw = g%geoLatBu(i-1,j )
485  xse = g%geoLonBu(i ,j-1) ; yse = g%geoLatBu(i ,j-1)
486  xsw = g%geoLonBu(i-1,j-1) ; ysw = g%geoLatBu(i-1,j-1)
487  ! This is a crude calculation that assume a geographic coordinate system
488  if (x<min(xne,xnw,xse,xsw) .or. x>max(xne,xnw,xse,xsw) .or. &
489  y<min(yne,ynw,yse,ysw) .or. y>max(yne,ynw,yse,ysw) ) then
490  return ! Avoid the more complicated calculation
491  endif
492  l0 = (x-xsw)*(yse-ysw) - (y-ysw)*(xse-xsw)
493  l1 = (x-xse)*(yne-yse) - (y-yse)*(xne-xse)
494  l2 = (x-xne)*(ynw-yne) - (y-yne)*(xnw-xne)
495  l3 = (x-xnw)*(ysw-ynw) - (y-ynw)*(xsw-xnw)
496 
497  p0 = sign(1., l0) ; if (l0 == 0.) p0=0.
498  p1 = sign(1., l1) ; if (l1 == 0.) p1=0.
499  p2 = sign(1., l2) ; if (l2 == 0.) p2=0.
500  p3 = sign(1., l3) ; if (l3 == 0.) p3=0.
501 
502  if ( (abs(p0)+abs(p2)) + (abs(p1)+abs(p3)) == abs((p0+p2) + (p1+p3)) ) then
503  ispointincell=.true.
504  endif
505 end function ispointincell
506 
507 !> Store an integer indicating which direction to work on first.
508 subroutine set_first_direction(G, y_first)
509  type(ocean_grid_type), intent(inout) :: g !< The ocean's grid structure
510  integer, intent(in) :: y_first !< The first direction to store
511 
512  g%first_direction = y_first
513 end subroutine set_first_direction
514 
515 !> Return global shape of horizontal grid
516 subroutine get_global_grid_size(G, niglobal, njglobal)
517  type(ocean_grid_type), intent(inout) :: g !< The horizontal grid type
518  integer, intent(out) :: niglobal !< i-index global size of grid
519  integer, intent(out) :: njglobal !< j-index global size of grid
520 
521  call get_global_shape(g%domain, niglobal, njglobal)
522 
523 end subroutine get_global_grid_size
524 
525 !> Allocate memory used by the ocean_grid_type and related structures.
526 subroutine allocate_metrics(G)
527  type(ocean_grid_type), intent(inout) :: G !< The horizontal grid type
528  integer :: isd, ied, jsd, jed, IsdB, IedB, JsdB, JedB, isg, ieg, jsg, jeg
529 
530  ! This subroutine allocates the lateral elements of the ocean_grid_type that
531  ! are always used and zeros them out.
532 
533  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
534  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
535  isg = g%isg ; ieg = g%ieg ; jsg = g%jsg ; jeg = g%jeg
536 
537  alloc_(g%dxT(isd:ied,jsd:jed)) ; g%dxT(:,:) = 0.0
538  alloc_(g%dxCu(isdb:iedb,jsd:jed)) ; g%dxCu(:,:) = 0.0
539  alloc_(g%dxCv(isd:ied,jsdb:jedb)) ; g%dxCv(:,:) = 0.0
540  alloc_(g%dxBu(isdb:iedb,jsdb:jedb)) ; g%dxBu(:,:) = 0.0
541  alloc_(g%IdxT(isd:ied,jsd:jed)) ; g%IdxT(:,:) = 0.0
542  alloc_(g%IdxCu(isdb:iedb,jsd:jed)) ; g%IdxCu(:,:) = 0.0
543  alloc_(g%IdxCv(isd:ied,jsdb:jedb)) ; g%IdxCv(:,:) = 0.0
544  alloc_(g%IdxBu(isdb:iedb,jsdb:jedb)) ; g%IdxBu(:,:) = 0.0
545 
546  alloc_(g%dyT(isd:ied,jsd:jed)) ; g%dyT(:,:) = 0.0
547  alloc_(g%dyCu(isdb:iedb,jsd:jed)) ; g%dyCu(:,:) = 0.0
548  alloc_(g%dyCv(isd:ied,jsdb:jedb)) ; g%dyCv(:,:) = 0.0
549  alloc_(g%dyBu(isdb:iedb,jsdb:jedb)) ; g%dyBu(:,:) = 0.0
550  alloc_(g%IdyT(isd:ied,jsd:jed)) ; g%IdyT(:,:) = 0.0
551  alloc_(g%IdyCu(isdb:iedb,jsd:jed)) ; g%IdyCu(:,:) = 0.0
552  alloc_(g%IdyCv(isd:ied,jsdb:jedb)) ; g%IdyCv(:,:) = 0.0
553  alloc_(g%IdyBu(isdb:iedb,jsdb:jedb)) ; g%IdyBu(:,:) = 0.0
554 
555  alloc_(g%areaT(isd:ied,jsd:jed)) ; g%areaT(:,:) = 0.0
556  alloc_(g%IareaT(isd:ied,jsd:jed)) ; g%IareaT(:,:) = 0.0
557  alloc_(g%areaBu(isdb:iedb,jsdb:jedb)) ; g%areaBu(:,:) = 0.0
558  alloc_(g%IareaBu(isdb:iedb,jsdb:jedb)) ; g%IareaBu(:,:) = 0.0
559 
560  alloc_(g%mask2dT(isd:ied,jsd:jed)) ; g%mask2dT(:,:) = 0.0
561  alloc_(g%mask2dCu(isdb:iedb,jsd:jed)) ; g%mask2dCu(:,:) = 0.0
562  alloc_(g%mask2dCv(isd:ied,jsdb:jedb)) ; g%mask2dCv(:,:) = 0.0
563  alloc_(g%mask2dBu(isdb:iedb,jsdb:jedb)) ; g%mask2dBu(:,:) = 0.0
564  alloc_(g%geoLatT(isd:ied,jsd:jed)) ; g%geoLatT(:,:) = 0.0
565  alloc_(g%geoLatCu(isdb:iedb,jsd:jed)) ; g%geoLatCu(:,:) = 0.0
566  alloc_(g%geoLatCv(isd:ied,jsdb:jedb)) ; g%geoLatCv(:,:) = 0.0
567  alloc_(g%geoLatBu(isdb:iedb,jsdb:jedb)) ; g%geoLatBu(:,:) = 0.0
568  alloc_(g%geoLonT(isd:ied,jsd:jed)) ; g%geoLonT(:,:) = 0.0
569  alloc_(g%geoLonCu(isdb:iedb,jsd:jed)) ; g%geoLonCu(:,:) = 0.0
570  alloc_(g%geoLonCv(isd:ied,jsdb:jedb)) ; g%geoLonCv(:,:) = 0.0
571  alloc_(g%geoLonBu(isdb:iedb,jsdb:jedb)) ; g%geoLonBu(:,:) = 0.0
572 
573  alloc_(g%dx_Cv(isd:ied,jsdb:jedb)) ; g%dx_Cv(:,:) = 0.0
574  alloc_(g%dy_Cu(isdb:iedb,jsd:jed)) ; g%dy_Cu(:,:) = 0.0
575 
576  alloc_(g%areaCu(isdb:iedb,jsd:jed)) ; g%areaCu(:,:) = 0.0
577  alloc_(g%areaCv(isd:ied,jsdb:jedb)) ; g%areaCv(:,:) = 0.0
578  alloc_(g%IareaCu(isdb:iedb,jsd:jed)) ; g%IareaCu(:,:) = 0.0
579  alloc_(g%IareaCv(isd:ied,jsdb:jedb)) ; g%IareaCv(:,:) = 0.0
580 
581  alloc_(g%bathyT(isd:ied, jsd:jed)) ; g%bathyT(:,:) = 0.0
582  alloc_(g%CoriolisBu(isdb:iedb, jsdb:jedb)) ; g%CoriolisBu(:,:) = 0.0
583  alloc_(g%dF_dx(isd:ied, jsd:jed)) ; g%dF_dx(:,:) = 0.0
584  alloc_(g%dF_dy(isd:ied, jsd:jed)) ; g%dF_dy(:,:) = 0.0
585 
586  alloc_(g%sin_rot(isd:ied,jsd:jed)) ; g%sin_rot(:,:) = 0.0
587  alloc_(g%cos_rot(isd:ied,jsd:jed)) ; g%cos_rot(:,:) = 1.0
588 
589  allocate(g%gridLonT(isg:ieg)) ; g%gridLonT(:) = 0.0
590  allocate(g%gridLonB(g%IsgB:g%IegB)) ; g%gridLonB(:) = 0.0
591  allocate(g%gridLatT(jsg:jeg)) ; g%gridLatT(:) = 0.0
592  allocate(g%gridLatB(g%JsgB:g%JegB)) ; g%gridLatB(:) = 0.0
593 
594 end subroutine allocate_metrics
595 
596 !> Release memory used by the ocean_grid_type and related structures.
597 subroutine mom_grid_end(G)
598  type(ocean_grid_type), intent(inout) :: g !< The horizontal grid type
599 
600  dealloc_(g%dxT) ; dealloc_(g%dxCu) ; dealloc_(g%dxCv) ; dealloc_(g%dxBu)
601  dealloc_(g%IdxT) ; dealloc_(g%IdxCu) ; dealloc_(g%IdxCv) ; dealloc_(g%IdxBu)
602 
603  dealloc_(g%dyT) ; dealloc_(g%dyCu) ; dealloc_(g%dyCv) ; dealloc_(g%dyBu)
604  dealloc_(g%IdyT) ; dealloc_(g%IdyCu) ; dealloc_(g%IdyCv) ; dealloc_(g%IdyBu)
605 
606  dealloc_(g%areaT) ; dealloc_(g%IareaT)
607  dealloc_(g%areaBu) ; dealloc_(g%IareaBu)
608  dealloc_(g%areaCu) ; dealloc_(g%IareaCu)
609  dealloc_(g%areaCv) ; dealloc_(g%IareaCv)
610 
611  dealloc_(g%mask2dT) ; dealloc_(g%mask2dCu)
612  dealloc_(g%mask2dCv) ; dealloc_(g%mask2dBu)
613 
614  dealloc_(g%geoLatT) ; dealloc_(g%geoLatCu)
615  dealloc_(g%geoLatCv) ; dealloc_(g%geoLatBu)
616  dealloc_(g%geoLonT) ; dealloc_(g%geoLonCu)
617  dealloc_(g%geoLonCv) ; dealloc_(g%geoLonBu)
618 
619  dealloc_(g%dx_Cv) ; dealloc_(g%dy_Cu)
620 
621  dealloc_(g%bathyT) ; dealloc_(g%CoriolisBu)
622  dealloc_(g%dF_dx) ; dealloc_(g%dF_dy)
623  dealloc_(g%sin_rot) ; dealloc_(g%cos_rot)
624 
625  if (g%bathymetry_at_vel) then
626  dealloc_(g%Dblock_u) ; dealloc_(g%Dopen_u)
627  dealloc_(g%Dblock_v) ; dealloc_(g%Dopen_v)
628  endif
629 
630  deallocate(g%gridLonT) ; deallocate(g%gridLatT)
631  deallocate(g%gridLonB) ; deallocate(g%gridLatB)
632 
633  deallocate(g%Domain%mpp_domain)
634  deallocate(g%Domain)
635 
636 end subroutine mom_grid_end
637 
638 !> \namespace mom_grid
639 !!
640 !! Grid metrics and their inverses are labelled according to their staggered location on a Arakawa C (or B) grid.
641 !! - Metrics centered on h- or T-points are labelled T, e.g. dxT is the distance across the cell in the x-direction.
642 !! - Metrics centered on u-points are labelled Cu (C-grid u location). e.g. dyCu is the y-distance between
643 !! two corners of a T-cell.
644 !! - Metrics centered on v-points are labelled Cv (C-grid v location). e.g. dyCv is the y-distance between two -points.
645 !! - Metrics centered on q-points are labelled Bu (B-grid u,v location). e.g. areaBu is the area centered on a q-point.
646 !!
647 !! \image html Grid_metrics.png "The labelling of distances (grid metrics) at various staggered
648 !! location on an T-cell and around a q-point."
649 !!
650 !! Areas centered at T-, u-, v- and q- points are `areaT`, `areaCu`, `areaCv` and `areaBu` respectively.
651 !!
652 !! The reciprocal of metrics are pre-calculated and also stored in the ocean_grid_type with a I prepended to the name.
653 !! For example, `1./areaT` is called `IareaT`, and `1./dyCv` is `IdyCv`.
654 !!
655 !! Geographic latitude and longitude (or model coordinates if not on a sphere) are stored in
656 !! `geoLatT`, `geoLonT` for T-points.
657 !! u-, v- and q- point coordinates are follow same pattern of replacing T with Cu, Cv and Bu respectively.
658 !!
659 !! Each location also has a 2D mask indicating whether the entire column is land or ocean.
660 !! `mask2dT` is 1 if the column is wet or 0 if the T-cell is land.
661 !! `mask2dCu` is 1 if both neighboring column are ocean, and 0 if either is land.
662 
663 end module mom_grid
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_unit_scaling::unit_scale_type
Describes various unit conversion factors.
Definition: MOM_unit_scaling.F90:14
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:16
mom_grid
Provides the ocean grid type.
Definition: MOM_grid.F90:2
mom_unit_scaling
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Definition: MOM_unit_scaling.F90:2
mom_domains::mom_domain_type
The MOM_domain_type contains information about the domain decompositoin.
Definition: MOM_domains.F90:104
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
mom_grid::ocean_grid_type
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:26