MOM6
MOM_grid_initialize.F90
1 !> Initializes horizontal grid
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 use mom_checksums, only : hchksum, bchksum
8 use mom_domains, only : pass_var, pass_vector, pe_here, root_pe, broadcast
9 use mom_domains, only : agrid, bgrid_ne, cgrid_ne, to_all, scalar_pair
10 use mom_domains, only : to_north, to_south, to_east, to_west
11 use mom_domains, only : mom_define_domain, mom_define_io_domain
12 use mom_domains, only : mom_domain_type
13 use mom_dyn_horgrid, only : dyn_horgrid_type, set_derived_dyn_horgrid
14 use mom_error_handler, only : mom_error, mom_mesg, fatal, is_root_pe
15 use mom_error_handler, only : calltree_enter, calltree_leave
17 use mom_io, only : mom_read_data, read_data, slasher, file_exists
18 use mom_io, only : corner, north_face, east_face
20 
21 use mpp_domains_mod, only : mpp_get_domain_extents, mpp_deallocate_domain
22 
23 implicit none ; private
24 
25 public set_grid_metrics, initialize_masks, adcroft_reciprocal
26 
27 ! A note on unit descriptions in comments: MOM6 uses units that can be rescaled for dimensional
28 ! consistency testing. These are noted in comments with units like Z, H, L, and T, along with
29 ! their mks counterparts with notation like "a velocity [Z T-1 ~> m s-1]". If the units
30 ! vary with the Boussinesq approximation, the Boussinesq variant is given first.
31 
32 !> Global positioning system (aka container for information to describe the grid)
33 type, public :: gps ; private
34  real :: len_lon !< The longitudinal or x-direction length of the domain.
35  real :: len_lat !< The latitudinal or y-direction length of the domain.
36  real :: west_lon !< The western longitude of the domain or the equivalent
37  !! starting value for the x-axis.
38  real :: south_lat !< The southern latitude of the domain or the equivalent
39  !! starting value for the y-axis.
40  real :: rad_earth !< The radius of the Earth [m].
41  real :: lat_enhance_factor !< The amount by which the meridional resolution
42  !! is enhanced within LAT_EQ_ENHANCE of the equator.
43  real :: lat_eq_enhance !< The latitude range to the north and south of the equator
44  !! over which the resolution is enhanced, in degrees.
45  logical :: isotropic !< If true, an isotropic grid on a sphere (also known as a Mercator grid)
46  !! is used. With an isotropic grid, the meridional extent of the domain
47  !! (LENLAT), the zonal extent (LENLON), and the number of grid points in each
48  !! direction are _not_ independent. In MOM the meridional extent is determined
49  !! to fit the zonal extent and the number of grid points, while grid is
50  !! perfectly isotropic.
51  logical :: equator_reference !< If true, the grid is defined to have the equator at the
52  !! nearest q- or h- grid point to (-LOWLAT*NJGLOBAL/LENLAT).
53  integer :: niglobal !< The number of i-points in the global grid computational domain
54  integer :: njglobal !< The number of j-points in the global grid computational domain
55 end type gps
56 
57 contains
58 
59 !> set_grid_metrics is used to set the primary values in the model's horizontal
60 !! grid. The bathymetry, land-sea mask and any restricted channel widths are
61 !! not known yet, so these are set later.
62 subroutine set_grid_metrics(G, param_file, US)
63  type(dyn_horgrid_type), intent(inout) :: g !< The dynamic horizontal grid type
64  type(param_file_type), intent(in) :: param_file !< Parameter file structure
65  type(unit_scale_type), optional, intent(in) :: us !< A dimensional unit scaling type
66  ! Local variables
67 ! This include declares and sets the variable "version".
68 #include "version_variable.h"
69  logical :: debug
70  character(len=256) :: config
71 
72  call calltree_enter("set_grid_metrics(), MOM_grid_initialize.F90")
73  call log_version(param_file, "MOM_grid_init", version, "")
74  call get_param(param_file, "MOM_grid_init", "GRID_CONFIG", config, &
75  "A character string that determines the method for "//&
76  "defining the horizontal grid. Current options are: \n"//&
77  " \t mosaic - read the grid from a mosaic (supergrid) \n"//&
78  " \t file set by GRID_FILE.\n"//&
79  " \t cartesian - use a (flat) Cartesian grid.\n"//&
80  " \t spherical - use a simple spherical grid.\n"//&
81  " \t mercator - use a Mercator spherical grid.", &
82  fail_if_missing=.true.)
83  call get_param(param_file, "MOM_grid_init", "DEBUG", debug, &
84  "If true, write out verbose debugging data.", &
85  default=.false., debuggingparam=.true.)
86 
87  ! These are defaults that may be changed in the next select block.
88  g%x_axis_units = "degrees_east" ; g%y_axis_units = "degrees_north"
89  select case (trim(config))
90  case ("mosaic"); call set_grid_metrics_from_mosaic(g, param_file)
91  case ("cartesian"); call set_grid_metrics_cartesian(g, param_file)
92  case ("spherical"); call set_grid_metrics_spherical(g, param_file)
93  case ("mercator"); call set_grid_metrics_mercator(g, param_file)
94  case ("file"); call mom_error(fatal, "MOM_grid_init: set_grid_metrics "//&
95  'GRID_CONFIG "file" is no longer a supported option. Use a '//&
96  'mosaic file ("mosaic") or one of the analytic forms instead.')
97  case default ; call mom_error(fatal, "MOM_grid_init: set_grid_metrics "//&
98  "Unrecognized grid configuration "//trim(config))
99  end select
100 
101  ! Calculate derived metrics (i.e. reciprocals and products)
102  call calltree_enter("set_derived_metrics(), MOM_grid_initialize.F90")
103  call set_derived_dyn_horgrid(g)
104  call calltree_leave("set_derived_metrics()")
105 
106  if (debug) call grid_metrics_chksum('MOM_grid_init/set_grid_metrics',g)
107 
108  call calltree_leave("set_grid_metrics()")
109 end subroutine set_grid_metrics
110 
111 ! ------------------------------------------------------------------------------
112 
113 !> grid_metrics_chksum performs a set of checksums on metrics on the grid for
114 !! debugging.
115 subroutine grid_metrics_chksum(parent, G)
116  character(len=*), intent(in) :: parent !< A string identifying the caller
117  type(dyn_horgrid_type), intent(in) :: G !< The dynamic horizontal grid type
118 
119  integer :: halo
120 
121  halo = min(g%ied-g%iec, g%jed-g%jec, 1)
122 
123  call hchksum_pair(trim(parent)//': d[xy]T', &
124  g%dxT, g%dyT, g%HI, haloshift=halo)
125 
126  call uvchksum(trim(parent)//': dxC[uv]', g%dxCu, g%dyCv, g%HI, haloshift=halo)
127 
128  call uvchksum(trim(parent)//': dxC[uv]', &
129  g%dyCu, g%dxCv, g%HI, haloshift=halo)
130 
131  call bchksum_pair(trim(parent)//': dxB[uv]', &
132  g%dxBu, g%dyBu, g%HI, haloshift=halo)
133 
134  call hchksum_pair(trim(parent)//': Id[xy]T', &
135  g%IdxT, g%IdyT, g%HI, haloshift=halo)
136 
137  call uvchksum(trim(parent)//': Id[xy]C[uv]', &
138  g%IdxCu, g%IdyCv, g%HI, haloshift=halo)
139 
140  call uvchksum(trim(parent)//': Id[xy]C[uv]', &
141  g%IdyCu, g%IdxCv, g%HI, haloshift=halo)
142 
143  call bchksum_pair(trim(parent)//': Id[xy]B[uv]', &
144  g%IdxBu, g%IdyBu, g%HI, haloshift=halo)
145 
146  call hchksum(g%areaT, trim(parent)//': areaT',g%HI, haloshift=halo)
147  call bchksum(g%areaBu, trim(parent)//': areaBu',g%HI, haloshift=halo)
148 
149  call hchksum(g%IareaT, trim(parent)//': IareaT',g%HI, haloshift=halo)
150  call bchksum(g%IareaBu, trim(parent)//': IareaBu',g%HI, haloshift=halo)
151 
152  call hchksum(g%geoLonT,trim(parent)//': geoLonT',g%HI, haloshift=halo)
153  call hchksum(g%geoLatT,trim(parent)//': geoLatT',g%HI, haloshift=halo)
154 
155  call bchksum(g%geoLonBu, trim(parent)//': geoLonBu',g%HI, haloshift=halo)
156  call bchksum(g%geoLatBu, trim(parent)//': geoLatBu',g%HI, haloshift=halo)
157 
158  call uvchksum(trim(parent)//': geoLonC[uv]', &
159  g%geoLonCu, g%geoLonCv, g%HI, haloshift=halo)
160 
161  call uvchksum(trim(parent)//': geoLatC[uv]', &
162  g%geoLatCu, g%geoLatCv, g%HI, haloshift=halo)
163 
164 end subroutine grid_metrics_chksum
165 
166 ! ------------------------------------------------------------------------------
167 
168 !> Sets the grid metrics from a mosaic file.
169 subroutine set_grid_metrics_from_mosaic(G, param_file)
170  type(dyn_horgrid_type), intent(inout) :: G !< The dynamic horizontal grid type
171  type(param_file_type), intent(in) :: param_file !< Parameter file structure
172  ! Local variables
173  real, dimension(G%isd :G%ied ,G%jsd :G%jed ) :: tempH1, tempH2, tempH3, tempH4
174  real, dimension(G%IsdB:G%IedB,G%JsdB:G%JedB) :: tempQ1, tempQ2, tempQ3, tempQ4
175  real, dimension(G%IsdB:G%IedB,G%jsd :G%jed ) :: tempE1, tempE2
176  real, dimension(G%isd :G%ied ,G%JsdB:G%JedB) :: tempN1, tempN2
177  ! These arrays are a holdover from earlier code in which the arrays in G were
178  ! macros and may have had reduced dimensions.
179  real, dimension(G%isd :G%ied ,G%jsd :G%jed ) :: dxT, dyT, areaT
180  real, dimension(G%IsdB:G%IedB,G%jsd :G%jed ) :: dxCu, dyCu
181  real, dimension(G%isd :G%ied ,G%JsdB:G%JedB) :: dxCv, dyCv
182  real, dimension(G%IsdB:G%IedB,G%JsdB:G%JedB) :: dxBu, dyBu, areaBu
183  ! This are symmetric arrays, corresponding to the data in the mosaic file
184  real, dimension(2*G%isd-2:2*G%ied+1,2*G%jsd-2:2*G%jed+1) :: tmpT
185  real, dimension(2*G%isd-3:2*G%ied+1,2*G%jsd-2:2*G%jed+1) :: tmpU
186  real, dimension(2*G%isd-2:2*G%ied+1,2*G%jsd-3:2*G%jed+1) :: tmpV
187  real, dimension(2*G%isd-3:2*G%ied+1,2*G%jsd-3:2*G%jed+1) :: tmpZ
188  real, dimension(:,:), allocatable :: tmpGlbl
189  character(len=200) :: filename, grid_file, inputdir
190  character(len=64) :: mdl = "MOM_grid_init set_grid_metrics_from_mosaic"
191  integer :: err=0, ni, nj, global_indices(4)
192  type(mom_domain_type) :: SGdom ! Supergrid domain
193  logical :: lon_bug ! If true use an older buggy answer in the tripolar longitude.
194  integer :: i, j, i2, j2
195  integer :: npei,npej
196  integer, dimension(:), allocatable :: exni,exnj
197  integer :: start(4), nread(4)
198 
199  call calltree_enter("set_grid_metrics_from_mosaic(), MOM_grid_initialize.F90")
200 
201  call get_param(param_file, mdl, "GRID_FILE", grid_file, &
202  "Name of the file from which to read horizontal grid data.", &
203  fail_if_missing=.true.)
204  call get_param(param_file, mdl, "USE_TRIPOLAR_GEOLONB_BUG", lon_bug, &
205  "If true, use older code that incorrectly sets the longitude "//&
206  "in some points along the tripolar fold to be off by 360 degrees.", &
207  default=.true.)
208  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
209  inputdir = slasher(inputdir)
210  filename = trim(adjustl(inputdir)) // trim(adjustl(grid_file))
211  call log_param(param_file, mdl, "INPUTDIR/GRID_FILE", filename)
212  if (.not.file_exists(filename)) &
213  call mom_error(fatal," set_grid_metrics_from_mosaic: Unable to open "//&
214  trim(filename))
215 
216  ! Initialize everything to 0.
217  dxcu(:,:) = 0.0 ; dycu(:,:) = 0.0
218  dxcv(:,:) = 0.0 ; dycv(:,:) = 0.0
219  dxbu(:,:) = 0.0 ; dybu(:,:) = 0.0 ; areabu(:,:) = 0.0
220 
221  !<MISSING CODE TO READ REFINEMENT LEVEL>
222  ni = 2*(g%iec-g%isc+1) ! i size of supergrid
223  nj = 2*(g%jec-g%jsc+1) ! j size of supergrid
224 
225  ! Define a domain for the supergrid (SGdom)
226  npei = g%domain%layout(1) ; npej = g%domain%layout(2)
227  allocate(exni(npei)) ; allocate(exnj(npej))
228  call mpp_get_domain_extents(g%domain%mpp_domain, exni, exnj)
229  allocate(sgdom%mpp_domain)
230  sgdom%nihalo = 2*g%domain%nihalo+1
231  sgdom%njhalo = 2*g%domain%njhalo+1
232  sgdom%niglobal = 2*g%domain%niglobal
233  sgdom%njglobal = 2*g%domain%njglobal
234  sgdom%layout(:) = g%domain%layout(:)
235  sgdom%io_layout(:) = g%domain%io_layout(:)
236  global_indices(1) = 1+sgdom%nihalo
237  global_indices(2) = sgdom%niglobal+sgdom%nihalo
238  global_indices(3) = 1+sgdom%njhalo
239  global_indices(4) = sgdom%njglobal+sgdom%njhalo
240  exni(:) = 2*exni(:) ; exnj(:) = 2*exnj(:)
241  if (associated(g%domain%maskmap)) then
242  call mom_define_domain(global_indices, sgdom%layout, sgdom%mpp_domain, &
243  xflags=g%domain%X_FLAGS, yflags=g%domain%Y_FLAGS, &
244  xhalo=sgdom%nihalo, yhalo=sgdom%njhalo, &
245  xextent=exni,yextent=exnj, &
246  symmetry=.true., name="MOM_MOSAIC", maskmap=g%domain%maskmap)
247  else
248  call mom_define_domain(global_indices, sgdom%layout, sgdom%mpp_domain, &
249  xflags=g%domain%X_FLAGS, yflags=g%domain%Y_FLAGS, &
250  xhalo=sgdom%nihalo, yhalo=sgdom%njhalo, &
251  xextent=exni,yextent=exnj, &
252  symmetry=.true., name="MOM_MOSAIC")
253  endif
254 
255  call mom_define_io_domain(sgdom%mpp_domain, sgdom%io_layout)
256  deallocate(exni)
257  deallocate(exnj)
258 
259  ! Read X from the supergrid
260  tmpz(:,:) = 999.
261  call mom_read_data(filename, 'x', tmpz, sgdom, position=corner)
262 
263  if (lon_bug) then
264  call pass_var(tmpz, sgdom, position=corner)
265  else
266  call pass_var(tmpz, sgdom, position=corner, inner_halo=0)
267  endif
268  call extrapolate_metric(tmpz, 2*(g%jsc-g%jsd)+2, missing=999.)
269  do j=g%jsd,g%jed ; do i=g%isd,g%ied ; i2 = 2*i ; j2 = 2*j
270  g%geoLonT(i,j) = tmpz(i2-1,j2-1)
271  enddo ; enddo
272  do j=g%JsdB,g%JedB ; do i=g%IsdB,g%IedB ; i2 = 2*i ; j2 = 2*j
273  g%geoLonBu(i,j) = tmpz(i2,j2)
274  enddo ; enddo
275  do j=g%jsd,g%jed ; do i=g%IsdB,g%IedB ; i2 = 2*i ; j2 = 2*j
276  g%geoLonCu(i,j) = tmpz(i2,j2-1)
277  enddo ; enddo
278  do j=g%JsdB,g%JedB ; do i=g%isd,g%ied ; i2 = 2*i ; j2 = 2*j
279  g%geoLonCv(i,j) = tmpz(i2-1,j2)
280  enddo ; enddo
281  ! For some reason, this messes up the solution...
282  ! call pass_var(G%geoLonBu, G%domain, position=CORNER)
283 
284  ! Read Y from the supergrid
285  tmpz(:,:) = 999.
286  call mom_read_data(filename, 'y', tmpz, sgdom, position=corner)
287 
288  call pass_var(tmpz, sgdom, position=corner)
289  call extrapolate_metric(tmpz, 2*(g%jsc-g%jsd)+2, missing=999.)
290  do j=g%jsd,g%jed ; do i=g%isd,g%ied ; i2 = 2*i ; j2 = 2*j
291  g%geoLatT(i,j) = tmpz(i2-1,j2-1)
292  enddo ; enddo
293  do j=g%JsdB,g%JedB ; do i=g%IsdB,g%IedB ; i2 = 2*i ; j2 = 2*j
294  g%geoLatBu(i,j) = tmpz(i2,j2)
295  enddo ; enddo
296  do j=g%jsd,g%jed ; do i=g%IsdB,g%IedB ; i2 = 2*i ; j2 = 2*j
297  g%geoLatCu(i,j) = tmpz(i2,j2-1)
298  enddo ; enddo
299  do j=g%JsdB,g%JedB ; do i=g%isd,g%ied ; i2 = 2*i ; j2 = 2*j
300  g%geoLatCv(i,j) = tmpz(i2-1,j2)
301  enddo ; enddo
302 
303  ! Read DX,DY from the supergrid
304  tmpu(:,:) = 0. ; tmpv(:,:) = 0.
305  call mom_read_data(filename,'dx',tmpv,sgdom,position=north_face)
306  call mom_read_data(filename,'dy',tmpu,sgdom,position=east_face)
307  call pass_vector(tmpu, tmpv, sgdom, to_all+scalar_pair, cgrid_ne)
308  call extrapolate_metric(tmpv, 2*(g%jsc-g%jsd)+2, missing=0.)
309  call extrapolate_metric(tmpu, 2*(g%jsc-g%jsd)+2, missing=0.)
310 
311  do j=g%jsd,g%jed ; do i=g%isd,g%ied ; i2 = 2*i ; j2 = 2*j
312  dxt(i,j) = tmpv(i2-1,j2-1) + tmpv(i2,j2-1)
313  dyt(i,j) = tmpu(i2-1,j2-1) + tmpu(i2-1,j2)
314  enddo ; enddo
315 
316  do j=g%jsd,g%jed ; do i=g%IsdB,g%IedB ; i2 = 2*i ; j2 = 2*j
317  dxcu(i,j) = tmpv(i2,j2-1) + tmpv(i2+1,j2-1)
318  dycu(i,j) = tmpu(i2,j2-1) + tmpu(i2,j2)
319  enddo ; enddo
320 
321  do j=g%JsdB,g%JedB ; do i=g%isd,g%ied ; i2 = 2*i ; j2 = 2*j
322  dxcv(i,j) = tmpv(i2-1,j2) + tmpv(i2,j2)
323  dycv(i,j) = tmpu(i2-1,j2) + tmpu(i2-1,j2+1)
324  enddo ; enddo
325 
326  do j=g%JsdB,g%JedB ; do i=g%IsdB,g%IedB ; i2 = 2*i ; j2 = 2*j
327  dxbu(i,j) = tmpv(i2,j2) + tmpv(i2+1,j2)
328  dybu(i,j) = tmpu(i2,j2) + tmpu(i2,j2+1)
329  enddo ; enddo
330 
331  ! Read AREA from the supergrid
332  tmpt(:,:) = 0.
333  call mom_read_data(filename, 'area', tmpt, sgdom)
334  call pass_var(tmpt, sgdom)
335  call extrapolate_metric(tmpt, 2*(g%jsc-g%jsd)+2, missing=0.)
336 
337  do j=g%jsd,g%jed ; do i=g%isd,g%ied ; i2 = 2*i ; j2 = 2*j
338  areat(i,j) = (tmpt(i2-1,j2-1) + tmpt(i2,j2)) + &
339  (tmpt(i2-1,j2) + tmpt(i2,j2-1))
340  enddo ; enddo
341  do j=g%JsdB,g%JedB ; do i=g%IsdB,g%IedB ; i2 = 2*i ; j2 = 2*j
342  areabu(i,j) = (tmpt(i2,j2) + tmpt(i2+1,j2+1)) + &
343  (tmpt(i2,j2+1) + tmpt(i2+1,j2))
344  enddo ; enddo
345 
346  ni=sgdom%niglobal
347  nj=sgdom%njglobal
348  call mpp_deallocate_domain(sgdom%mpp_domain)
349  deallocate(sgdom%mpp_domain)
350 
351  call pass_vector(dycu, dxcv, g%Domain, to_all+scalar_pair, cgrid_ne)
352  call pass_vector(dxcu, dycv, g%Domain, to_all+scalar_pair, cgrid_ne)
353  call pass_vector(dxbu, dybu, g%Domain, to_all+scalar_pair, bgrid_ne)
354  call pass_var(areat, g%Domain)
355  call pass_var(areabu, g%Domain, position=corner)
356 
357  do i=g%isd,g%ied ; do j=g%jsd,g%jed
358  g%dxT(i,j) = dxt(i,j) ; g%dyT(i,j) = dyt(i,j) ; g%areaT(i,j) = areat(i,j)
359  enddo ; enddo
360  do i=g%IsdB,g%IedB ; do j=g%jsd,g%jed
361  g%dxCu(i,j) = dxcu(i,j) ; g%dyCu(i,j) = dycu(i,j)
362  enddo ; enddo
363  do i=g%isd,g%ied ; do j=g%JsdB,g%JedB
364  g%dxCv(i,j) = dxcv(i,j) ; g%dyCv(i,j) = dycv(i,j)
365  enddo ; enddo
366  do i=g%IsdB,g%IedB ; do j=g%JsdB,g%JedB
367  g%dxBu(i,j) = dxbu(i,j) ; g%dyBu(i,j) = dybu(i,j) ; g%areaBu(i,j) = areabu(i,j)
368  enddo ; enddo
369 
370  ! Construct axes for diagnostic output (only necessary because "ferret" uses
371  ! broken convention for interpretting netCDF files).
372  start(:) = 1 ; nread(:) = 1
373  start(2) = 2 ; nread(1) = ni+1 ; nread(2) = 2
374  allocate( tmpglbl(ni+1,2) )
375  if (is_root_pe()) &
376  call read_data(filename, "x", tmpglbl, start, nread, no_domain=.true.)
377  call broadcast(tmpglbl, 2*(ni+1), root_pe())
378 
379  ! I don't know why the second axis is 1 or 2 here. -RWH
380  do i=g%isg,g%ieg
381  g%gridLonT(i) = tmpglbl(2*(i-g%isg)+2,2)
382  enddo
383  ! Note that the dynamic grid always uses symmetric memory for the global
384  ! arrays G%gridLatB and G%gridLonB.
385  do i=g%isg-1,g%ieg
386  g%gridLonB(i) = tmpglbl(2*(i-g%isg)+3,1)
387  enddo
388  deallocate( tmpglbl )
389 
390  allocate( tmpglbl(1, nj+1) )
391  start(:) = 1 ; nread(:) = 1
392  start(1) = int(ni/4)+1 ; nread(2) = nj+1
393  if (is_root_pe()) &
394  call read_data(filename, "y", tmpglbl, start, nread, no_domain=.true.)
395  call broadcast(tmpglbl, nj+1, root_pe())
396 
397  do j=g%jsg,g%jeg
398  g%gridLatT(j) = tmpglbl(1,2*(j-g%jsg)+2)
399  enddo
400  do j=g%jsg-1,g%jeg
401  g%gridLatB(j) = tmpglbl(1,2*(j-g%jsg)+3)
402  enddo
403  deallocate( tmpglbl )
404 
405  call calltree_leave("set_grid_metrics_from_mosaic()")
406 end subroutine set_grid_metrics_from_mosaic
407 
408 
409 ! ------------------------------------------------------------------------------
410 
411 !> Calculate the values of the metric terms for a Cartesian grid that
412 !! might be used and save them in arrays.
413 !!
414 !! Within this subroutine, the x- and y- grid spacings and their
415 !! inverses and the cell areas centered on h, q, u, and v points are
416 !! calculated, as are the geographic locations of each of these 4
417 !! sets of points.
418 subroutine set_grid_metrics_cartesian(G, param_file)
419  type(dyn_horgrid_type), intent(inout) :: G !< The dynamic horizontal grid type
420  type(param_file_type), intent(in) :: param_file !< Parameter file structure
421  ! Local variables
422  integer :: i, j, isd, ied, jsd, jed, IsdB, IedB, JsdB, JedB, I1off, J1off
423  integer :: niglobal, njglobal
424  real :: grid_latT(G%jsd:G%jed), grid_latB(G%JsdB:G%JedB)
425  real :: grid_lonT(G%isd:G%ied), grid_lonB(G%IsdB:G%IedB)
426  real :: dx_everywhere, dy_everywhere ! Grid spacings in m.
427  real :: I_dx, I_dy ! Inverse grid spacings in m.
428  real :: PI
429  character(len=80) :: units_temp
430  character(len=48) :: mdl = "MOM_grid_init set_grid_metrics_cartesian"
431 
432  niglobal = g%Domain%niglobal ; njglobal = g%Domain%njglobal
433  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
434  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
435  i1off = g%idg_offset ; j1off = g%jdg_offset
436 
437  call calltree_enter("set_grid_metrics_cartesian(), MOM_grid_initialize.F90")
438 
439  pi = 4.0*atan(1.0)
440 
441  call get_param(param_file, mdl, "AXIS_UNITS", units_temp, &
442  "The units for the Cartesian axes. Valid entries are: \n"//&
443  " \t degrees - degrees of latitude and longitude \n"//&
444  " \t m - meters \n \t k - kilometers", default="degrees")
445  call get_param(param_file, mdl, "SOUTHLAT", g%south_lat, &
446  "The southern latitude of the domain or the equivalent "//&
447  "starting value for the y-axis.", units=units_temp, &
448  fail_if_missing=.true.)
449  call get_param(param_file, mdl, "LENLAT", g%len_lat, &
450  "The latitudinal or y-direction length of the domain.", &
451  units=units_temp, fail_if_missing=.true.)
452  call get_param(param_file, mdl, "WESTLON", g%west_lon, &
453  "The western longitude of the domain or the equivalent "//&
454  "starting value for the x-axis.", units=units_temp, &
455  default=0.0)
456  call get_param(param_file, mdl, "LENLON", g%len_lon, &
457  "The longitudinal or x-direction length of the domain.", &
458  units=units_temp, fail_if_missing=.true.)
459  call get_param(param_file, mdl, "RAD_EARTH", g%Rad_Earth, &
460  "The radius of the Earth.", units="m", default=6.378e6)
461 
462  if (units_temp(1:1) == 'k') then
463  g%x_axis_units = "kilometers" ; g%y_axis_units = "kilometers"
464  elseif (units_temp(1:1) == 'm') then
465  g%x_axis_units = "meters" ; g%y_axis_units = "meters"
466  endif
467  call log_param(param_file, mdl, "explicit AXIS_UNITS", g%x_axis_units)
468 
469  ! Note that the dynamic grid always uses symmetric memory for the global
470  ! arrays G%gridLatB and G%gridLonB.
471  do j=g%jsg-1,g%jeg
472  g%gridLatB(j) = g%south_lat + g%len_lat*real(j-(g%jsg-1))/real(njglobal)
473  enddo
474  do j=g%jsg,g%jeg
475  g%gridLatT(j) = g%south_lat + g%len_lat*(real(j-g%jsg)+0.5)/real(njglobal)
476  enddo
477  do i=g%isg-1,g%ieg
478  g%gridLonB(i) = g%west_lon + g%len_lon*real(i-(g%isg-1))/real(niglobal)
479  enddo
480  do i=g%isg,g%ieg
481  g%gridLonT(i) = g%west_lon + g%len_lon*(real(i-g%isg)+0.5)/real(niglobal)
482  enddo
483 
484  do j=jsdb,jedb
485  grid_latb(j) = g%south_lat + g%len_lat*real(j+j1off-(g%jsg-1))/real(njglobal)
486  enddo
487  do j=jsd,jed
488  grid_latt(j) = g%south_lat + g%len_lat*(real(j+j1off-g%jsg)+0.5)/real(njglobal)
489  enddo
490  do i=isdb,iedb
491  grid_lonb(i) = g%west_lon + g%len_lon*real(i+i1off-(g%isg-1))/real(niglobal)
492  enddo
493  do i=isd,ied
494  grid_lont(i) = g%west_lon + g%len_lon*(real(i+i1off-g%isg)+0.5)/real(niglobal)
495  enddo
496 
497  if (units_temp(1:1) == 'k') then ! Axes are measured in km.
498  dx_everywhere = 1000.0 * g%len_lon / (real(niglobal))
499  dy_everywhere = 1000.0 * g%len_lat / (real(njglobal))
500  elseif (units_temp(1:1) == 'm') then ! Axes are measured in m.
501  dx_everywhere = g%len_lon / (real(niglobal))
502  dy_everywhere = g%len_lat / (real(njglobal))
503  else ! Axes are measured in degrees of latitude and longitude.
504  dx_everywhere = g%Rad_Earth * g%len_lon * pi / (180.0 * niglobal)
505  dy_everywhere = g%Rad_Earth * g%len_lat * pi / (180.0 * njglobal)
506  endif
507 
508  i_dx = 1.0 / dx_everywhere ; i_dy = 1.0 / dy_everywhere
509 
510  do j=jsdb,jedb ; do i=isdb,iedb
511  g%geoLonBu(i,j) = grid_lonb(i) ; g%geoLatBu(i,j) = grid_latb(j)
512 
513  g%dxBu(i,j) = dx_everywhere ; g%IdxBu(i,j) = i_dx
514  g%dyBu(i,j) = dy_everywhere ; g%IdyBu(i,j) = i_dy
515  g%areaBu(i,j) = dx_everywhere * dy_everywhere ; g%IareaBu(i,j) = i_dx * i_dy
516  enddo ; enddo
517 
518  do j=jsd,jed ; do i=isd,ied
519  g%geoLonT(i,j) = grid_lont(i) ; g%geoLatT(i,j) = grid_latt(j)
520  g%dxT(i,j) = dx_everywhere ; g%IdxT(i,j) = i_dx
521  g%dyT(i,j) = dy_everywhere ; g%IdyT(i,j) = i_dy
522  g%areaT(i,j) = dx_everywhere * dy_everywhere ; g%IareaT(i,j) = i_dx * i_dy
523  enddo ; enddo
524 
525  do j=jsd,jed ; do i=isdb,iedb
526  g%geoLonCu(i,j) = grid_lonb(i) ; g%geoLatCu(i,j) = grid_latt(j)
527 
528  g%dxCu(i,j) = dx_everywhere ; g%IdxCu(i,j) = i_dx
529  g%dyCu(i,j) = dy_everywhere ; g%IdyCu(i,j) = i_dy
530  enddo ; enddo
531 
532  do j=jsdb,jedb ; do i=isd,ied
533  g%geoLonCv(i,j) = grid_lont(i) ; g%geoLatCv(i,j) = grid_latb(j)
534 
535  g%dxCv(i,j) = dx_everywhere ; g%IdxCv(i,j) = i_dx
536  g%dyCv(i,j) = dy_everywhere ; g%IdyCv(i,j) = i_dy
537  enddo ; enddo
538 
539  call calltree_leave("set_grid_metrics_cartesian()")
540 end subroutine set_grid_metrics_cartesian
541 
542 ! ------------------------------------------------------------------------------
543 
544 !> Calculate the values of the metric terms that might be used
545 !! and save them in arrays.
546 !!
547 !! Within this subroutine, the x- and y- grid spacings and their
548 !! inverses and the cell areas centered on h, q, u, and v points are
549 !! calculated, as are the geographic locations of each of these 4
550 !! sets of points.
551 subroutine set_grid_metrics_spherical(G, param_file)
552  type(dyn_horgrid_type), intent(inout) :: G !< The dynamic horizontal grid type
553  type(param_file_type), intent(in) :: param_file !< Parameter file structure
554  ! Local variables
555  real :: PI, PI_180! PI = 3.1415926... as 4*atan(1)
556  integer :: i, j, isd, ied, jsd, jed
557  integer :: is, ie, js, je, Isq, Ieq, Jsq, Jeq, IsdB, IedB, JsdB, JedB
558  integer :: i_offset, j_offset
559  real :: grid_latT(G%jsd:G%jed), grid_latB(G%JsdB:G%JedB)
560  real :: grid_lonT(G%isd:G%ied), grid_lonB(G%IsdB:G%IedB)
561  real :: dLon,dLat,latitude,longitude,dL_di
562  character(len=48) :: mdl = "MOM_grid_init set_grid_metrics_spherical"
563 
564  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
565  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
566  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
567  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
568  i_offset = g%idg_offset ; j_offset = g%jdg_offset
569 
570  call calltree_enter("set_grid_metrics_spherical(), MOM_grid_initialize.F90")
571 
572 ! Calculate the values of the metric terms that might be used
573 ! and save them in arrays.
574  pi = 4.0*atan(1.0); pi_180 = atan(1.0)/45.
575 
576  call get_param(param_file, mdl, "SOUTHLAT", g%south_lat, &
577  "The southern latitude of the domain.", units="degrees", &
578  fail_if_missing=.true.)
579  call get_param(param_file, mdl, "LENLAT", g%len_lat, &
580  "The latitudinal length of the domain.", units="degrees", &
581  fail_if_missing=.true.)
582  call get_param(param_file, mdl, "WESTLON", g%west_lon, &
583  "The western longitude of the domain.", units="degrees", &
584  default=0.0)
585  call get_param(param_file, mdl, "LENLON", g%len_lon, &
586  "The longitudinal length of the domain.", units="degrees", &
587  fail_if_missing=.true.)
588  call get_param(param_file, mdl, "RAD_EARTH", g%Rad_Earth, &
589  "The radius of the Earth.", units="m", default=6.378e6)
590 
591  dlon = g%len_lon/g%Domain%niglobal
592  dlat = g%len_lat/g%Domain%njglobal
593 
594  ! Note that the dynamic grid always uses symmetric memory for the global
595  ! arrays G%gridLatB and G%gridLonB.
596  do j=g%jsg-1,g%jeg
597  latitude = g%south_lat + dlat*(real(j-(g%jsg-1)))
598  g%gridLatB(j) = min(max(latitude,-90.),90.)
599  enddo
600  do j=g%jsg,g%jeg
601  latitude = g%south_lat + dlat*(real(j-g%jsg)+0.5)
602  g%gridLatT(j) = min(max(latitude,-90.),90.)
603  enddo
604  do i=g%isg-1,g%ieg
605  g%gridLonB(i) = g%west_lon + dlon*(real(i-(g%isg-1)))
606  enddo
607  do i=g%isg,g%ieg
608  g%gridLonT(i) = g%west_lon + dlon*(real(i-g%isg)+0.5)
609  enddo
610 
611  do j=jsdb,jedb
612  latitude = g%south_lat + dlat* real(j+j_offset-(g%jsg-1))
613  grid_latb(j) = min(max(latitude,-90.),90.)
614  enddo
615  do j=jsd,jed
616  latitude = g%south_lat + dlat*(real(j+j_offset-g%jsg)+0.5)
617  grid_latt(j) = min(max(latitude,-90.),90.)
618  enddo
619  do i=isdb,iedb
620  grid_lonb(i) = g%west_lon + dlon*real(i+i_offset-(g%isg-1))
621  enddo
622  do i=isd,ied
623  grid_lont(i) = g%west_lon + dlon*(real(i+i_offset-g%isg)+0.5)
624  enddo
625 
626  dl_di = (g%len_lon * 4.0*atan(1.0)) / (180.0 * g%Domain%niglobal)
627  do j=jsdb,jedb ; do i=isdb,iedb
628  g%geoLonBu(i,j) = grid_lonb(i)
629  g%geoLatBu(i,j) = grid_latb(j)
630 
631  ! The following line is needed to reproduce the solution from
632  ! set_grid_metrics_mercator when used to generate a simple spherical grid.
633  g%dxBu(i,j) = g%Rad_Earth * cos( g%geoLatBu(i,j)*pi_180 ) * dl_di
634 ! G%dxBu(I,J) = G%Rad_Earth * dLon*PI_180 * COS( G%geoLatBu(I,J)*PI_180 )
635  g%dyBu(i,j) = g%Rad_Earth * dlat*pi_180
636  g%areaBu(i,j) = g%dxBu(i,j) * g%dyBu(i,j)
637  enddo ; enddo
638 
639  do j=jsdb,jedb ; do i=isd,ied
640  g%geoLonCv(i,j) = grid_lont(i)
641  g%geoLatCv(i,j) = grid_latb(j)
642 
643  ! The following line is needed to reproduce the solution from
644  ! set_grid_metrics_mercator when used to generate a simple spherical grid.
645  g%dxCv(i,j) = g%Rad_Earth * cos( g%geoLatCv(i,j)*pi_180 ) * dl_di
646 ! G%dxCv(i,J) = G%Rad_Earth * (dLon*PI_180) * COS( G%geoLatCv(i,J)*PI_180 )
647  g%dyCv(i,j) = g%Rad_Earth * dlat*pi_180
648  enddo ; enddo
649 
650  do j=jsd,jed ; do i=isdb,iedb
651  g%geoLonCu(i,j) = grid_lonb(i)
652  g%geoLatCu(i,j) = grid_latt(j)
653 
654  ! The following line is needed to reproduce the solution from
655  ! set_grid_metrics_mercator when used to generate a simple spherical grid.
656  g%dxCu(i,j) = g%Rad_Earth * cos( g%geoLatCu(i,j)*pi_180 ) * dl_di
657 ! G%dxCu(I,j) = G%Rad_Earth * dLon*PI_180 * COS( latitude )
658  g%dyCu(i,j) = g%Rad_Earth * dlat*pi_180
659  enddo ; enddo
660 
661  do j=jsd,jed ; do i=isd,ied
662  g%geoLonT(i,j) = grid_lont(i)
663  g%geoLatT(i,j) = grid_latt(j)
664 
665  ! The following line is needed to reproduce the solution from
666  ! set_grid_metrics_mercator when used to generate a simple spherical grid.
667  g%dxT(i,j) = g%Rad_Earth * cos( g%geoLatT(i,j)*pi_180 ) * dl_di
668 ! G%dxT(i,j) = G%Rad_Earth * dLon*PI_180 * COS( latitude )
669  g%dyT(i,j) = g%Rad_Earth * dlat*pi_180
670 
671 ! latitude = G%geoLatCv(i,J)*PI_180 ! In radians
672 ! dL_di = G%geoLatCv(i,max(jsd,J-1))*PI_180 ! In radians
673 ! G%areaT(i,j) = Rad_Earth**2*dLon*dLat*ABS(SIN(latitude)-SIN(dL_di))
674  g%areaT(i,j) = g%dxT(i,j) * g%dyT(i,j)
675  enddo ; enddo
676 
677  call calltree_leave("set_grid_metrics_spherical()")
678 end subroutine set_grid_metrics_spherical
679 
680 !> Calculate the values of the metric terms that might be used
681 !! and save them in arrays.
682 !!
683 !! Within this subroutine, the x- and y- grid spacings and their
684 !! inverses and the cell areas centered on h, q, u, and v points are
685 !! calculated, as are the geographic locations of each of these 4
686 !! sets of points.
687 subroutine set_grid_metrics_mercator(G, param_file)
688  type(dyn_horgrid_type), intent(inout) :: G !< The dynamic horizontal grid type
689  type(param_file_type), intent(in) :: param_file !< Parameter file structure
690  ! Local variables
691  integer :: i, j, isd, ied, jsd, jed
692  integer :: I_off, J_off
693  type(gps) :: GP
694  character(len=128) :: warnmesg
695  character(len=48) :: mdl = "MOM_grid_init set_grid_metrics_mercator"
696  real :: PI, PI_2! PI = 3.1415926... as 4*atan(1), PI_2 = (PI) /2.0
697  real :: y_q, y_h, jd, x_q, x_h, id
698  real, dimension(G%isd:G%ied,G%jsd:G%jed) :: &
699  xh, yh ! Latitude and longitude of h points in radians.
700  real, dimension(G%IsdB:G%IedB,G%jsd:G%jed) :: &
701  xu, yu ! Latitude and longitude of u points in radians.
702  real, dimension(G%isd:G%ied,G%JsdB:G%JedB) :: &
703  xv, yv ! Latitude and longitude of v points in radians.
704  real, dimension(G%IsdB:G%IedB,G%JsdB:G%JedB) :: &
705  xq, yq ! Latitude and longitude of q points in radians.
706  real :: fnRef ! fnRef is the value of Int_dj_dy or
707  ! Int_dj_dy at a latitude or longitude that is
708  real :: jRef, iRef ! being set to be at grid index jRef or iRef.
709  integer :: itt1, itt2
710  logical :: debug = .false., simple_area = .true.
711  integer :: is, ie, js, je, Isq, Ieq, Jsq, Jeq, IsdB, IedB, JsdB, JedB
712 
713  ! All of the metric terms should be defined over the domain from
714  ! isd to ied. Outside of the physical domain, both the metrics
715  ! and their inverses may be set to zero.
716  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
717  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
718  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
719  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
720  i_off = g%idg_offset ; j_off = g%jdg_offset
721 
722  gp%niglobal = g%Domain%niglobal
723  gp%njglobal = g%Domain%njglobal
724 
725  call calltree_enter("set_grid_metrics_mercator(), MOM_grid_initialize.F90")
726 
727  ! Calculate the values of the metric terms that might be used
728  ! and save them in arrays.
729  pi = 4.0*atan(1.0) ; pi_2 = 0.5*pi
730 
731  call get_param(param_file, mdl, "SOUTHLAT", gp%south_lat, &
732  "The southern latitude of the domain.", units="degrees", &
733  fail_if_missing=.true.)
734  call get_param(param_file, mdl, "LENLAT", gp%len_lat, &
735  "The latitudinal length of the domain.", units="degrees", &
736  fail_if_missing=.true.)
737  call get_param(param_file, mdl, "WESTLON", gp%west_lon, &
738  "The western longitude of the domain.", units="degrees", &
739  default=0.0)
740  call get_param(param_file, mdl, "LENLON", gp%len_lon, &
741  "The longitudinal length of the domain.", units="degrees", &
742  fail_if_missing=.true.)
743  call get_param(param_file, mdl, "RAD_EARTH", gp%Rad_Earth, &
744  "The radius of the Earth.", units="m", default=6.378e6)
745  g%south_lat = gp%south_lat ; g%len_lat = gp%len_lat
746  g%west_lon = gp%west_lon ; g%len_lon = gp%len_lon
747  g%Rad_Earth = gp%Rad_Earth
748  call get_param(param_file, mdl, "ISOTROPIC", gp%isotropic, &
749  "If true, an isotropic grid on a sphere (also known as "//&
750  "a Mercator grid) is used. With an isotropic grid, the "//&
751  "meridional extent of the domain (LENLAT), the zonal "//&
752  "extent (LENLON), and the number of grid points in each "//&
753  "direction are _not_ independent. In MOM the meridional "//&
754  "extent is determined to fit the zonal extent and the "//&
755  "number of grid points, while grid is perfectly isotropic.", &
756  default=.false.)
757  call get_param(param_file, mdl, "EQUATOR_REFERENCE", gp%equator_reference, &
758  "If true, the grid is defined to have the equator at the "//&
759  "nearest q- or h- grid point to (-LOWLAT*NJGLOBAL/LENLAT).", &
760  default=.true.)
761  call get_param(param_file, mdl, "LAT_ENHANCE_FACTOR", gp%Lat_enhance_factor, &
762  "The amount by which the meridional resolution is "//&
763  "enhanced within LAT_EQ_ENHANCE of the equator.", &
764  units="nondim", default=1.0)
765  call get_param(param_file, mdl, "LAT_EQ_ENHANCE", gp%Lat_eq_enhance, &
766  "The latitude range to the north and south of the equator "//&
767  "over which the resolution is enhanced.", units="degrees", &
768  default=0.0)
769 
770  ! With an isotropic grid, the north-south extent of the domain,
771  ! the east-west extent, and the number of grid points in each
772  ! direction are _not_ independent. Here the north-south extent
773  ! will be determined to fit the east-west extent and the number of
774  ! grid points. The grid is perfectly isotropic.
775  if (gp%equator_reference) then
776  ! With the following expression, the equator will always be placed
777  ! on either h or q points, in a position consistent with the ratio
778  ! GP%south_lat to GP%len_lat.
779  jref = (g%jsg-1) + 0.5*floor(gp%njglobal*((-1.0*gp%south_lat*2.0)/gp%len_lat)+0.5)
780  fnref = int_dj_dy(0.0, gp)
781  else
782  ! The following line sets the reference latitude GP%south_lat at j=js-1 (or -2?)
783  jref = (g%jsg-1)
784  fnref = int_dj_dy((gp%south_lat*pi/180.0), gp)
785  endif
786 
787  ! These calculations no longer depend on the the order in which they
788  ! are performed because they all use the same (poor) starting guess and
789  ! iterate to convergence.
790  ! Note that the dynamic grid always uses symmetric memory for the global
791  ! arrays G%gridLatB and G%gridLonB.
792  do j=g%jsg-1,g%jeg
793  jd = fnref + (j - jref)
794  y_q = find_root(int_dj_dy, dy_dj, gp, jd, 0.0, -1.0*pi_2, pi_2, itt2)
795  g%gridLatB(j) = y_q*180.0/pi
796  ! if (is_root_pe()) &
797  ! write(*, '("J, y_q = ",I4,ES14.4," itts = ",I4)') j, y_q, itt2
798  enddo
799  do j=g%jsg,g%jeg
800  jd = fnref + (j - jref) - 0.5
801  y_h = find_root(int_dj_dy, dy_dj, gp, jd, 0.0, -1.0*pi_2, pi_2, itt1)
802  g%gridLatT(j) = y_h*180.0/pi
803  ! if (is_root_pe()) &
804  ! write(*, '("j, y_h = ",I4,ES14.4," itts = ",I4)') j, y_h, itt1
805  enddo
806  do j=jsdb+j_off,jedb+j_off
807  jd = fnref + (j - jref)
808  y_q = find_root(int_dj_dy, dy_dj, gp, jd, 0.0, -1.0*pi_2, pi_2, itt2)
809  do i=isdb,iedb ; yq(i,j-j_off) = y_q ; enddo
810  do i=isd,ied ; yv(i,j-j_off) = y_q ; enddo
811  enddo
812  do j=jsd+j_off,jed+j_off
813  jd = fnref + (j - jref) - 0.5
814  y_h = find_root(int_dj_dy, dy_dj, gp, jd, 0.0, -1.0*pi_2, pi_2, itt1)
815  if ((j >= jsd+j_off) .and. (j <= jed+j_off)) then
816  do i=isd,ied ; yh(i,j-j_off) = y_h ; enddo
817  do i=isdb,iedb ; yu(i,j-j_off) = y_h ; enddo
818  endif
819  enddo
820 
821  ! Determine the longitudes of the various points.
822 
823  ! These two lines place the western edge of the domain at GP%west_lon.
824  iref = (g%isg-1) + gp%niglobal
825  fnref = int_di_dx(((gp%west_lon+gp%len_lon)*pi/180.0), gp)
826 
827  ! These calculations no longer depend on the the order in which they
828  ! are performed because they all use the same (poor) starting guess and
829  ! iterate to convergence.
830  do i=g%isg-1,g%ieg
831  id = fnref + (i - iref)
832  x_q = find_root(int_di_dx, dx_di, gp, id, 0.0, -4.0*pi, 4.0*pi, itt2)
833  g%gridLonB(i) = x_q*180.0/pi
834  enddo
835  do i=g%isg,g%ieg
836  id = fnref + (i - iref) - 0.5
837  x_h = find_root(int_di_dx, dx_di, gp, id, 0.0, -4.0*pi, 4.0*pi, itt1)
838  g%gridLonT(i) = x_h*180.0/pi
839  enddo
840  do i=isdb+i_off,iedb+i_off
841  id = fnref + (i - iref)
842  x_q = find_root(int_di_dx, dx_di, gp, id, 0.0, -4.0*pi, 4.0*pi, itt2)
843  do j=jsdb,jedb ; xq(i-i_off,j) = x_q ; enddo
844  do j=jsd,jed ; xu(i-i_off,j) = x_q ; enddo
845  enddo
846  do i=isd+i_off,ied+i_off
847  id = fnref + (i - iref) - 0.5
848  x_h = find_root(int_di_dx, dx_di, gp, id, 0.0, -4.0*pi, 4.0*pi, itt1)
849  do j=jsd,jed ; xh(i-i_off,j) = x_h ; enddo
850  do j=jsdb,jedb ; xv(i-i_off,j) = x_h ; enddo
851  enddo
852 
853  do j=jsdb,jedb ; do i=isdb,iedb
854  g%geoLonBu(i,j) = xq(i,j)*180.0/pi
855  g%geoLatBu(i,j) = yq(i,j)*180.0/pi
856  g%dxBu(i,j) = ds_di(xq(i,j), yq(i,j), gp)
857  g%dyBu(i,j) = ds_dj(xq(i,j), yq(i,j), gp)
858 
859  g%areaBu(i,j) = g%dxBu(i,j) * g%dyBu(i,j)
860  g%IareaBu(i,j) = 1.0 / g%areaBu(i,j)
861  enddo ; enddo
862 
863  do j=jsd,jed ; do i=isd,ied
864  g%geoLonT(i,j) = xh(i,j)*180.0/pi
865  g%geoLatT(i,j) = yh(i,j)*180.0/pi
866  g%dxT(i,j) = ds_di(xh(i,j), yh(i,j), gp)
867  g%dyT(i,j) = ds_dj(xh(i,j), yh(i,j), gp)
868 
869  g%areaT(i,j) = g%dxT(i,j)*g%dyT(i,j)
870  g%IareaT(i,j) = 1.0 / g%areaT(i,j)
871  enddo ; enddo
872 
873  do j=jsd,jed ; do i=isdb,iedb
874  g%geoLonCu(i,j) = xu(i,j)*180.0/pi
875  g%geoLatCu(i,j) = yu(i,j)*180.0/pi
876  g%dxCu(i,j) = ds_di(xu(i,j), yu(i,j), gp)
877  g%dyCu(i,j) = ds_dj(xu(i,j), yu(i,j), gp)
878  enddo ; enddo
879 
880  do j=jsdb,jedb ; do i=isd,ied
881  g%geoLonCv(i,j) = xv(i,j)*180.0/pi
882  g%geoLatCv(i,j) = yv(i,j)*180.0/pi
883  g%dxCv(i,j) = ds_di(xv(i,j), yv(i,j), gp)
884  g%dyCv(i,j) = ds_dj(xv(i,j), yv(i,j), gp)
885  enddo ; enddo
886 
887  if (.not.simple_area) then
888  do j=jsdb+1,jed ; do i=isdb+1,ied
889  g%areaT(i,j) = gp%Rad_Earth**2 * &
890  (dl(xq(i-1,j-1),xq(i-1,j),yq(i-1,j-1),yq(i-1,j)) + &
891  (dl(xq(i-1,j),xq(i,j),yq(i-1,j),yq(i,j)) + &
892  (dl(xq(i,j),xq(i,j-1),yq(i,j),yq(i,j-1)) + &
893  dl(xq(i,j-1),xq(i-1,j-1),yq(i,j-1),yq(i-1,j-1)))))
894  enddo ;enddo
895  if ((isdb == isd) .or. (jsdb == jsq)) then
896  ! Fill in row and column 1 to calculate the area in the southernmost
897  ! and westernmost land cells when we are not using symmetric memory.
898  ! The pass_var call updates these values if they are not land cells.
899  g%areaT(isd+1,jsd) = g%areaT(isd+1,jsd+1)
900  do j=jsd,jed ; g%areaT(isd,j) = g%areaT(isd+1,j) ; enddo
901  do i=isd,ied ; g%areaT(i,jsd) = g%areaT(i,jsd+1) ; enddo
902  ! Now replace the data in the halos, if value values exist.
903  call pass_var(g%areaT,g%Domain)
904  endif
905  do j=jsd,jed ; do i=isd,ied
906  g%IareaT(i,j) = 1.0 / g%areaT(i,j)
907  enddo ; enddo
908  endif
909 
910  call calltree_leave("set_grid_metrics_mercator()")
911 end subroutine set_grid_metrics_mercator
912 
913 
914 !> This function returns the grid spacing in the logical x direction.
915 function ds_di(x, y, GP)
916  real, intent(in) :: x !< The longitude in question
917  real, intent(in) :: y !< The latitude in question
918  type(gps), intent(in) :: gp !< A structure of grid parameters
919  real :: ds_di
920  ! Local variables
921 
922  ds_di = gp%Rad_Earth * cos(y) * dx_di(x,gp)
923  ! In general, this might be...
924  ! ds_di = GP%Rad_Earth * sqrt( cos(y)*cos(y) * dx_di(x,y,GP)*dx_di(x,y,GP) + &
925  ! dy_di(x,y,GP)*dy_di(x,y,GP))
926 end function ds_di
927 
928 !> This function returns the grid spacing in the logical y direction.
929 function ds_dj(x, y, GP)
930  real, intent(in) :: x !< The longitude in question
931  real, intent(in) :: y !< The latitude in question
932  type(gps), intent(in) :: gp !< A structure of grid parameters
933  ! Local variables
934  real :: ds_dj
935 
936  ds_dj = gp%Rad_Earth * dy_dj(y,gp)
937  ! In general, this might be...
938  ! ds_dj = GP%Rad_Earth * sqrt( cos(y)*cos(y) * dx_dj(x,y,GP)*dx_dj(x,y,GP) + &
939  ! dy_dj(x,y,GP)*dy_dj(x,y,GP))
940 end function ds_dj
941 
942 !> This function returns the contribution from the line integral along one of the four sides of a
943 !! cell face to the area of a cell, assuming that the sides follow a linear path in latitude and
944 !! longitude (i.e., on a Mercator grid).
945 function dl(x1, x2, y1, y2)
946  real, intent(in) :: x1 !< Segment starting longitude, in degrees E.
947  real, intent(in) :: x2 !< Segment ending longitude, in degrees E.
948  real, intent(in) :: y1 !< Segment ending latitude, in degrees N.
949  real, intent(in) :: y2 !< Segment ending latitude, in degrees N.
950  ! Local variables
951  real :: dl
952  real :: r, dy
953 
954  dy = y2 - y1
955 
956  if (abs(dy) > 2.5e-8) then
957  r = ((1.0 - cos(dy))*cos(y1) + sin(dy)*sin(y1)) / dy
958  else
959  r = (0.5*dy*cos(y1) + sin(y1))
960  endif
961  dl = r * (x2 - x1)
962 
963 end function dl
964 
965 !> This subroutine finds and returns the value of y at which the monotonically increasing
966 !! function fn takes the value fnval, also returning in ittmax the number of iterations of
967 !! Newton's method that were used to polish the root.
968 function find_root( fn, dy_df, GP, fnval, y1, ymin, ymax, ittmax)
969  real :: find_root !< The value of y where fn(y) = fnval that will be returned
970  real, external :: fn !< The external function whose root is being sought
971  real, external :: dy_df !< The inverse of the derivative of that function
972  type(gps), intent(in) :: gp !< A structure of grid parameters
973  real, intent(in) :: fnval !< The value of fn being sought
974  real, intent(in) :: y1 !< A first guess for y
975  real, intent(in) :: ymin !< The minimum permitted value of y
976  real, intent(in) :: ymax !< The maximum permitted value of y
977  integer, intent(out) :: ittmax !< The number of iterations used to polish the root
978  ! Local variables
979  real :: y, y_next
980  real :: ybot, ytop, fnbot, fntop
981  integer :: itt
982  character(len=256) :: warnmesg
983 
984  real :: dy_dfn, dy, fny
985 
986 ! Bracket the root. Do not use the bounding values because the value at the
987 ! function at the bounds could be infinite, as is the case for the Mercator
988 ! grid recursion relation. (I.e., this is a search on an open interval.)
989  ybot = y1
990  fnbot = fn(ybot,gp) - fnval ; itt = 0
991  do while (fnbot > 0.0)
992  if ((ybot - 2.0*dy_df(ybot,gp)) < (0.5*(ybot+ymin))) then
993  ! Go twice as far as the secant method would normally go.
994  ybot = ybot - 2.0*dy_df(ybot,gp)
995  else ! But stay within the open interval!
996  ybot = 0.5*(ybot+ymin) ; itt = itt + 1
997  endif
998  fnbot = fn(ybot,gp) - fnval
999 
1000  if ((itt > 50) .and. (fnbot > 0.0)) then
1001  write(warnmesg, '("PE ",I2," unable to find bottom bound for grid function. &
1002  &x = ",ES10.4,", xmax = ",ES10.4,", fn = ",ES10.4,", dfn_dx = ",ES10.4,&
1003  &", seeking fn = ",ES10.4," - fn = ",ES10.4,".")') &
1004  pe_here(),ybot,ymin,fn(ybot,gp),dy_df(ybot,gp),fnval, fnbot
1005  call mom_error(fatal,warnmesg)
1006  endif
1007  enddo
1008 
1009  ytop = y1
1010  fntop = fn(ytop,gp) - fnval ; itt = 0
1011  do while (fntop < 0.0)
1012  if ((ytop + 2.0*dy_df(ytop,gp)) < (0.5*(ytop+ymax))) then
1013  ! Go twice as far as the secant method would normally go.
1014  ytop = ytop + 2.0*dy_df(ytop,gp)
1015  else ! But stay within the open interval!
1016  ytop = 0.5*(ytop+ymax) ; itt = itt + 1
1017  endif
1018  fntop = fn(ytop,gp) - fnval
1019 
1020  if ((itt > 50) .and. (fntop < 0.0)) then
1021  write(warnmesg, '("PE ",I2," unable to find top bound for grid function. &
1022  &x = ",ES10.4,", xmax = ",ES10.4,", fn = ",ES10.4,", dfn_dx = ",ES10.4, &
1023  &", seeking fn = ",ES10.4," - fn = ",ES10.4,".")') &
1024  pe_here(),ytop,ymax,fn(ytop,gp),dy_df(ytop,gp),fnval,fntop
1025  call mom_error(fatal,warnmesg)
1026  endif
1027  enddo
1028 
1029  ! Find the root using a bracketed variant of Newton's method, starting
1030  ! with a false-positon method first guess.
1031  if ((fntop < 0.0) .or. (fnbot > 0.0) .or. (ytop < ybot)) then
1032  write(warnmesg, '("PE ",I2," find_root failed to bracket function. y = ",&
1033  &2ES10.4,", fn = ",2ES10.4,".")') pe_here(),ybot,ytop,fnbot,fntop
1034  call mom_error(fatal, warnmesg)
1035  endif
1036 
1037  if (fntop == 0.0) then ; y = ytop ; fny = fntop
1038  elseif (fnbot == 0.0) then ; y = ybot ; fny = fnbot
1039  else
1040  y = (ybot*fntop - ytop*fnbot) / (fntop - fnbot)
1041  fny = fn(y,gp) - fnval
1042  if (fny < 0.0) then ; fnbot = fny ; ybot = y
1043  else ; fntop = fny ; ytop = y ; endif
1044  endif
1045 
1046  do itt=1,50
1047  dy_dfn = dy_df(y,gp)
1048 
1049  dy = -1.0* fny * dy_dfn
1050  y_next = y + dy
1051  if ((y_next >= ytop) .or. (y_next <= ybot)) then
1052  ! The Newton's method estimate has escaped bracketing, so use the
1053  ! false-position method instead. The complicated test is to properly
1054  ! handle the case where the iteration is down to roundoff level differences.
1055  y_next = y
1056  if (abs(fntop - fnbot) > epsilon(y) * (abs(fntop) + abs(fnbot))) &
1057  y_next = (ybot*fntop - ytop*fnbot) / (fntop - fnbot)
1058  endif
1059 
1060  dy = y_next - y
1061  if (abs(dy) < (2.0*epsilon(y)*(abs(y) + abs(y_next)) + 1.0e-20)) then
1062  y = y_next ; exit
1063  endif
1064  y = y_next
1065 
1066  fny = fn(y,gp) - fnval
1067  if (fny > 0.0) then ; ytop = y ; fntop = fny
1068  elseif (fny < 0.0) then ; ybot = y ; fnbot = fny
1069  else ; exit ; endif
1070 
1071  enddo
1072  if (abs(y) < 1e-12) y = 0.0
1073 
1074  ittmax = itt
1075  find_root = y
1076 end function find_root
1077 
1078 !> This function calculates and returns the value of dx/di, where x is the
1079 !! longitude in Radians, and i is the integral north-south grid index.
1080 function dx_di(x, GP)
1081  real, intent(in) :: x !< The longitude in question
1082  type(gps), intent(in) :: gp !< A structure of grid parameters
1083  real :: dx_di
1084 
1085  dx_di = (gp%len_lon * 4.0*atan(1.0)) / (180.0 * gp%niglobal)
1086 
1087 end function dx_di
1088 
1089 !> This function calculates and returns the integral of the inverse
1090 !! of dx/di to the point x, in radians.
1091 function int_di_dx(x, GP)
1092  real, intent(in) :: x !< The longitude in question
1093  type(gps), intent(in) :: gp !< A structure of grid parameters
1094  real :: int_di_dx
1095 
1096  int_di_dx = x * ((180.0 * gp%niglobal) / (gp%len_lon * 4.0*atan(1.0)))
1097 
1098 end function int_di_dx
1099 
1100 !> This subroutine calculates and returns the value of dy/dj, where y is the
1101 !! latitude in Radians, and j is the integral north-south grid index.
1102 function dy_dj(y, GP)
1103  real, intent(in) :: y !< The latitude in question
1104  type(gps), intent(in) :: gp !< A structure of grid parameters
1105  real :: dy_dj
1106  ! Local variables
1107  real :: pi ! 3.1415926... calculated as 4*atan(1)
1108  real :: c0 ! The constant that converts the nominal y-spacing in
1109  ! gridpoints to the nominal spacing in Radians.
1110  real :: y_eq_enhance ! The latitude in radians within which the resolution
1111  ! is enhanced.
1112  pi = 4.0*atan(1.0)
1113  if (gp%isotropic) then
1114  c0 = (gp%len_lon * pi) / (180.0 * gp%niglobal)
1115  y_eq_enhance = pi*abs(gp%lat_eq_enhance)/180.0
1116  if (abs(y) < y_eq_enhance) then
1117  dy_dj = c0 * (cos(y) / (1.0 + 0.5*cos(y) * (gp%lat_enhance_factor - 1.0) * &
1118  (1.0+cos(pi*y/y_eq_enhance)) ))
1119  else
1120  dy_dj = c0 * cos(y)
1121  endif
1122  else
1123  c0 = (gp%len_lat * pi) / (180.0 * gp%njglobal)
1124  dy_dj = c0
1125  endif
1126 
1127 end function dy_dj
1128 
1129 !> This subroutine calculates and returns the integral of the inverse
1130 !! of dy/dj to the point y, in radians.
1131 function int_dj_dy(y, GP)
1132  real, intent(in) :: y !< The latitude in question
1133  type(gps), intent(in) :: gp !< A structure of grid parameters
1134  real :: int_dj_dy
1135  ! Local variables
1136  real :: i_c0 = 0.0 ! The inverse of the constant that converts the
1137  ! nominal spacing in gridpoints to the nominal
1138  ! spacing in Radians.
1139  real :: pi ! 3.1415926... calculated as 4*atan(1)
1140  real :: y_eq_enhance ! The latitude in radians from
1141  ! from the equator within which the
1142  ! meridional grid spacing is enhanced by
1143  ! a factor of GP%lat_enhance_factor.
1144  real :: r
1145 
1146  pi = 4.0*atan(1.0)
1147  if (gp%isotropic) then
1148  i_c0 = (180.0 * gp%niglobal) / (gp%len_lon * pi)
1149  y_eq_enhance = pi*abs(gp%lat_eq_enhance)/180.0
1150 
1151  if (y >= 0.0) then
1152  r = i_c0 * log((1.0 + sin(y))/cos(y))
1153  else
1154  r = -1.0 * i_c0 * log((1.0 - sin(y))/cos(y))
1155  endif
1156 
1157  if (y >= y_eq_enhance) then
1158  r = r + i_c0*0.5*(gp%lat_enhance_factor - 1.0)*y_eq_enhance
1159  elseif (y <= -y_eq_enhance) then
1160  r = r - i_c0*0.5*(gp%lat_enhance_factor - 1.0)*y_eq_enhance
1161  else
1162  r = r + i_c0*0.5*(gp%lat_enhance_factor - 1.0) * &
1163  (y + (y_eq_enhance/pi)*sin(pi*y/y_eq_enhance))
1164  endif
1165  else
1166  i_c0 = (180.0 * gp%njglobal) / (gp%len_lat * pi)
1167  r = i_c0 * y
1168  endif
1169 
1170  int_dj_dy = r
1171 end function int_dj_dy
1172 
1173 !> Extrapolates missing metric data into all the halo regions.
1174 subroutine extrapolate_metric(var, jh, missing)
1175  real, dimension(:,:), intent(inout) :: var !< The array in which to fill in halos
1176  integer, intent(in) :: jh !< The size of the halos to be filled
1177  real, optional, intent(in) :: missing !< The missing data fill value, 0 by default.
1178  ! Local variables
1179  real :: badval
1180  integer :: i,j
1181 
1182  badval = 0.0 ; if (present(missing)) badval = missing
1183 
1184  ! Fill in southern halo by extrapolating from the computational domain
1185  do j=lbound(var,2)+jh,lbound(var,2),-1 ; do i=lbound(var,1),ubound(var,1)
1186  if (var(i,j)==badval) var(i,j) = 2.0*var(i,j+1)-var(i,j+2)
1187  enddo ; enddo
1188 
1189  ! Fill in northern halo by extrapolating from the computational domain
1190  do j=ubound(var,2)-jh,ubound(var,2) ; do i=lbound(var,1),ubound(var,1)
1191  if (var(i,j)==badval) var(i,j) = 2.0*var(i,j-1)-var(i,j-2)
1192  enddo ; enddo
1193 
1194  ! Fill in western halo by extrapolating from the computational domain
1195  do j=lbound(var,2),ubound(var,2) ; do i=lbound(var,1)+jh,lbound(var,1),-1
1196  if (var(i,j)==badval) var(i,j) = 2.0*var(i+1,j)-var(i+2,j)
1197  enddo ; enddo
1198 
1199  ! Fill in eastern halo by extrapolating from the computational domain
1200  do j=lbound(var,2),ubound(var,2) ; do i=ubound(var,1)-jh,ubound(var,1)
1201  if (var(i,j)==badval) var(i,j) = 2.0*var(i-1,j)-var(i-2,j)
1202  enddo ; enddo
1203 
1204 end subroutine extrapolate_metric
1205 
1206 !> This function implements Adcroft's rule for reciprocals, namely that
1207 !! Adcroft_Inv(x) = 1/x for |x|>0 or 0 for x=0.
1208 function adcroft_reciprocal(val) result(I_val)
1209  real, intent(in) :: val !< The value being inverted.
1210  real :: i_val !< The Adcroft reciprocal of val.
1211 
1212  i_val = 0.0
1213  if (val /= 0.0) i_val = 1.0/val
1214 end function adcroft_reciprocal
1215 
1216 !> Initializes the grid masks and any metrics that come with masks already applied.
1217 !!
1218 !! Initialize_masks sets mask2dT, mask2dCu, mask2dCv, and mask2dBu to mask out
1219 !! flow over any points which are shallower than Dmin and permit an
1220 !! appropriate treatment of the boundary conditions. mask2dCu and mask2dCv
1221 !! are 0.0 at any points adjacent to a land point. mask2dBu is 0.0 at
1222 !! any land or boundary point. For points in the interior, mask2dCu,
1223 !! mask2dCv, and mask2dBu are all 1.0.
1224 subroutine initialize_masks(G, PF, US)
1225  type(dyn_horgrid_type), intent(inout) :: g !< The dynamic horizontal grid type
1226  type(param_file_type), intent(in) :: pf !< Parameter file structure
1227  type(unit_scale_type), optional, intent(in) :: us !< A dimensional unit scaling type
1228  ! Local variables
1229  real :: m_to_z_scale ! A unit conversion factor from m to Z.
1230  real :: dmin ! The depth for masking in the same units as G%bathyT [Z ~> m].
1231  real :: min_depth ! The minimum ocean depth in the same units as G%bathyT [Z ~> m].
1232  real :: mask_depth ! The depth shallower than which to mask a point as land [Z ~> m].
1233  character(len=40) :: mdl = "MOM_grid_init initialize_masks"
1234  integer :: i, j
1235 
1236  call calltree_enter("initialize_masks(), MOM_grid_initialize.F90")
1237  m_to_z_scale = 1.0 ; if (present(us)) m_to_z_scale = us%m_to_Z
1238  call get_param(pf, mdl, "MINIMUM_DEPTH", min_depth, &
1239  "If MASKING_DEPTH is unspecified, then anything shallower than "//&
1240  "MINIMUM_DEPTH is assumed to be land and all fluxes are masked out. "//&
1241  "If MASKING_DEPTH is specified, then all depths shallower than "//&
1242  "MINIMUM_DEPTH but deeper than MASKING_DEPTH are rounded to MINIMUM_DEPTH.", &
1243  units="m", default=0.0, scale=m_to_z_scale)
1244  call get_param(pf, mdl, "MASKING_DEPTH", mask_depth, &
1245  "The depth below which to mask points as land points, for which all "//&
1246  "fluxes are zeroed out. MASKING_DEPTH is ignored if negative.", &
1247  units="m", default=-9999.0, scale=m_to_z_scale)
1248 
1249  dmin = min_depth
1250  if (mask_depth>=0.) dmin = mask_depth
1251 
1252  g%mask2dCu(:,:) = 0.0 ; g%mask2dCv(:,:) = 0.0 ; g%mask2dBu(:,:) = 0.0
1253 
1254  ! Construct the h-point or T-point mask
1255  do j=g%jsd,g%jed ; do i=g%isd,g%ied
1256  if (g%bathyT(i,j) <= dmin) then
1257  g%mask2dT(i,j) = 0.0
1258  else
1259  g%mask2dT(i,j) = 1.0
1260  endif
1261  enddo ; enddo
1262 
1263  do j=g%jsd,g%jed ; do i=g%isd,g%ied-1
1264  if ((g%bathyT(i,j) <= dmin) .or. (g%bathyT(i+1,j) <= dmin)) then
1265  g%mask2dCu(i,j) = 0.0
1266  else
1267  g%mask2dCu(i,j) = 1.0
1268  endif
1269  enddo ; enddo
1270 
1271  do j=g%jsd,g%jed-1 ; do i=g%isd,g%ied
1272  if ((g%bathyT(i,j) <= dmin) .or. (g%bathyT(i,j+1) <= dmin)) then
1273  g%mask2dCv(i,j) = 0.0
1274  else
1275  g%mask2dCv(i,j) = 1.0
1276  endif
1277  enddo ; enddo
1278 
1279  do j=g%jsd,g%jed-1 ; do i=g%isd,g%ied-1
1280  if ((g%bathyT(i+1,j) <= dmin) .or. (g%bathyT(i+1,j+1) <= dmin) .or. &
1281  (g%bathyT(i,j) <= dmin) .or. (g%bathyT(i,j+1) <= dmin)) then
1282  g%mask2dBu(i,j) = 0.0
1283  else
1284  g%mask2dBu(i,j) = 1.0
1285  endif
1286  enddo ; enddo
1287 
1288  call pass_var(g%mask2dBu, g%Domain, position=corner)
1289  call pass_vector(g%mask2dCu, g%mask2dCv, g%Domain, to_all+scalar_pair, cgrid_ne)
1290 
1291  do j=g%jsd,g%jed ; do i=g%IsdB,g%IedB
1292  g%dy_Cu(i,j) = g%mask2dCu(i,j) * g%dyCu(i,j)
1293  g%areaCu(i,j) = g%dxCu(i,j) * g%dy_Cu(i,j)
1294  g%IareaCu(i,j) = g%mask2dCu(i,j) * adcroft_reciprocal(g%areaCu(i,j))
1295  enddo ; enddo
1296 
1297  do j=g%JsdB,g%JedB ; do i=g%isd,g%ied
1298  g%dx_Cv(i,j) = g%mask2dCv(i,j) * g%dxCv(i,j)
1299  g%areaCv(i,j) = g%dyCv(i,j) * g%dx_Cv(i,j)
1300  g%IareaCv(i,j) = g%mask2dCv(i,j) * adcroft_reciprocal(g%areaCv(i,j))
1301  enddo ; enddo
1302 
1303  call calltree_leave("initialize_masks()")
1304 end subroutine initialize_masks
1305 
1306 !> \namespace mom_grid_initialize
1307 !!
1308 !! The metric terms have the form Dzp, IDzp, or DXDYp, where z can
1309 !! be X or Y, and p can be q, u, v, or h. z describes the direction
1310 !! of the metric, while p describes the location. IDzp is the
1311 !! inverse of Dzp, while DXDYp is the product of DXp and DYp except
1312 !! that areaT is calculated analytically from the latitudes and
1313 !! longitudes of the surrounding q points.
1314 !!
1315 !! On a sphere, a variety of grids can be implemented by defining
1316 !! analytic expressions for dx_di, dy_dj (where x and y are latitude
1317 !! and longitude, and i and j are grid indices) and the expressions
1318 !! for the integrals of their inverses in the four subroutines
1319 !! dy_dj, Int_dj_dy, dx_di, and Int_di_dx.
1320 !!
1321 !! initialize_masks sets up land masks based on the depth field.
1322 !! The one argument is the minimum ocean depth. Depths that are
1323 !! less than this are interpreted as land points.
1324 
1325 end module mom_grid_initialize
mom_checksums::bchksum
Checksums an array (2d or 3d) staggered at corner points.
Definition: MOM_checksums.F90:53
mom_file_parser::log_version
An overloaded interface to log version information about modules.
Definition: MOM_file_parser.F90:109
mom_dyn_horgrid
Contains a shareable dynamic type for describing horizontal grids and metric data and utilty routines...
Definition: MOM_dyn_horgrid.F90:3
mom_checksums::uvchksum
Checksums a pair velocity arrays (2d or 3d) staggered at C-grid locations.
Definition: MOM_checksums.F90:28
mom_domains::pass_var
Do a halo update on an array.
Definition: MOM_domains.F90:49
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_grid_initialize::gps
Global positioning system (aka container for information to describe the grid)
Definition: MOM_grid_initialize.F90:33
mom_io
This module contains I/O framework code.
Definition: MOM_io.F90:2
mom_unit_scaling::unit_scale_type
Describes various unit conversion factors.
Definition: MOM_unit_scaling.F90:14
mom_domains::pass_vector
Do a halo update on a pair of arrays representing the two components of a vector.
Definition: MOM_domains.F90:54
mom_checksums
Routines to calculate checksums of various array and vector types.
Definition: MOM_checksums.F90:2
mom_domains
Describes the decomposed MOM domain and has routines for communications across PEs.
Definition: MOM_domains.F90:2
mom_io::mom_read_data
Read a data field from a file.
Definition: MOM_io.F90:74
mom_file_parser
The MOM6 facility to parse input files for runtime parameters.
Definition: MOM_file_parser.F90:2
mom_checksums::bchksum_pair
Checksums a pair of arrays (2d or 3d) staggered at corner points.
Definition: MOM_checksums.F90:43
mom_unit_scaling
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Definition: MOM_unit_scaling.F90:2
mom_grid_initialize
Initializes horizontal grid.
Definition: MOM_grid_initialize.F90:2
mom_io::file_exists
Indicate whether a file exists, perhaps with domain decomposition.
Definition: MOM_io.F90:68
mom_checksums::hchksum_pair
Checksums a pair of arrays (2d or 3d) staggered at tracer points.
Definition: MOM_checksums.F90:23
mom_checksums::hchksum
Checksums an array (2d or 3d) staggered at tracer points.
Definition: MOM_checksums.F90:48
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
mom_dyn_horgrid::dyn_horgrid_type
Describes the horizontal ocean grid with only dynamic memory arrays.
Definition: MOM_dyn_horgrid.F90:22