MOM6
MOM_regridding.F90
1 !> Generates vertical grids as part of the ALE algorithm
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 use mom_error_handler, only : mom_error, fatal, warning
8 use mom_io, only : file_exists, field_exists, field_size, mom_read_data
9 use mom_io, only : slasher
11 use mom_variables, only : ocean_grid_type, thermo_var_ptrs
14 use mom_string_functions, only : uppercase, extractword, extract_integer, extract_real
15 
16 use mom_remapping, only : remapping_cs
18 use regrid_consts, only : coordinatemode, default_coordinate_mode
19 use regrid_consts, only : regridding_layer, regridding_zstar
20 use regrid_consts, only : regridding_rho, regridding_sigma
21 use regrid_consts, only : regridding_arbitrary, regridding_sigma_shelf_zstar
22 use regrid_consts, only : regridding_hycom1, regridding_slight, regridding_adaptive
23 use regrid_interp, only : interp_cs_type, set_interp_scheme, set_interp_extrap
24 
25 use coord_zlike, only : init_coord_zlike, zlike_cs, set_zlike_params, build_zstar_column, end_coord_zlike
26 use coord_sigma, only : init_coord_sigma, sigma_cs, set_sigma_params, build_sigma_column, end_coord_sigma
27 use coord_rho, only : init_coord_rho, rho_cs, set_rho_params, build_rho_column, end_coord_rho
28 use coord_rho, only : old_inflate_layers_1d
29 use coord_hycom, only : init_coord_hycom, hycom_cs, set_hycom_params, build_hycom1_column, end_coord_hycom
30 use coord_slight, only : init_coord_slight, slight_cs, set_slight_params, build_slight_column, end_coord_slight
31 use coord_adapt, only : init_coord_adapt, adapt_cs, set_adapt_params, build_adapt_column, end_coord_adapt
32 
33 use netcdf ! Used by check_grid_def()
34 
35 implicit none ; private
36 
37 #include <MOM_memory.h>
38 
39 ! A note on unit descriptions in comments: MOM6 uses units that can be rescaled for dimensional
40 ! consistency testing. These are noted in comments with units like Z, H, L, and T, along with
41 ! their mks counterparts with notation like "a velocity [Z T-1 ~> m s-1]". If the units
42 ! vary with the Boussinesq approximation, the Boussinesq variant is given first.
43 
44 !> Regridding control structure
45 type, public :: regridding_cs ; private
46 
47  !> This array is set by function setCoordinateResolution()
48  !! It contains the "resolution" or delta coordinate of the target
49  !! coordinate. It has the units of the target coordinate, e.g.
50  !! [Z ~> m] for z*, non-dimensional for sigma, etc.
51  real, dimension(:), allocatable :: coordinateresolution
52 
53  !> This is a scaling factor that restores coordinateResolution to values in
54  !! the natural units for output.
55  real :: coord_scale = 1.0
56 
57  !> This array is set by function set_target_densities()
58  !! This array is the nominal coordinate of interfaces and is the
59  !! running sum of coordinateResolution, in [R ~> kg m-3]. i.e.
60  !! target_density(k+1) = coordinateResolution(k) + coordinateResolution(k)
61  !! It is only used in "rho", "SLight" or "Hycom" mode.
62  real, dimension(:), allocatable :: target_density
63 
64  !> A flag to indicate that the target_density arrays has been filled with data.
65  logical :: target_density_set = .false.
66 
67  !> This array is set by function set_regrid_max_depths()
68  !! It specifies the maximum depth that every interface is allowed to take [H ~> m or kg m-2].
69  real, dimension(:), allocatable :: max_interface_depths
70 
71  !> This array is set by function set_regrid_max_thickness()
72  !! It specifies the maximum depth that every interface is allowed to take [H ~> m or kg m-2].
73  real, dimension(:), allocatable :: max_layer_thickness
74 
75  integer :: nk !< Number of layers/levels in generated grid
76 
77  !> Indicates which grid to use in the vertical (z*, sigma, target interface
78  !! densities)
79  integer :: regridding_scheme
80 
81  !> Interpolation control structure
82  type(interp_cs_type) :: interp_cs
83 
84  !> Minimum thickness allowed when building the new grid through regridding [H ~> m or kg m-2].
85  real :: min_thickness
86 
87  !> Reference pressure for potential density calculations [R L2 T-2 ~> Pa]
88  real :: ref_pressure = 2.e7
89 
90  !> Weight given to old coordinate when blending between new and old grids [nondim]
91  !! Used only below depth_of_time_filter_shallow, with a cubic variation
92  !! from zero to full effect between depth_of_time_filter_shallow and
93  !! depth_of_time_filter_deep.
94  real :: old_grid_weight = 0.
95 
96  !> Depth above which no time-filtering of grid is applied [H ~> m or kg m-2]
97  real :: depth_of_time_filter_shallow = 0.
98 
99  !> Depth below which time-filtering of grid is applied at full effect [H ~> m or kg m-2]
100  real :: depth_of_time_filter_deep = 0.
101 
102  !> Fraction (between 0 and 1) of compressibility to add to potential density
103  !! profiles when interpolating for target grid positions. [nondim]
104  real :: compressibility_fraction = 0.
105 
106  !> If true, each interface is given a maximum depth based on a rescaling of
107  !! the indexing of coordinateResolution.
108  logical :: set_maximum_depths = .false.
109 
110  !> A scaling factor (> 1) of the rate at which the coordinateResolution list
111  !! is traversed to set the minimum depth of interfaces.
112  real :: max_depth_index_scale = 2.0
113 
114  !> If true, integrate for interface positions from the top downward.
115  !! If false, integrate from the bottom upward, as does the rest of the model.
116  logical :: integrate_downward_for_e = .true.
117 
118  !> If true, use the order of arithmetic and expressions that recover the remapping answers from 2018.
119  !! If false, use more robust forms of the same remapping expressions.
120  logical :: remap_answers_2018 = .true.
121 
122  type(zlike_cs), pointer :: zlike_cs => null() !< Control structure for z-like coordinate generator
123  type(sigma_cs), pointer :: sigma_cs => null() !< Control structure for sigma coordinate generator
124  type(rho_cs), pointer :: rho_cs => null() !< Control structure for rho coordinate generator
125  type(hycom_cs), pointer :: hycom_cs => null() !< Control structure for hybrid coordinate generator
126  type(slight_cs), pointer :: slight_cs => null() !< Control structure for Slight-coordinate generator
127  type(adapt_cs), pointer :: adapt_cs => null() !< Control structure for adaptive coordinate generator
128 
129 end type
130 
131 ! The following routines are visible to the outside world
132 public initialize_regridding, end_regridding, regridding_main
133 public inflate_vanished_layers_old, check_remapping_grid, check_grid_column
134 public set_regrid_params, get_regrid_size
135 public uniformresolution, setcoordinateresolution
136 public build_rho_column
137 public set_target_densities_from_gv, set_target_densities
138 public set_regrid_max_depths, set_regrid_max_thickness
139 public getcoordinateresolution, getcoordinateinterfaces
140 public getcoordinateunits, getcoordinateshortname, getstaticthickness
141 public default_coordinate_mode
142 public get_zlike_cs, get_sigma_cs, get_rho_cs
143 
144 !> Documentation for coordinate options
145 character(len=*), parameter, public :: regriddingcoordinatemodedoc = &
146  " LAYER - Isopycnal or stacked shallow water layers\n"//&
147  " ZSTAR, Z* - stretched geopotential z*\n"//&
148  " SIGMA_SHELF_ZSTAR - stretched geopotential z* ignoring shelf\n"//&
149  " SIGMA - terrain following coordinates\n"//&
150  " RHO - continuous isopycnal\n"//&
151  " HYCOM1 - HyCOM-like hybrid coordinate\n"//&
152  " SLIGHT - stretched coordinates above continuous isopycnal\n"//&
153  " ADAPTIVE - optimize for smooth neutral density surfaces"
154 
155 !> Documentation for regridding interpolation schemes
156 character(len=*), parameter, public :: regriddinginterpschemedoc = &
157  " P1M_H2 (2nd-order accurate)\n"//&
158  " P1M_H4 (2nd-order accurate)\n"//&
159  " P1M_IH4 (2nd-order accurate)\n"//&
160  " PLM (2nd-order accurate)\n"//&
161  " PPM_H4 (3rd-order accurate)\n"//&
162  " PPM_IH4 (3rd-order accurate)\n"//&
163  " P3M_IH4IH3 (4th-order accurate)\n"//&
164  " P3M_IH6IH5 (4th-order accurate)\n"//&
165  " PQM_IH4IH3 (4th-order accurate)\n"//&
166  " PQM_IH6IH5 (5th-order accurate)"
167 
168 !> Default interpolation scheme
169 character(len=*), parameter, public :: regriddingdefaultinterpscheme = "P1M_H2"
170 !> Default mode for boundary extrapolation
171 logical, parameter, public :: regriddingdefaultboundaryextrapolation = .false.
172 !> Default minimum thickness for some coordinate generation modes
173 real, parameter, public :: regriddingdefaultminthickness = 1.e-3
174 
175 #undef __DO_SAFETY_CHECKS__
176 
177 contains
178 
179 !> Initialization and configures a regridding control structure based on customizable run-time parameters
180 subroutine initialize_regridding(CS, GV, US, max_depth, param_file, mdl, coord_mode, param_prefix, param_suffix)
181  type(regridding_cs), intent(inout) :: cs !< Regridding control structure
182  type(verticalgrid_type), intent(in) :: gv !< Ocean vertical grid structure
183  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
184  real, intent(in) :: max_depth !< The maximum depth of the ocean [Z ~> m].
185  type(param_file_type), intent(in) :: param_file !< Parameter file
186  character(len=*), intent(in) :: mdl !< Name of calling module.
187  character(len=*), intent(in) :: coord_mode !< Coordinate mode
188  character(len=*), intent(in) :: param_prefix !< String to prefix to parameter names.
189  !! If empty, causes main model parameters to be used.
190  character(len=*), intent(in) :: param_suffix !< String to append to parameter names.
191 
192  ! Local variables
193  integer :: ke ! Number of levels
194  character(len=80) :: string, string2, varname ! Temporary strings
195  character(len=40) :: coord_units, param_name, coord_res_param ! Temporary strings
196  character(len=200) :: inputdir, filename
197  character(len=320) :: message ! Temporary strings
198  character(len=12) :: expected_units ! Temporary strings
199  logical :: tmplogical, fix_haloclines, set_max, do_sum, main_parameters
200  logical :: coord_is_state_dependent, ierr
201  logical :: default_2018_answers, remap_answers_2018
202  real :: filt_len, strat_tol, index_scale, tmpreal, p_ref
203  real :: maximum_depth ! The maximum depth of the ocean [m] (not in Z).
204  real :: dz_fixed_sfc, rho_avg_depth, nlay_sfc_int
205  real :: adapttimeratio, adaptzoom, adaptzoomcoeff, adaptbuoycoeff, adaptalpha
206  real :: adaptdrho0 ! Reference density difference for stratification-dependent diffusion. [R ~> kg m-3]
207  integer :: nz_fixed_sfc, k, nzf(4)
208  real, dimension(:), allocatable :: dz ! Resolution (thickness) in units of coordinate, which may be [m]
209  ! or [Z ~> m] or [H ~> m or kg m-2] or [R ~> kg m-3] or other units.
210  real, dimension(:), allocatable :: h_max ! Maximum layer thicknesses [H ~> m or kg m-2]
211  real, dimension(:), allocatable :: z_max ! Maximum interface depths [H ~> m or kg m-2] or other
212  ! units depending on the coordinate
213  real, dimension(:), allocatable :: dz_max ! Thicknesses used to find maximum interface depths
214  ! [H ~> m or kg m-2] or other units
215  real, dimension(:), allocatable :: rho_target ! Target density used in HYBRID mode [kg m-3]
216  ! Thicknesses [m] that give level centers corresponding to table 2 of WOA09
217  real, dimension(40) :: woa09_dz = (/ 5., 10., 10., 15., 22.5, 25., 25., 25., &
218  37.5, 50., 50., 75., 100., 100., 100., 100., &
219  100., 100., 100., 100., 100., 100., 100., 175., &
220  250., 375., 500., 500., 500., 500., 500., 500., &
221  500., 500., 500., 500., 500., 500., 500., 500. /)
222 
223  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
224  inputdir = slasher(inputdir)
225 
226  main_parameters=.false.
227  if (len_trim(param_prefix)==0) main_parameters=.true.
228  if (main_parameters .and. len_trim(param_suffix)>0) call mom_error(fatal,trim(mdl)//', initialize_regridding: '// &
229  'Suffix provided without prefix for parameter names!')
230 
231  cs%nk = 0
232  cs%regridding_scheme = coordinatemode(coord_mode)
233  coord_is_state_dependent = state_dependent(coord_mode)
234  maximum_depth = us%Z_to_m*max_depth
235 
236  if (main_parameters) then
237  ! Read coordinate units parameter (main model = REGRIDDING_COORDINATE_UNITS)
238  call get_param(param_file, mdl, "REGRIDDING_COORDINATE_UNITS", coord_units, &
239  "Units of the regridding coordinate.", default=coordinateunits(coord_mode))
240  else
241  coord_units=coordinateunits(coord_mode)
242  endif
243 
244  if (coord_is_state_dependent) then
245  if (main_parameters) then
246  param_name = "INTERPOLATION_SCHEME"
247  string2 = regriddingdefaultinterpscheme
248  else
249  param_name = trim(param_prefix)//"_INTERP_SCHEME_"//trim(param_suffix)
250  string2 = 'PPM_H4' ! Default for diagnostics
251  endif
252  call get_param(param_file, mdl, "INTERPOLATION_SCHEME", string, &
253  "This sets the interpolation scheme to use to "//&
254  "determine the new grid. These parameters are "//&
255  "only relevant when REGRIDDING_COORDINATE_MODE is "//&
256  "set to a function of state. Otherwise, it is not "//&
257  "used. It can be one of the following schemes: "//&
258  trim(regriddinginterpschemedoc), default=trim(string2))
259  call set_regrid_params(cs, interp_scheme=string)
260 
261  call get_param(param_file, mdl, "DEFAULT_2018_ANSWERS", default_2018_answers, &
262  "This sets the default value for the various _2018_ANSWERS parameters.", &
263  default=.false.)
264  call get_param(param_file, mdl, "REMAPPING_2018_ANSWERS", remap_answers_2018, &
265  "If true, use the order of arithmetic and expressions that recover the "//&
266  "answers from the end of 2018. Otherwise, use updated and more robust "//&
267  "forms of the same expressions.", default=default_2018_answers)
268  call set_regrid_params(cs, remap_answers_2018=remap_answers_2018)
269  endif
270 
271  if (main_parameters .and. coord_is_state_dependent) then
272  call get_param(param_file, mdl, "BOUNDARY_EXTRAPOLATION", tmplogical, &
273  "When defined, a proper high-order reconstruction "//&
274  "scheme is used within boundary cells rather "//&
275  "than PCM. E.g., if PPM is used for remapping, a "//&
276  "PPM reconstruction will also be used within "//&
277  "boundary cells.", default=regriddingdefaultboundaryextrapolation)
278  call set_regrid_params(cs, boundary_extrapolation=tmplogical)
279  else
280  call set_regrid_params(cs, boundary_extrapolation=.false.)
281  endif
282 
283  ! Read coordinate configuration parameter (main model = ALE_COORDINATE_CONFIG)
284  if (main_parameters) then
285  param_name = "ALE_COORDINATE_CONFIG"
286  coord_res_param = "ALE_RESOLUTION"
287  string2 = 'UNIFORM'
288  else
289  param_name = trim(param_prefix)//"_DEF_"//trim(param_suffix)
290  coord_res_param = trim(param_prefix)//"_RES_"//trim(param_suffix)
291  string2 = 'UNIFORM'
292  if (maximum_depth>3000.) string2='WOA09' ! For convenience
293  endif
294  call get_param(param_file, mdl, param_name, string, &
295  "Determines how to specify the coordinate "//&
296  "resolution. Valid options are:\n"//&
297  " PARAM - use the vector-parameter "//trim(coord_res_param)//"\n"//&
298  " UNIFORM[:N] - uniformly distributed\n"//&
299  " FILE:string - read from a file. The string specifies\n"//&
300  " the filename and variable name, separated\n"//&
301  " by a comma or space, e.g. FILE:lev.nc,dz\n"//&
302  " or FILE:lev.nc,interfaces=zw\n"//&
303  " WOA09[:N] - the WOA09 vertical grid (approximately)\n"//&
304  " FNC1:string - FNC1:dz_min,H_total,power,precision\n"//&
305  " HYBRID:string - read from a file. The string specifies\n"//&
306  " the filename and two variable names, separated\n"//&
307  " by a comma or space, for sigma-2 and dz. e.g.\n"//&
308  " HYBRID:vgrid.nc,sigma2,dz",&
309  default=trim(string2))
310  message = "The distribution of vertical resolution for the target\n"//&
311  "grid used for Eulerian-like coordinates. For example,\n"//&
312  "in z-coordinate mode, the parameter is a list of level\n"//&
313  "thicknesses (in m). In sigma-coordinate mode, the list\n"//&
314  "is of non-dimensional fractions of the water column."
315  if (index(trim(string),'UNIFORM')==1) then
316  if (len_trim(string)==7) then
317  ke = gv%ke ! Use model nk by default
318  tmpreal = maximum_depth
319  elseif (index(trim(string),'UNIFORM:')==1 .and. len_trim(string)>8) then
320  ! Format is "UNIFORM:N" or "UNIFORM:N,dz"
321  ke = extract_integer(string(9:len_trim(string)),'',1)
322  tmpreal = extract_real(string(9:len_trim(string)),',',2,missing_value=maximum_depth)
323  else
324  call mom_error(fatal,trim(mdl)//', initialize_regridding: '// &
325  'Unable to interpret "'//trim(string)//'".')
326  endif
327  allocate(dz(ke))
328  dz(:) = uniformresolution(ke, coord_mode, tmpreal, &
329  us%R_to_kg_m3*(gv%Rlay(1) + 0.5*(gv%Rlay(1)-gv%Rlay(min(2,ke)))), &
330  us%R_to_kg_m3*(gv%Rlay(ke) + 0.5*(gv%Rlay(ke)-gv%Rlay(max(ke-1,1)))) )
331  if (main_parameters) call log_param(param_file, mdl, "!"//coord_res_param, dz, &
332  trim(message), units=trim(coord_units))
333  elseif (trim(string)=='PARAM') then
334  ! Read coordinate resolution (main model = ALE_RESOLUTION)
335  ke = gv%ke ! Use model nk by default
336  allocate(dz(ke))
337  call get_param(param_file, mdl, coord_res_param, dz, &
338  trim(message), units=trim(coord_units), fail_if_missing=.true.)
339  elseif (index(trim(string),'FILE:')==1) then
340  ! FILE:filename,var_name is assumed to be reading level thickness variables
341  ! FILE:filename,interfaces=var_name reads positions
342  if (string(6:6)=='.' .or. string(6:6)=='/') then
343  ! If we specified "FILE:./xyz" or "FILE:/xyz" then we have a relative or absolute path
344  filename = trim( extractword(trim(string(6:80)), 1) )
345  else
346  ! Otherwise assume we should look for the file in INPUTDIR
347  filename = trim(inputdir) // trim( extractword(trim(string(6:80)), 1) )
348  endif
349  if (.not. file_exists(filename)) call mom_error(fatal,trim(mdl)//", initialize_regridding: "// &
350  "Specified file not found: Looking for '"//trim(filename)//"' ("//trim(string)//")")
351 
352  varname = trim( extractword(trim(string(6:)), 2) )
353  if (len_trim(varname)==0) then
354  if (field_exists(filename,'dz')) then; varname = 'dz'
355  elseif (field_exists(filename,'dsigma')) then; varname = 'dsigma'
356  elseif (field_exists(filename,'ztest')) then; varname = 'ztest'
357  else ; call mom_error(fatal,trim(mdl)//", initialize_regridding: "// &
358  "Coordinate variable not specified and none could be guessed.")
359  endif
360  endif
361  ! This check fails when the variable is a dimension variable! -AJA
362  !if (.not. field_exists(fileName,trim(varName))) call MOM_error(FATAL,trim(mdl)//", initialize_regridding: "// &
363  ! "Specified field not found: Looking for '"//trim(varName)//"' ("//trim(string)//")")
364  if (cs%regridding_scheme == regridding_sigma) then
365  expected_units = 'nondim'
366  elseif (cs%regridding_scheme == regridding_rho) then
367  expected_units = 'kg m-3'
368  else
369  expected_units = 'meters'
370  endif
371  if (index(trim(varname),'interfaces=')==1) then
372  varname=trim(varname(12:))
373  call check_grid_def(filename, varname, expected_units, message, ierr)
374  if (ierr) call mom_error(fatal,trim(mdl)//", initialize_regridding: "//&
375  "Unsupported format in grid definition '"//trim(filename)//"'. Error message "//trim(message))
376  call field_size(trim(filename), trim(varname), nzf)
377  ke = nzf(1)-1
378  if (ke < 1) call mom_error(fatal, trim(mdl)//" initialize_regridding via Var "//trim(varname)//&
379  "in FILE "//trim(filename)//" requires at least 2 target interface values.")
380  if (cs%regridding_scheme == regridding_rho) then
381  allocate(rho_target(ke+1))
382  call mom_read_data(trim(filename), trim(varname), rho_target)
383  else
384  allocate(dz(ke))
385  allocate(z_max(ke+1))
386  call mom_read_data(trim(filename), trim(varname), z_max)
387  dz(:) = abs(z_max(1:ke) - z_max(2:ke+1))
388  deallocate(z_max)
389  endif
390  else
391  ! Assume reading resolution
392  call field_size(trim(filename), trim(varname), nzf)
393  ke = nzf(1)
394  allocate(dz(ke))
395  call mom_read_data(trim(filename), trim(varname), dz)
396  endif
397  if (main_parameters .and. (ke/=gv%ke)) then
398  call mom_error(fatal,trim(mdl)//', initialize_regridding: '// &
399  'Mismatch in number of model levels and "'//trim(string)//'".')
400  endif
401  if (main_parameters) call log_param(param_file, mdl, "!"//coord_res_param, dz, &
402  trim(message), units=coordinateunits(coord_mode))
403  elseif (index(trim(string),'FNC1:')==1) then
404  ke = gv%ke; allocate(dz(ke))
405  call dz_function1( trim(string(6:)), dz )
406  if (main_parameters) call log_param(param_file, mdl, "!"//coord_res_param, dz, &
407  trim(message), units=coordinateunits(coord_mode))
408  elseif (index(trim(string),'RFNC1:')==1) then
409  ! Function used for set target interface densities
410  ke = rho_function1( trim(string(7:)), rho_target )
411  elseif (index(trim(string),'HYBRID:')==1) then
412  ke = gv%ke; allocate(dz(ke))
413  ! The following assumes the FILE: syntax of above but without "FILE:" in the string
414  allocate(rho_target(ke+1))
415  filename = trim( extractword(trim(string(8:)), 1) )
416  if (filename(1:1)/='.' .and. filename(1:1)/='/') filename = trim(inputdir) // trim( filename )
417  if (.not. file_exists(filename)) call mom_error(fatal,trim(mdl)//", initialize_regridding: HYBRID "// &
418  "Specified file not found: Looking for '"//trim(filename)//"' ("//trim(string)//")")
419  varname = trim( extractword(trim(string(8:)), 2) )
420  if (.not. field_exists(filename,varname)) call mom_error(fatal,trim(mdl)//", initialize_regridding: HYBRID "// &
421  "Specified field not found: Looking for '"//trim(varname)//"' ("//trim(string)//")")
422  call mom_read_data(trim(filename), trim(varname), rho_target)
423  varname = trim( extractword(trim(string(8:)), 3) )
424  if (varname(1:5) == 'FNC1:') then ! Use FNC1 to calculate dz
425  call dz_function1( trim(string((index(trim(string),'FNC1:')+5):)), dz )
426  else ! Read dz from file
427  if (.not. field_exists(filename,varname)) call mom_error(fatal,trim(mdl)//", initialize_regridding: HYBRID "// &
428  "Specified field not found: Looking for '"//trim(varname)//"' ("//trim(string)//")")
429  call mom_read_data(trim(filename), trim(varname), dz)
430  endif
431  if (main_parameters) then
432  call log_param(param_file, mdl, "!"//coord_res_param, dz, &
433  trim(message), units=coordinateunits(coord_mode))
434  call log_param(param_file, mdl, "!TARGET_DENSITIES", rho_target, &
435  'HYBRID target densities for interfaces', units=coordinateunits(coord_mode))
436  endif
437  elseif (index(trim(string),'WOA09')==1) then
438  if (len_trim(string)==5) then
439  tmpreal = 0. ; ke = 0
440  do while (tmpreal<maximum_depth)
441  ke = ke + 1
442  tmpreal = tmpreal + woa09_dz(ke)
443  enddo
444  elseif (index(trim(string),'WOA09:')==1) then
445  if (len_trim(string)==6) call mom_error(fatal,trim(mdl)//', initialize_regridding: '// &
446  'Expected string of form "WOA09:N" but got "'//trim(string)//'".')
447  ke = extract_integer(string(7:len_trim(string)),'',1)
448  endif
449  if (ke>40 .or. ke<1) call mom_error(fatal,trim(mdl)//', initialize_regridding: '// &
450  'For "WOA05:N" N must 0<N<41 but got "'//trim(string)//'".')
451  allocate(dz(ke))
452  dz(1:ke) = woa09_dz(1:ke)
453  else
454  call mom_error(fatal,trim(mdl)//", initialize_regridding: "// &
455  "Unrecognized coordinate configuration"//trim(string))
456  endif
457 
458  if (main_parameters) then
459  ! This is a work around to apparently needed to work with the from_Z initialization... ???
460  if (coordinatemode(coord_mode) == regridding_zstar .or. &
461  coordinatemode(coord_mode) == regridding_hycom1 .or. &
462  coordinatemode(coord_mode) == regridding_slight .or. &
463  coordinatemode(coord_mode) == regridding_adaptive) then
464  ! Adjust target grid to be consistent with maximum_depth
465  tmpreal = sum( dz(:) )
466  if (tmpreal < maximum_depth) then
467  dz(ke) = dz(ke) + ( maximum_depth - tmpreal )
468  elseif (tmpreal > maximum_depth) then
469  if ( dz(ke) + ( maximum_depth - tmpreal ) > 0. ) then
470  dz(ke) = dz(ke) + ( maximum_depth - tmpreal )
471  else
472  call mom_error(fatal,trim(mdl)//", initialize_regridding: "// &
473  "MAXIMUM_DEPTH was too shallow to adjust bottom layer of DZ!"//trim(string))
474  endif
475  endif
476  endif
477  endif
478 
479  cs%nk=ke
480 
481  ! Target resolution (for fixed coordinates)
482  allocate( cs%coordinateResolution(cs%nk) ); cs%coordinateResolution(:) = -1.e30
483  if (state_dependent(cs%regridding_scheme)) then
484  ! Target values
485  allocate( cs%target_density(cs%nk+1) ); cs%target_density(:) = -1.e30*us%kg_m3_to_R
486  endif
487 
488  if (allocated(dz)) then
489  if (coordinatemode(coord_mode) == regridding_sigma) then
490  call setcoordinateresolution(dz, cs, scale=1.0)
491  elseif (coordinatemode(coord_mode) == regridding_rho) then
492  call setcoordinateresolution(dz, cs, scale=us%kg_m3_to_R)
493  cs%coord_scale = us%R_to_kg_m3
494  elseif (coordinatemode(coord_mode) == regridding_adaptive) then
495  call setcoordinateresolution(dz, cs, scale=gv%m_to_H)
496  cs%coord_scale = gv%H_to_m
497  else
498  call setcoordinateresolution(dz, cs, scale=us%m_to_Z)
499  cs%coord_scale = us%Z_to_m
500  endif
501  endif
502 
503  if (allocated(rho_target)) then
504  call set_target_densities(cs, us%kg_m3_to_R*rho_target)
505  deallocate(rho_target)
506 
507  ! \todo This line looks like it would overwrite the target densities set just above?
508  elseif (coordinatemode(coord_mode) == regridding_rho) then
509  call set_target_densities_from_gv(gv, us, cs)
510  call log_param(param_file, mdl, "!TARGET_DENSITIES", us%R_to_kg_m3*cs%target_density(:), &
511  'RHO target densities for interfaces', units=coordinateunits(coord_mode))
512  endif
513 
514  ! initialise coordinate-specific control structure
515  call initcoord(cs, gv, us, coord_mode)
516 
517  if (main_parameters .and. coord_is_state_dependent) then
518  call get_param(param_file, mdl, "P_REF", p_ref, &
519  "The pressure that is used for calculating the coordinate "//&
520  "density. (1 Pa = 1e4 dbar, so 2e7 is commonly used.) "//&
521  "This is only used if USE_EOS and ENABLE_THERMODYNAMICS are true.", &
522  units="Pa", default=2.0e7, scale=us%kg_m3_to_R*us%m_s_to_L_T**2)
523  call get_param(param_file, mdl, "REGRID_COMPRESSIBILITY_FRACTION", tmpreal, &
524  "When interpolating potential density profiles we can add "//&
525  "some artificial compressibility solely to make homogeneous "//&
526  "regions appear stratified.", units="nondim", default=0.)
527  call set_regrid_params(cs, compress_fraction=tmpreal, ref_pressure=p_ref)
528  endif
529 
530  if (main_parameters) then
531  call get_param(param_file, mdl, "MIN_THICKNESS", tmpreal, &
532  "When regridding, this is the minimum layer "//&
533  "thickness allowed.", units="m", scale=gv%m_to_H, &
534  default=regriddingdefaultminthickness )
535  call set_regrid_params(cs, min_thickness=tmpreal)
536  else
537  call set_regrid_params(cs, min_thickness=0.)
538  endif
539 
540  if (coordinatemode(coord_mode) == regridding_slight) then
541  ! Set SLight-specific regridding parameters.
542  call get_param(param_file, mdl, "SLIGHT_DZ_SURFACE", dz_fixed_sfc, &
543  "The nominal thickness of fixed thickness near-surface "//&
544  "layers with the SLight coordinate.", units="m", default=1.0, scale=gv%m_to_H)
545  call get_param(param_file, mdl, "SLIGHT_NZ_SURFACE_FIXED", nz_fixed_sfc, &
546  "The number of fixed-depth surface layers with the SLight "//&
547  "coordinate.", units="nondimensional", default=2)
548  call get_param(param_file, mdl, "SLIGHT_SURFACE_AVG_DEPTH", rho_avg_depth, &
549  "The thickness of the surface region over which to average "//&
550  "when calculating the density to use to define the interior "//&
551  "with the SLight coordinate.", units="m", default=1.0, scale=gv%m_to_H)
552  call get_param(param_file, mdl, "SLIGHT_NLAY_TO_INTERIOR", nlay_sfc_int, &
553  "The number of layers to offset the surface density when "//&
554  "defining where the interior ocean starts with SLight.", &
555  units="nondimensional", default=2.0)
556  call get_param(param_file, mdl, "SLIGHT_FIX_HALOCLINES", fix_haloclines, &
557  "If true, identify regions above the reference pressure "//&
558  "where the reference pressure systematically underestimates "//&
559  "the stratification and use this in the definition of the "//&
560  "interior with the SLight coordinate.", default=.false.)
561 
562  call set_regrid_params(cs, dz_min_surface=dz_fixed_sfc, &
563  nz_fixed_surface=nz_fixed_sfc, rho_ml_avg_depth=rho_avg_depth, &
564  nlay_ml_to_interior=nlay_sfc_int, fix_haloclines=fix_haloclines)
565  if (fix_haloclines) then
566  ! Set additional parameters related to SLIGHT_FIX_HALOCLINES.
567  call get_param(param_file, mdl, "HALOCLINE_FILTER_LENGTH", filt_len, &
568  "A length scale over which to smooth the temperature and "//&
569  "salinity before identifying erroneously unstable haloclines.", &
570  units="m", default=2.0, scale=gv%m_to_H)
571  call get_param(param_file, mdl, "HALOCLINE_STRAT_TOL", strat_tol, &
572  "A tolerance for the ratio of the stratification of the "//&
573  "apparent coordinate stratification to the actual value "//&
574  "that is used to identify erroneously unstable haloclines. "//&
575  "This ratio is 1 when they are equal, and sensible values "//&
576  "are between 0 and 0.5.", units="nondimensional", default=0.2)
577  call set_regrid_params(cs, halocline_filt_len=filt_len, &
578  halocline_strat_tol=strat_tol)
579  endif
580 
581  endif
582 
583  if (coordinatemode(coord_mode) == regridding_adaptive) then
584  call get_param(param_file, mdl, "ADAPT_TIME_RATIO", adapttimeratio, &
585  "Ratio of ALE timestep to grid timescale.", units="nondim", default=1.0e-1)
586  call get_param(param_file, mdl, "ADAPT_ZOOM_DEPTH", adaptzoom, &
587  "Depth of near-surface zooming region.", units="m", default=200.0, scale=gv%m_to_H)
588  call get_param(param_file, mdl, "ADAPT_ZOOM_COEFF", adaptzoomcoeff, &
589  "Coefficient of near-surface zooming diffusivity.", units="nondim", default=0.2)
590  call get_param(param_file, mdl, "ADAPT_BUOY_COEFF", adaptbuoycoeff, &
591  "Coefficient of buoyancy diffusivity.", units="nondim", default=0.8)
592  call get_param(param_file, mdl, "ADAPT_ALPHA", adaptalpha, &
593  "Scaling on optimization tendency.", units="nondim", default=1.0)
594  call get_param(param_file, mdl, "ADAPT_DO_MIN_DEPTH", tmplogical, &
595  "If true, make a HyCOM-like mixed layer by preventing interfaces "//&
596  "from being shallower than the depths specified by the regridding coordinate.", &
597  default=.false.)
598  call get_param(param_file, mdl, "ADAPT_DRHO0", adaptdrho0, &
599  "Reference density difference for stratification-dependent diffusion.", &
600  units="kg m-3", default=0.5, scale=us%kg_m3_to_R)
601 
602  call set_regrid_params(cs, adapttimeratio=adapttimeratio, adaptzoom=adaptzoom, &
603  adaptzoomcoeff=adaptzoomcoeff, adaptbuoycoeff=adaptbuoycoeff, adaptalpha=adaptalpha, &
604  adaptdomin=tmplogical, adaptdrho0=adaptdrho0)
605  endif
606 
607  if (main_parameters .and. coord_is_state_dependent) then
608  call get_param(param_file, mdl, "MAXIMUM_INT_DEPTH_CONFIG", string, &
609  "Determines how to specify the maximum interface depths.\n"//&
610  "Valid options are:\n"//&
611  " NONE - there are no maximum interface depths\n"//&
612  " PARAM - use the vector-parameter MAXIMUM_INTERFACE_DEPTHS\n"//&
613  " FILE:string - read from a file. The string specifies\n"//&
614  " the filename and variable name, separated\n"//&
615  " by a comma or space, e.g. FILE:lev.nc,Z\n"//&
616  " FNC1:string - FNC1:dz_min,H_total,power,precision",&
617  default='NONE')
618  message = "The list of maximum depths for each interface."
619  allocate(z_max(ke+1))
620  allocate(dz_max(ke))
621  if ( trim(string) == "NONE") then
622  ! Do nothing.
623  elseif ( trim(string) == "PARAM") then
624  call get_param(param_file, mdl, "MAXIMUM_INTERFACE_DEPTHS", z_max, &
625  trim(message), units="m", scale=gv%m_to_H, fail_if_missing=.true.)
626  call set_regrid_max_depths(cs, z_max)
627  elseif (index(trim(string),'FILE:')==1) then
628  if (string(6:6)=='.' .or. string(6:6)=='/') then
629  ! If we specified "FILE:./xyz" or "FILE:/xyz" then we have a relative or absolute path
630  filename = trim( extractword(trim(string(6:80)), 1) )
631  else
632  ! Otherwise assume we should look for the file in INPUTDIR
633  filename = trim(inputdir) // trim( extractword(trim(string(6:80)), 1) )
634  endif
635  if (.not. file_exists(filename)) call mom_error(fatal,trim(mdl)//", initialize_regridding: "// &
636  "Specified file not found: Looking for '"//trim(filename)//"' ("//trim(string)//")")
637 
638  do_sum = .false.
639  varname = trim( extractword(trim(string(6:)), 2) )
640  if (.not. field_exists(filename,varname)) call mom_error(fatal,trim(mdl)//", initialize_regridding: "// &
641  "Specified field not found: Looking for '"//trim(varname)//"' ("//trim(string)//")")
642  if (len_trim(varname)==0) then
643  if (field_exists(filename,'z_max')) then; varname = 'z_max'
644  elseif (field_exists(filename,'dz')) then; varname = 'dz' ; do_sum = .true.
645  elseif (field_exists(filename,'dz_max')) then; varname = 'dz_max' ; do_sum = .true.
646  else ; call mom_error(fatal,trim(mdl)//", initialize_regridding: "// &
647  "MAXIMUM_INT_DEPTHS variable not specified and none could be guessed.")
648  endif
649  endif
650  if (do_sum) then
651  call mom_read_data(trim(filename), trim(varname), dz_max)
652  z_max(1) = 0.0 ; do k=1,ke ; z_max(k+1) = z_max(k) + dz_max(k) ; enddo
653  else
654  call mom_read_data(trim(filename), trim(varname), z_max)
655  endif
656  call log_param(param_file, mdl, "!MAXIMUM_INT_DEPTHS", z_max, &
657  trim(message), units=coordinateunits(coord_mode))
658  call set_regrid_max_depths(cs, z_max, gv%m_to_H)
659  elseif (index(trim(string),'FNC1:')==1) then
660  call dz_function1( trim(string(6:)), dz_max )
661  if ((coordinatemode(coord_mode) == regridding_slight) .and. &
662  (dz_fixed_sfc > 0.0)) then
663  do k=1,nz_fixed_sfc ; dz_max(k) = dz_fixed_sfc ; enddo
664  endif
665  z_max(1) = 0.0 ; do k=1,ke ; z_max(k+1) = z_max(k) + dz_max(k) ; enddo
666  call log_param(param_file, mdl, "!MAXIMUM_INT_DEPTHS", z_max, &
667  trim(message), units=coordinateunits(coord_mode))
668  call set_regrid_max_depths(cs, z_max, gv%m_to_H)
669  else
670  call mom_error(fatal,trim(mdl)//", initialize_regridding: "// &
671  "Unrecognized MAXIMUM_INT_DEPTH_CONFIG "//trim(string))
672  endif
673  deallocate(z_max)
674  deallocate(dz_max)
675 
676  ! Optionally specify maximum thicknesses for each layer, enforced by moving
677  ! the interface below a layer downward.
678  call get_param(param_file, mdl, "MAX_LAYER_THICKNESS_CONFIG", string, &
679  "Determines how to specify the maximum layer thicknesses.\n"//&
680  "Valid options are:\n"//&
681  " NONE - there are no maximum layer thicknesses\n"//&
682  " PARAM - use the vector-parameter MAX_LAYER_THICKNESS\n"//&
683  " FILE:string - read from a file. The string specifies\n"//&
684  " the filename and variable name, separated\n"//&
685  " by a comma or space, e.g. FILE:lev.nc,Z\n"//&
686  " FNC1:string - FNC1:dz_min,H_total,power,precision",&
687  default='NONE')
688  message = "The list of maximum thickness for each layer."
689  allocate(h_max(ke))
690  if ( trim(string) == "NONE") then
691  ! Do nothing.
692  elseif ( trim(string) == "PARAM") then
693  call get_param(param_file, mdl, "MAX_LAYER_THICKNESS", h_max, &
694  trim(message), units="m", fail_if_missing=.true., scale=gv%m_to_H)
695  call set_regrid_max_thickness(cs, h_max)
696  elseif (index(trim(string),'FILE:')==1) then
697  if (string(6:6)=='.' .or. string(6:6)=='/') then
698  ! If we specified "FILE:./xyz" or "FILE:/xyz" then we have a relative or absolute path
699  filename = trim( extractword(trim(string(6:80)), 1) )
700  else
701  ! Otherwise assume we should look for the file in INPUTDIR
702  filename = trim(inputdir) // trim( extractword(trim(string(6:80)), 1) )
703  endif
704  if (.not. file_exists(filename)) call mom_error(fatal,trim(mdl)//", initialize_regridding: "// &
705  "Specified file not found: Looking for '"//trim(filename)//"' ("//trim(string)//")")
706 
707  varname = trim( extractword(trim(string(6:)), 2) )
708  if (.not. field_exists(filename,varname)) call mom_error(fatal,trim(mdl)//", initialize_regridding: "// &
709  "Specified field not found: Looking for '"//trim(varname)//"' ("//trim(string)//")")
710  if (len_trim(varname)==0) then
711  if (field_exists(filename,'h_max')) then; varname = 'h_max'
712  elseif (field_exists(filename,'dz_max')) then; varname = 'dz_max'
713  else ; call mom_error(fatal,trim(mdl)//", initialize_regridding: "// &
714  "MAXIMUM_INT_DEPTHS variable not specified and none could be guessed.")
715  endif
716  endif
717  call mom_read_data(trim(filename), trim(varname), h_max)
718  call log_param(param_file, mdl, "!MAX_LAYER_THICKNESS", h_max, &
719  trim(message), units=coordinateunits(coord_mode))
720  call set_regrid_max_thickness(cs, h_max, gv%m_to_H)
721  elseif (index(trim(string),'FNC1:')==1) then
722  call dz_function1( trim(string(6:)), h_max )
723  call log_param(param_file, mdl, "!MAX_LAYER_THICKNESS", h_max, &
724  trim(message), units=coordinateunits(coord_mode))
725  call set_regrid_max_thickness(cs, h_max, gv%m_to_H)
726  else
727  call mom_error(fatal,trim(mdl)//", initialize_regridding: "// &
728  "Unrecognized MAX_LAYER_THICKNESS_CONFIG "//trim(string))
729  endif
730  deallocate(h_max)
731  endif
732 
733  if (allocated(dz)) deallocate(dz)
734 end subroutine initialize_regridding
735 
736 !> Do some basic checks on the vertical grid definition file, variable
737 subroutine check_grid_def(filename, varname, expected_units, msg, ierr)
738  character(len=*), intent(in) :: filename !< File name
739  character(len=*), intent(in) :: varname !< Variable name
740  character(len=*), intent(in) :: expected_units !< Expected units of variable
741  character(len=*), intent(inout) :: msg !< Message to use for errors
742  logical, intent(out) :: ierr !< True if an error occurs
743  ! Local variables
744  character (len=200) :: units, long_name
745  integer :: ncid, status, intid, vid
746  integer :: i
747 
748  ierr = .false.
749  status = nf90_open(trim(filename), nf90_nowrite, ncid)
750  if (status /= nf90_noerr) then
751  ierr = .true.
752  msg = 'File not found: '//trim(filename)
753  return
754  endif
755 
756  status = nf90_inq_varid(ncid, trim(varname), vid)
757  if (status /= nf90_noerr) then
758  ierr = .true.
759  msg = 'Var not found: '//trim(varname)
760  return
761  endif
762 
763  status = nf90_get_att(ncid, vid, "units", units)
764  if (status /= nf90_noerr) then
765  ierr = .true.
766  msg = 'Attribute not found: units'
767  return
768  endif
769  ! NF90_GET_ATT can return attributes with null characters, which TRIM will not truncate.
770  ! This loop replaces any null characters with a space so that the following check between
771  ! the read units and the expected units will pass
772  do i=1,len_trim(units)
773  if (units(i:i) == char(0)) units(i:i) = " "
774  enddo
775 
776  if (trim(units) /= trim(expected_units)) then
777  if (trim(expected_units) == "meters") then
778  if (trim(units) /= "m") then
779  ierr = .true.
780  endif
781  else
782  ierr = .true.
783  endif
784  endif
785 
786  if (ierr) then
787  msg = 'Units incorrect: '//trim(units)//' /= '//trim(expected_units)
788  endif
789 
790 end subroutine check_grid_def
791 
792 !> Deallocation of regridding memory
793 subroutine end_regridding(CS)
794  type(regridding_cs), intent(inout) :: cs !< Regridding control structure
795 
796  if (associated(cs%zlike_CS)) call end_coord_zlike(cs%zlike_CS)
797  if (associated(cs%sigma_CS)) call end_coord_sigma(cs%sigma_CS)
798  if (associated(cs%rho_CS)) call end_coord_rho(cs%rho_CS)
799  if (associated(cs%hycom_CS)) call end_coord_hycom(cs%hycom_CS)
800  if (associated(cs%slight_CS)) call end_coord_slight(cs%slight_CS)
801  if (associated(cs%adapt_CS)) call end_coord_adapt(cs%adapt_CS)
802 
803  deallocate( cs%coordinateResolution )
804  if (allocated(cs%target_density)) deallocate( cs%target_density )
805  if (allocated(cs%max_interface_depths) ) deallocate( cs%max_interface_depths )
806  if (allocated(cs%max_layer_thickness) ) deallocate( cs%max_layer_thickness )
807 
808 end subroutine end_regridding
809 
810 !------------------------------------------------------------------------------
811 !> Dispatching regridding routine for orchestrating regridding & remapping
812 subroutine regridding_main( remapCS, CS, G, GV, h, tv, h_new, dzInterface, frac_shelf_h, conv_adjust)
813 !------------------------------------------------------------------------------
814 ! This routine takes care of (1) building a new grid and (2) remapping between
815 ! the old grid and the new grid. The creation of the new grid can be based
816 ! on z coordinates, target interface densities, sigma coordinates or any
817 ! arbitrary coordinate system.
818 ! The MOM6 interface positions are always calculated from the bottom up by
819 ! accumulating the layer thicknesses starting at z=-G%bathyT. z increases
820 ! upwards (decreasing k-index).
821 ! The new grid is defined by the change in position of those interfaces in z
822 ! dzInterface = zNew - zOld.
823 ! Thus, if the regridding inflates the top layer, hNew(1) > hOld(1), then the
824 ! second interface moves downward, zNew(2) < zOld(2), and dzInterface(2) < 0.
825 ! hNew(k) = hOld(k) - dzInterface(k+1) + dzInterface(k)
826 ! IMPORTANT NOTE:
827 ! This is the converse of the sign convention used in the remapping code!
828 !------------------------------------------------------------------------------
829 
830  ! Arguments
831  type(remapping_cs), intent(in) :: remapcs !< Remapping parameters and options
832  type(regridding_cs), intent(in) :: cs !< Regridding control structure
833  type(ocean_grid_type), intent(in) :: g !< Ocean grid structure
834  type(verticalgrid_type), intent(in) :: gv !< Ocean vertical grid structure
835  real, dimension(SZI_(G),SZJ_(G), SZK_(GV)), intent(inout) :: h !< Current 3D grid obtained after
836  !! the last time step
837  type(thermo_var_ptrs), intent(inout) :: tv !< Thermodynamical variables (T, S, ...)
838  real, dimension(SZI_(G),SZJ_(G), CS%nk), intent(inout) :: h_new !< New 3D grid consistent with target coordinate
839  real, dimension(SZI_(G),SZJ_(G), CS%nk+1), intent(inout) :: dzinterface !< The change in position of each interface
840  real, dimension(:,:), optional, pointer :: frac_shelf_h !< Fractional ice shelf coverage
841  logical, optional, intent(in ) :: conv_adjust !< If true, do convective adjustment
842  ! Local variables
843  real :: trickgnucompiler
844  logical :: use_ice_shelf
845  logical :: do_convective_adjustment
846 
847  do_convective_adjustment = .true.
848  if (present(conv_adjust)) do_convective_adjustment = conv_adjust
849 
850  use_ice_shelf = .false.
851  if (present(frac_shelf_h)) then
852  if (associated(frac_shelf_h)) use_ice_shelf = .true.
853  endif
854 
855  select case ( cs%regridding_scheme )
856 
857  case ( regridding_zstar )
858  if (use_ice_shelf) then
859  call build_zstar_grid( cs, g, gv, h, dzinterface, frac_shelf_h )
860  else
861  call build_zstar_grid( cs, g, gv, h, dzinterface )
862  endif
863  call calc_h_new_by_dz(cs, g, gv, h, dzinterface, h_new)
864 
865  case ( regridding_sigma_shelf_zstar)
866  call build_zstar_grid( cs, g, gv, h, dzinterface )
867  call calc_h_new_by_dz(cs, g, gv, h, dzinterface, h_new)
868 
869  case ( regridding_sigma )
870  call build_sigma_grid( cs, g, gv, h, dzinterface )
871  call calc_h_new_by_dz(cs, g, gv, h, dzinterface, h_new)
872 
873  case ( regridding_rho )
874  if (do_convective_adjustment) call convective_adjustment(g, gv, h, tv)
875  call build_rho_grid( g, gv, g%US, h, tv, dzinterface, remapcs, cs )
876  call calc_h_new_by_dz(cs, g, gv, h, dzinterface, h_new)
877 
878  case ( regridding_arbitrary )
879  call build_grid_arbitrary( g, gv, h, dzinterface, trickgnucompiler, cs )
880  call calc_h_new_by_dz(cs, g, gv, h, dzinterface, h_new)
881 
882  case ( regridding_hycom1 )
883  call build_grid_hycom1( g, gv, g%US, h, tv, h_new, dzinterface, cs )
884 
885  case ( regridding_slight )
886  call build_grid_slight( g, gv, g%US, h, tv, dzinterface, cs )
887  call calc_h_new_by_dz(cs, g, gv, h, dzinterface, h_new)
888 
889  case ( regridding_adaptive )
890  call build_grid_adaptive(g, gv, g%US, h, tv, dzinterface, remapcs, cs)
891  call calc_h_new_by_dz(cs, g, gv, h, dzinterface, h_new)
892 
893  case default
894  call mom_error(fatal,'MOM_regridding, regridding_main: '//&
895  'Unknown regridding scheme selected!')
896 
897  end select ! type of grid
898 
899 #ifdef __DO_SAFETY_CHECKS__
900  call check_remapping_grid(g, gv, h, dzinterface,'in regridding_main')
901 #endif
902 
903 end subroutine regridding_main
904 
905 !> Calculates h_new from h + delta_k dzInterface
906 subroutine calc_h_new_by_dz(CS, G, GV, h, dzInterface, h_new)
907  type(regridding_cs), intent(in) :: CS !< Regridding control structure
908  type(ocean_grid_type), intent(in) :: G !< Grid structure
909  type(verticalgrid_type), intent(in) :: GV !< Ocean vertical grid structure
910  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: h !< Old layer thicknesses (arbitrary units)
911  real, dimension(SZI_(G),SZJ_(G),CS%nk+1), intent(in) :: dzInterface !< Change in interface positions (same as h)
912  real, dimension(SZI_(G),SZJ_(G),CS%nk), intent(inout) :: h_new !< New layer thicknesses (same as h)
913  ! Local variables
914  integer :: i, j, k, nki
915 
916  nki = min(cs%nk, gv%ke)
917 
918  !$OMP parallel do default(shared)
919  do j = g%jsc-1,g%jec+1
920  do i = g%isc-1,g%iec+1
921  if (g%mask2dT(i,j)>0.) then
922  do k=1,nki
923  h_new(i,j,k) = max( 0., h(i,j,k) + ( dzinterface(i,j,k) - dzinterface(i,j,k+1) ) )
924  enddo
925  if (cs%nk > gv%ke) then
926  do k=nki+1, cs%nk
927  h_new(i,j,k) = max( 0., dzinterface(i,j,k) - dzinterface(i,j,k+1) )
928  enddo
929  endif
930  else
931  h_new(i,j,1:nki) = h(i,j,1:nki)
932  if (cs%nk > gv%ke) h_new(i,j,nki+1:cs%nk) = 0.
933  ! On land points, why are we keeping the original h rather than setting to zero? -AJA
934  endif
935  enddo
936  enddo
937 
938 end subroutine calc_h_new_by_dz
939 
940 !> Check that the total thickness of two grids match
941 subroutine check_remapping_grid( G, GV, h, dzInterface, msg )
942  type(ocean_grid_type), intent(in) :: g !< Grid structure
943  type(verticalgrid_type), intent(in) :: gv !< Ocean vertical grid structure
944  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
945  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)+1), intent(in) :: dzinterface !< Change in interface positions
946  !! [H ~> m or kg m-2]
947  character(len=*), intent(in) :: msg !< Message to append to errors
948  ! Local variables
949  integer :: i, j
950 
951  !$OMP parallel do default(shared)
952  do j = g%jsc-1,g%jec+1 ; do i = g%isc-1,g%iec+1
953  if (g%mask2dT(i,j)>0.) call check_grid_column( gv%ke, gv%Z_to_H*g%bathyT(i,j), h(i,j,:), dzinterface(i,j,:), msg )
954  enddo ; enddo
955 
956 end subroutine check_remapping_grid
957 
958 !> Check that the total thickness of new and old grids are consistent
959 subroutine check_grid_column( nk, depth, h, dzInterface, msg )
960  integer, intent(in) :: nk !< Number of cells
961  real, intent(in) :: depth !< Depth of bottom [Z ~> m] or arbitrary units
962  real, dimension(nk), intent(in) :: h !< Cell thicknesses [Z ~> m] or arbitrary units
963  real, dimension(nk+1), intent(in) :: dzinterface !< Change in interface positions (same units as h)
964  character(len=*), intent(in) :: msg !< Message to append to errors
965  ! Local variables
966  integer :: k
967  real :: eps, total_h_old, total_h_new, h_new, z_old, z_new
968 
969  eps =1. ; eps = epsilon(eps)
970 
971  ! Total thickness of grid h
972  total_h_old = 0.
973  do k = 1,nk
974  total_h_old = total_h_old + h(k)
975  enddo
976 
977  ! Integrate upwards for the interfaces consistent with the rest of MOM6
978  z_old = - depth
979  if (depth == 0.) z_old = - total_h_old
980  total_h_new = 0.
981  do k = nk,1,-1
982  z_old = z_old + h(k) ! Old interface position above layer k
983  z_new = z_old + dzinterface(k) ! New interface position based on dzInterface
984  h_new = h(k) + ( dzinterface(k) - dzinterface(k+1) ) ! New thickness
985  if (h_new<0.) then
986  write(0,*) 'k,h,hnew=',k,h(k),h_new
987  write(0,*) 'dzI(k+1),dzI(k)=',dzinterface(k+1),dzinterface(k)
988  call mom_error( fatal, 'MOM_regridding, check_grid_column: '//&
989  'Negative layer thickness implied by re-gridding, '//trim(msg))
990  endif
991  total_h_new = total_h_new + h_new
992 
993  enddo
994 
995  ! Conservation by implied h_new
996  if (abs(total_h_new-total_h_old)>real(nk-1)*0.5*(total_h_old+total_h_new)*eps) then
997  write(0,*) 'nk=',nk
998  do k = 1,nk
999  write(0,*) 'k,h,hnew=',k,h(k),h(k)+(dzinterface(k)-dzinterface(k+1))
1000  enddo
1001  write(0,*) 'Hold,Hnew,Hnew-Hold=',total_h_old,total_h_new,total_h_new-total_h_old
1002  write(0,*) 'eps,(n)/2*eps*H=',eps,real(nk-1)*0.5*(total_h_old+total_h_new)*eps
1003  call mom_error( fatal, 'MOM_regridding, check_grid_column: '//&
1004  'Re-gridding did NOT conserve total thickness to within roundoff '//trim(msg))
1005  endif
1006 
1007  ! Check that the top and bottom are intentionally moving
1008  if (dzinterface(1) /= 0.) call mom_error( fatal, &
1009  'MOM_regridding, check_grid_column: Non-zero dzInterface at surface! '//trim(msg))
1010  if (dzinterface(nk+1) /= 0.) call mom_error( fatal, &
1011  'MOM_regridding, check_grid_column: Non-zero dzInterface at bottom! '//trim(msg))
1012 
1013 end subroutine check_grid_column
1014 
1015 !> Returns the change in interface position motion after filtering and
1016 !! assuming the top and bottom interfaces do not move. The filtering is
1017 !! a function of depth, and is applied as the integrated average filtering
1018 !! over the trajectory of the interface. By design, this code can not give
1019 !! tangled interfaces provided that z_old and z_new are not already tangled.
1020 subroutine filtered_grid_motion( CS, nk, z_old, z_new, dz_g )
1021  type(regridding_cs), intent(in) :: CS !< Regridding control structure
1022  integer, intent(in) :: nk !< Number of cells in source grid
1023  real, dimension(nk+1), intent(in) :: z_old !< Old grid position [H ~> m or kg m-2]
1024  real, dimension(CS%nk+1), intent(in) :: z_new !< New grid position [H ~> m or kg m-2]
1025  real, dimension(CS%nk+1), intent(inout) :: dz_g !< Change in interface positions [H ~> m or kg m-2]
1026  ! Local variables
1027  real :: sgn ! The sign convention for downward.
1028  real :: dz_tgt, zr1, z_old_k
1029  real :: Aq, Bq, dz0, z0, F0
1030  real :: zs, zd, dzwt, Idzwt
1031  real :: wtd, Iwtd
1032  real :: Int_zs, Int_zd, dInt_zs_zd
1033 ! For debugging:
1034  real, dimension(nk+1) :: z_act
1035 ! real, dimension(nk+1) :: ddz_g_s, ddz_g_d
1036  logical :: debug = .false.
1037  integer :: k
1038 
1039  if ((z_old(nk+1) - z_old(1)) * (z_new(cs%nk+1) - z_new(1)) < 0.0) then
1040  call mom_error(fatal, "filtered_grid_motion: z_old and z_new use different sign conventions.")
1041  elseif ((z_old(nk+1) - z_old(1)) * (z_new(cs%nk+1) - z_new(1)) == 0.0) then
1042  ! This is a massless column, so do nothing and return.
1043  do k=1,cs%nk+1 ; dz_g(k) = 0.0 ; enddo ; return
1044  elseif ((z_old(nk+1) - z_old(1)) + (z_new(cs%nk+1) - z_new(1)) > 0.0) then
1045  sgn = 1.0
1046  else
1047  sgn = -1.0
1048  endif
1049 
1050  if (debug) then
1051  do k=2,cs%nk+1
1052  if (sgn*(z_new(k)-z_new(k-1)) < -5e-16*(abs(z_new(k))+abs(z_new(k-1))) ) &
1053  call mom_error(fatal, "filtered_grid_motion: z_new is tangled.")
1054  enddo
1055  do k=2,nk+1
1056  if (sgn*(z_old(k)-z_old(k-1)) < -5e-16*(abs(z_old(k))+abs(z_old(k-1))) ) &
1057  call mom_error(fatal, "filtered_grid_motion: z_old is tangled.")
1058  enddo
1059  ! ddz_g_s(:) = 0.0 ; ddz_g_d(:) = 0.0
1060  endif
1061 
1062  zs = cs%depth_of_time_filter_shallow
1063  zd = cs%depth_of_time_filter_deep
1064  wtd = 1.0 - cs%old_grid_weight
1065  iwtd = 1.0 / wtd
1066 
1067  dzwt = (zd - zs)
1068  idzwt = 0.0 ; if (abs(zd - zs) > 0.0) idzwt = 1.0 / (zd - zs)
1069  dint_zs_zd = 0.5*(1.0 + iwtd) * (zd - zs)
1070  aq = 0.5*(iwtd - 1.0)
1071 
1072  dz_g(1) = 0.0
1073  z_old_k = z_old(1)
1074  do k = 2,cs%nk+1
1075  if (k<=nk+1) z_old_k = z_old(k) ! This allows for virtual z_old interface at bottom of the model
1076  ! zr1 is positive and increases with depth, and dz_tgt is positive downward.
1077  dz_tgt = sgn*(z_new(k) - z_old_k)
1078  zr1 = sgn*(z_old_k - z_old(1))
1079 
1080  ! First, handle the two simple and common cases that do not pass through
1081  ! the adjustment rate transition zone.
1082  if ((zr1 > zd) .and. (zr1 + wtd * dz_tgt > zd)) then
1083  dz_g(k) = sgn * wtd * dz_tgt
1084  elseif ((zr1 < zs) .and. (zr1 + dz_tgt < zs)) then
1085  dz_g(k) = sgn * dz_tgt
1086  else
1087  ! Find the new value by inverting the equation
1088  ! integral(0 to dz_new) Iwt(z) dz = dz_tgt
1089  ! This is trivial where Iwt is a constant, and agrees with the two limits above.
1090 
1091  ! Take test values at the transition points to figure out which segment
1092  ! the new value will be found in.
1093  if (zr1 >= zd) then
1094  int_zd = iwtd*(zd - zr1)
1095  int_zs = int_zd - dint_zs_zd
1096  elseif (zr1 <= zs) then
1097  int_zs = (zs - zr1)
1098  int_zd = dint_zs_zd + (zs - zr1)
1099  else
1100 ! Int_zd = (zd - zr1) * (Iwtd + 0.5*(1.0 - Iwtd) * (zd - zr1) / (zd - zs))
1101  int_zd = (zd - zr1) * (iwtd*(0.5*(zd+zr1) - zs) + 0.5*(zd - zr1)) * idzwt
1102  int_zs = (zs - zr1) * (0.5*iwtd * ((zr1 - zs)) + (zd - 0.5*(zr1+zs))) * idzwt
1103  ! It has been verified that Int_zs = Int_zd - dInt_zs_zd to within roundoff.
1104  endif
1105 
1106  if (dz_tgt >= int_zd) then ! The new location is in the deep, slow region.
1107  dz_g(k) = sgn * ((zd-zr1) + wtd*(dz_tgt - int_zd))
1108  elseif (dz_tgt <= int_zs) then ! The new location is in the shallow region.
1109  dz_g(k) = sgn * ((zs-zr1) + (dz_tgt - int_zs))
1110  else ! We need to solve a quadratic equation for z_new.
1111  ! For accuracy, do the integral from the starting depth or the nearest
1112  ! edge of the transition region. The results with each choice are
1113  ! mathematically equivalent, but differ in roundoff, and this choice
1114  ! should minimize the likelihood of inadvertently overlapping interfaces.
1115  if (zr1 <= zs) then ; dz0 = zs-zr1 ; z0 = zs ; f0 = dz_tgt - int_zs
1116  elseif (zr1 >= zd) then ; dz0 = zd-zr1 ; z0 = zd ; f0 = dz_tgt - int_zd
1117  else ; dz0 = 0.0 ; z0 = zr1 ; f0 = dz_tgt ; endif
1118 
1119  bq = (dzwt + 2.0*aq*(z0-zs))
1120  ! Solve the quadratic: Aq*(zn-z0)**2 + Bq*(zn-z0) - F0*dzwt = 0
1121  ! Note that b>=0, and the two terms in the standard form cancel for the right root.
1122  dz_g(k) = sgn * (dz0 + 2.0*f0*dzwt / (bq + sqrt(bq**2 + 4.0*aq*f0*dzwt) ))
1123 
1124 ! if (debug) then
1125 ! dz0 = zs-zr1 ; z0 = zs ; F0 = dz_tgt - Int_zs ; Bq = (dzwt + 2.0*Aq*(z0-zs))
1126 ! ddz_g_s(k) = sgn * (dz0 + 2.0*F0*dzwt / (Bq + sqrt(Bq**2 + 4.0*Aq*F0*dzwt) )) - dz_g(k)
1127 ! dz0 = zd-zr1 ; z0 = zd ; F0 = dz_tgt - Int_zd ; Bq = (dzwt + 2.0*Aq*(z0-zs))
1128 ! ddz_g_d(k) = sgn * (dz0 + 2.0*F0*dzwt / (Bq + sqrt(Bq**2 + 4.0*Aq*F0*dzwt) )) - dz_g(k)
1129 !
1130 ! if (abs(ddz_g_s(k)) > 1e-12*(abs(dz_g(k)) + abs(dz_g(k)+ddz_g_s(k)))) &
1131 ! call MOM_error(WARNING, "filtered_grid_motion: Expect z_output to be tangled (sc).")
1132 ! if (abs(ddz_g_d(k) - ddz_g_s(k)) > 1e-12*(abs(dz_g(k)+ddz_g_d(k)) + abs(dz_g(k)+ddz_g_s(k)))) &
1133 ! call MOM_error(WARNING, "filtered_grid_motion: Expect z_output to be tangled.")
1134 ! endif
1135  endif
1136 
1137  endif
1138  enddo
1139  !dz_g(CS%nk+1) = 0.0
1140 
1141  if (debug) then
1142  z_old_k = z_old(1)
1143  do k=1,cs%nk+1
1144  if (k<=nk+1) z_old_k = z_old(k) ! This allows for virtual z_old interface at bottom of the model
1145  z_act(k) = z_old_k + dz_g(k)
1146  enddo
1147  do k=2,cs%nk+1
1148  if (sgn*((z_act(k))-z_act(k-1)) < -1e-15*(abs(z_act(k))+abs(z_act(k-1))) ) &
1149  call mom_error(fatal, "filtered_grid_motion: z_output is tangled.")
1150  enddo
1151  endif
1152 
1153 end subroutine filtered_grid_motion
1154 
1155 !> Builds a z*-ccordinate grid with partial steps (Adcroft and Campin, 2004).
1156 !! z* is defined as
1157 !! z* = (z-eta)/(H+eta)*H s.t. z*=0 when z=eta and z*=-H when z=-H .
1158 subroutine build_zstar_grid( CS, G, GV, h, dzInterface, frac_shelf_h)
1160  ! Arguments
1161  type(regridding_cs), intent(in) :: CS !< Regridding control structure
1162  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
1163  type(verticalgrid_type), intent(in) :: GV !< ocean vertical grid structure
1164  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
1165  real, dimension(SZI_(G),SZJ_(G), CS%nk+1), intent(inout) :: dzInterface !< The change in interface depth
1166  !! [H ~> m or kg m-2].
1167  real, dimension(:,:), optional, pointer :: frac_shelf_h !< Fractional ice shelf coverage [nondim].
1168  ! Local variables
1169  real :: nominalDepth, totalThickness, dh ! Depths and thicknesses [H ~> m or kg m-2]
1170  real, dimension(SZK_(GV)+1) :: zOld, zNew ! Coordinate interface heights [H ~> m or kg m-2]
1171  integer :: i, j, k, nz
1172  logical :: ice_shelf
1173 
1174  nz = gv%ke
1175  ice_shelf = .false.
1176  if (present(frac_shelf_h)) then
1177  if (associated(frac_shelf_h)) ice_shelf = .true.
1178  endif
1179 
1180 !$OMP parallel do default(none) shared(G,GV,dzInterface,CS,nz,h,frac_shelf_h,ice_shelf) &
1181 !$OMP private(nominalDepth,totalThickness,zNew,dh,zOld)
1182  do j = g%jsc-1,g%jec+1
1183  do i = g%isc-1,g%iec+1
1184 
1185  if (g%mask2dT(i,j)==0.) then
1186  dzinterface(i,j,:) = 0.
1187  cycle
1188  endif
1189 
1190  ! Local depth (G%bathyT is positive)
1191  nominaldepth = g%bathyT(i,j)*gv%Z_to_H
1192 
1193  ! Determine water column thickness
1194  totalthickness = 0.0
1195  do k = 1,nz
1196  totalthickness = totalthickness + h(i,j,k)
1197  enddo
1198 
1199  zold(nz+1) = - nominaldepth
1200  do k = nz,1,-1
1201  zold(k) = zold(k+1) + h(i,j,k)
1202  enddo
1203 
1204  if (ice_shelf) then
1205  if (frac_shelf_h(i,j) > 0.) then ! under ice shelf
1206  call build_zstar_column(cs%zlike_CS, nominaldepth, totalthickness, znew, &
1207  z_rigid_top = totalthickness-nominaldepth, &
1208  eta_orig=zold(1), zscale=gv%Z_to_H)
1209  else
1210  call build_zstar_column(cs%zlike_CS, nominaldepth, totalthickness, &
1211  znew, zscale=gv%Z_to_H)
1212  endif
1213  else
1214  call build_zstar_column(cs%zlike_CS, nominaldepth, totalthickness, &
1215  znew, zscale=gv%Z_to_H)
1216  endif
1217 
1218  ! Calculate the final change in grid position after blending new and old grids
1219  call filtered_grid_motion( cs, nz, zold, znew, dzinterface(i,j,:) )
1220 
1221 #ifdef __DO_SAFETY_CHECKS__
1222  dh=max(nominaldepth,totalthickness)
1223  if (abs(znew(1)-zold(1))>(nz-1)*0.5*epsilon(dh)*dh) then
1224  write(0,*) 'min_thickness=',cs%min_thickness
1225  write(0,*) 'nominalDepth=',nominaldepth,'totalThickness=',totalthickness
1226  write(0,*) 'dzInterface(1) = ',dzinterface(i,j,1),epsilon(dh),nz
1227  do k=1,nz+1
1228  write(0,*) k,zold(k),znew(k)
1229  enddo
1230  do k=1,nz
1231  write(0,*) k,h(i,j,k),znew(k)-znew(k+1),cs%coordinateResolution(k)
1232  enddo
1233  call mom_error( fatal, &
1234  'MOM_regridding, build_zstar_grid(): top surface has moved!!!' )
1235  endif
1236 #endif
1237 
1238  call adjust_interface_motion( cs, nz, h(i,j,:), dzinterface(i,j,:) )
1239 
1240  enddo
1241  enddo
1242 
1243 end subroutine build_zstar_grid
1244 
1245 !------------------------------------------------------------------------------
1246 ! Build sigma grid
1247 !> This routine builds a grid based on terrain-following coordinates.
1248 subroutine build_sigma_grid( CS, G, GV, h, dzInterface )
1249 !------------------------------------------------------------------------------
1250 ! This routine builds a grid based on terrain-following coordinates.
1251 ! The module parameter coordinateResolution(:) determines the resolution in
1252 ! sigma coordinate, dSigma(:). sigma-coordinates are defined by
1253 ! sigma = (eta-z)/(H+eta) s.t. sigma=0 at z=eta and sigma=1 at z=-H .
1254 !------------------------------------------------------------------------------
1255 
1256  ! Arguments
1257  type(regridding_cs), intent(in) :: CS !< Regridding control structure
1258  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
1259  type(verticalgrid_type), intent(in) :: GV !< ocean vertical grid structure
1260  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
1261  real, dimension(SZI_(G),SZJ_(G), CS%nk+1), intent(inout) :: dzInterface !< The change in interface depth
1262  !! [H ~> m or kg m-2]
1263 
1264  ! Local variables
1265  integer :: i, j, k
1266  integer :: nz
1267  real :: nominalDepth, totalThickness, dh
1268  real, dimension(SZK_(GV)+1) :: zOld, zNew
1269 
1270  nz = gv%ke
1271 
1272  do i = g%isc-1,g%iec+1
1273  do j = g%jsc-1,g%jec+1
1274 
1275  if (g%mask2dT(i,j)==0.) then
1276  dzinterface(i,j,:) = 0.
1277  cycle
1278  endif
1279 
1280  ! The rest of the model defines grids integrating up from the bottom
1281  nominaldepth = g%bathyT(i,j)*gv%Z_to_H
1282 
1283  ! Determine water column height
1284  totalthickness = 0.0
1285  do k = 1,nz
1286  totalthickness = totalthickness + h(i,j,k)
1287  enddo
1288 
1289  call build_sigma_column(cs%sigma_CS, nominaldepth, totalthickness, znew)
1290 
1291  ! Calculate the final change in grid position after blending new and old grids
1292  zold(nz+1) = -nominaldepth
1293  do k = nz,1,-1
1294  zold(k) = zold(k+1) + h(i, j, k)
1295  enddo
1296 
1297  call filtered_grid_motion( cs, nz, zold, znew, dzinterface(i,j,:) )
1298 
1299 #ifdef __DO_SAFETY_CHECKS__
1300  dh=max(nominaldepth,totalthickness)
1301  if (abs(znew(1)-zold(1))>(cs%nk-1)*0.5*epsilon(dh)*dh) then
1302  write(0,*) 'min_thickness=',cs%min_thickness
1303  write(0,*) 'nominalDepth=',nominaldepth,'totalThickness=',totalthickness
1304  write(0,*) 'dzInterface(1) = ',dzinterface(i,j,1),epsilon(dh),nz,cs%nk
1305  do k=1,nz+1
1306  write(0,*) k,zold(k),znew(k)
1307  enddo
1308  do k=1,cs%nk
1309  write(0,*) k,h(i,j,k),znew(k)-znew(k+1),totalthickness*cs%coordinateResolution(k),cs%coordinateResolution(k)
1310  enddo
1311  call mom_error( fatal, &
1312  'MOM_regridding, build_sigma_grid: top surface has moved!!!' )
1313  endif
1314  dzinterface(i,j,1) = 0.
1315  dzinterface(i,j,cs%nk+1) = 0.
1316 #endif
1317 
1318  enddo
1319  enddo
1320 
1321 end subroutine build_sigma_grid
1322 
1323 !------------------------------------------------------------------------------
1324 ! Build grid based on target interface densities
1325 !------------------------------------------------------------------------------
1326 !> This routine builds a new grid based on a given set of target interface densities.
1327 subroutine build_rho_grid( G, GV, US, h, tv, dzInterface, remapCS, CS )
1328 !------------------------------------------------------------------------------
1329 ! This routine builds a new grid based on a given set of target interface
1330 ! densities (these target densities are computed by taking the mean value
1331 ! of given layer densities). The algorithn operates as follows within each
1332 ! column:
1333 ! 1. Given T & S within each layer, the layer densities are computed.
1334 ! 2. Based on these layer densities, a global density profile is reconstructed
1335 ! (this profile is monotonically increasing and may be discontinuous)
1336 ! 3. The new grid interfaces are determined based on the target interface
1337 ! densities.
1338 ! 4. T & S are remapped onto the new grid.
1339 ! 5. Return to step 1 until convergence or until the maximum number of
1340 ! iterations is reached, whichever comes first.
1341 !------------------------------------------------------------------------------
1342 
1343  ! Arguments
1344  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
1345  type(verticalgrid_type), intent(in) :: GV !< Ocean vertical grid structure
1346  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
1347  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
1348  type(thermo_var_ptrs), intent(in) :: tv !< Thermodynamics structure
1349  real, dimension(SZI_(G),SZJ_(G), SZK_(GV)+1), intent(inout) :: dzInterface !< The change in interface depth
1350  !! [H ~> m or kg m-2]
1351  type(remapping_cs), intent(in) :: remapCS !< The remapping control structure
1352  type(regridding_cs), intent(in) :: CS !< Regridding control structure
1353 
1354  ! Local variables
1355  integer :: nz
1356  integer :: i, j, k
1357  real :: nominalDepth ! Depth of the bottom of the ocean, positive downward [H ~> m or kg m-2]
1358  real, dimension(SZK_(GV)+1) :: zOld, zNew ! Old and new interface heights [H ~> m or kg m-2]
1359  real :: h_neglect, h_neglect_edge ! Negligible thicknesses [H ~> m or kg m-2]
1360 #ifdef __DO_SAFETY_CHECKS__
1361  real :: totalThickness
1362  real :: dh
1363 #endif
1364 
1365  if (.not.cs%remap_answers_2018) then
1366  h_neglect = gv%H_subroundoff ; h_neglect_edge = gv%H_subroundoff
1367  elseif (gv%Boussinesq) then
1368  h_neglect = gv%m_to_H*1.0e-30 ; h_neglect_edge = gv%m_to_H*1.0e-10
1369  else
1370  h_neglect = gv%kg_m2_to_H*1.0e-30 ; h_neglect_edge = gv%kg_m2_to_H*1.0e-10
1371  endif
1372 
1373  nz = gv%ke
1374 
1375  if (.not.cs%target_density_set) call mom_error(fatal, "build_rho_grid: "//&
1376  "Target densities must be set before build_rho_grid is called.")
1377 
1378  ! Build grid based on target interface densities
1379  do j = g%jsc-1,g%jec+1
1380  do i = g%isc-1,g%iec+1
1381 
1382  if (g%mask2dT(i,j)==0.) then
1383  dzinterface(i,j,:) = 0.
1384  cycle
1385  endif
1386 
1387 
1388  ! Local depth (G%bathyT is positive)
1389  nominaldepth = g%bathyT(i,j)*gv%Z_to_H
1390 
1391  call build_rho_column(cs%rho_CS, nz, nominaldepth, h(i, j, :), &
1392  tv%T(i, j, :), tv%S(i, j, :), tv%eqn_of_state, znew, &
1393  h_neglect=h_neglect, h_neglect_edge=h_neglect_edge)
1394 
1395  if (cs%integrate_downward_for_e) then
1396  zold(1) = 0.
1397  do k = 1,nz
1398  zold(k+1) = zold(k) - h(i,j,k)
1399  enddo
1400  else
1401  ! The rest of the model defines grids integrating up from the bottom
1402  zold(nz+1) = - nominaldepth
1403  do k = nz,1,-1
1404  zold(k) = zold(k+1) + h(i,j,k)
1405  enddo
1406  endif
1407 
1408  ! Calculate the final change in grid position after blending new and old grids
1409  call filtered_grid_motion( cs, nz, zold, znew, dzinterface(i,j,:) )
1410 
1411 #ifdef __DO_SAFETY_CHECKS__
1412  do k = 2,nz
1413  if (znew(k) > zold(1)) then
1414  write(0,*) 'zOld=',zold
1415  write(0,*) 'zNew=',znew
1416  call mom_error( fatal, 'MOM_regridding, build_rho_grid: '//&
1417  'interior interface above surface!' )
1418  endif
1419  if (znew(k) > znew(k-1)) then
1420  write(0,*) 'zOld=',zold
1421  write(0,*) 'zNew=',znew
1422  call mom_error( fatal, 'MOM_regridding, build_rho_grid: '//&
1423  'interior interfaces cross!' )
1424  endif
1425  enddo
1426 
1427  totalthickness = 0.0
1428  do k = 1,nz
1429  totalthickness = totalthickness + h(i,j,k)
1430  enddo
1431 
1432  dh=max(nominaldepth,totalthickness)
1433  if (abs(znew(1)-zold(1))>(nz-1)*0.5*epsilon(dh)*dh) then
1434  write(0,*) 'min_thickness=',cs%min_thickness
1435  write(0,*) 'nominalDepth=',nominaldepth,'totalThickness=',totalthickness
1436  write(0,*) 'zNew(1)-zOld(1) = ',znew(1)-zold(1),epsilon(dh),nz
1437  do k=1,nz+1
1438  write(0,*) k,zold(k),znew(k)
1439  enddo
1440  do k=1,nz
1441  write(0,*) k,h(i,j,k),znew(k)-znew(k+1)
1442  enddo
1443  call mom_error( fatal, &
1444  'MOM_regridding, build_rho_grid: top surface has moved!!!' )
1445  endif
1446 #endif
1447 
1448  enddo ! end loop on i
1449  enddo ! end loop on j
1450 
1451 end subroutine build_rho_grid
1452 
1453 !> Builds a simple HyCOM-like grid with the deepest location of potential
1454 !! density interpolated from the column profile and a clipping of depth for
1455 !! each interface to a fixed z* or p* grid. This should probably be (optionally?)
1456 !! changed to find the nearest location of the target density.
1457 !! \remark { Based on Bleck, 2002: An oceanice general circulation model framed in
1458 !! hybrid isopycnic-Cartesian coordinates, Ocean Modelling 37, 55-88.
1459 !! http://dx.doi.org/10.1016/S1463-5003(01)00012-9 }
1460 subroutine build_grid_hycom1( G, GV, US, h, tv, h_new, dzInterface, CS )
1461  type(ocean_grid_type), intent(in) :: G !< Grid structure
1462  type(verticalgrid_type), intent(in) :: GV !< Ocean vertical grid structure
1463  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
1464  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: h !< Existing model thickness [H ~> m or kg m-2]
1465  type(thermo_var_ptrs), intent(in) :: tv !< Thermodynamics structure
1466  type(regridding_cs), intent(in) :: CS !< Regridding control structure
1467  real, dimension(SZI_(G),SZJ_(G),CS%nk), intent(inout) :: h_new !< New layer thicknesses [H ~> m or kg m-2]
1468  real, dimension(SZI_(G),SZJ_(G),CS%nk+1), intent(inout) :: dzInterface !< Changes in interface position
1469 
1470  ! Local variables
1471  real, dimension(SZK_(GV)+1) :: z_col ! Source interface positions relative to the surface [H ~> m or kg m-2]
1472  real, dimension(CS%nk+1) :: z_col_new ! New interface positions relative to the surface [H ~> m or kg m-2]
1473  real, dimension(SZK_(GV)+1) :: dz_col ! The realized change in z_col [H ~> m or kg m-2]
1474  real, dimension(SZK_(GV)) :: p_col ! Layer center pressure [R L2 T-2 ~> Pa]
1475  real :: ref_pres ! The reference pressure [R L2 T-2 ~> Pa]
1476  integer :: i, j, k, nki
1477  real :: depth
1478  real :: h_neglect, h_neglect_edge
1479 
1480  if (.not.cs%remap_answers_2018) then
1481  h_neglect = gv%H_subroundoff ; h_neglect_edge = gv%H_subroundoff
1482  elseif (gv%Boussinesq) then
1483  h_neglect = gv%m_to_H*1.0e-30 ; h_neglect_edge = gv%m_to_H*1.0e-10
1484  else
1485  h_neglect = gv%kg_m2_to_H*1.0e-30 ; h_neglect_edge = gv%kg_m2_to_H*1.0e-10
1486  endif
1487 
1488  if (.not.cs%target_density_set) call mom_error(fatal, "build_grid_HyCOM1 : "//&
1489  "Target densities must be set before build_grid_HyCOM1 is called.")
1490 
1491  nki = min(gv%ke, cs%nk)
1492 
1493  ! Build grid based on target interface densities
1494  do j = g%jsc-1,g%jec+1 ; do i = g%isc-1,g%iec+1
1495  if (g%mask2dT(i,j)>0.) then
1496 
1497  depth = g%bathyT(i,j) * gv%Z_to_H
1498 
1499  z_col(1) = 0. ! Work downward rather than bottom up
1500  do k = 1, gv%ke
1501  z_col(k+1) = z_col(k) + h(i,j,k)
1502  p_col(k) = tv%P_Ref + cs%compressibility_fraction * &
1503  ( 0.5 * ( z_col(k) + z_col(k+1) ) * (gv%H_to_RZ*gv%g_Earth) - tv%P_Ref )
1504  enddo
1505 
1506  call build_hycom1_column(cs%hycom_CS, tv%eqn_of_state, gv%ke, depth, &
1507  h(i,j,:), tv%T(i,j,:), tv%S(i,j,:), p_col, &
1508  z_col, z_col_new, zscale=gv%Z_to_H, &
1509  h_neglect=h_neglect, h_neglect_edge=h_neglect_edge)
1510 
1511  ! Calculate the final change in grid position after blending new and old grids
1512  call filtered_grid_motion( cs, gv%ke, z_col, z_col_new, dz_col )
1513 
1514  ! This adjusts things robust to round-off errors
1515  dz_col(:) = -dz_col(:)
1516  call adjust_interface_motion( cs, gv%ke, h(i,j,:), dz_col(:) )
1517 
1518  dzinterface(i,j,1:nki+1) = dz_col(1:nki+1)
1519  if (nki<cs%nk) dzinterface(i,j,nki+2:cs%nk+1) = 0.
1520 
1521  else ! on land
1522  dzinterface(i,j,:) = 0.
1523  endif ! mask2dT
1524  enddo ; enddo ! i,j
1525 
1526  call calc_h_new_by_dz(cs, g, gv, h, dzinterface, h_new)
1527 
1528 end subroutine build_grid_hycom1
1529 
1530 !> This subroutine builds an adaptive grid that follows density surfaces where
1531 !! possible, subject to constraints on the smoothness of interface heights.
1532 subroutine build_grid_adaptive(G, GV, US, h, tv, dzInterface, remapCS, CS)
1533  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1534  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
1535  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
1536  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
1537  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various
1538  !! thermodynamic variables
1539  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)+1), intent(inout) :: dzInterface !< The change in interface depth
1540  !! [H ~> m or kg m-2]
1541  type(remapping_cs), intent(in) :: remapCS !< The remapping control structure
1542  type(regridding_cs), intent(in) :: CS !< Regridding control structure
1543 
1544  ! local variables
1545  integer :: i, j, k, nz ! indices and dimension lengths
1546  ! temperature, salinity and pressure on interfaces
1547  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)+1) :: tInt, sInt
1548  ! current interface positions and after tendency term is applied
1549  ! positive downward
1550  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)+1) :: zInt ! Interface depths [H ~> m or kg m-2]
1551  real, dimension(SZK_(GV)+1) :: zNext ! New interface depths [H ~> m or kg m-2]
1552 
1553  nz = gv%ke
1554 
1555  ! position surface at z = 0.
1556  zint(:,:,1) = 0.
1557 
1558  ! work on interior interfaces
1559  do k = 2, nz ; do j = g%jsc-2,g%jec+2 ; do i = g%isc-2,g%iec+2
1560  tint(i,j,k) = 0.5 * (tv%T(i,j,k-1) + tv%T(i,j,k))
1561  sint(i,j,k) = 0.5 * (tv%S(i,j,k-1) + tv%S(i,j,k))
1562  zint(i,j,k) = zint(i,j,k-1) + h(i,j,k-1) ! zInt in [H]
1563  enddo ; enddo ; enddo
1564 
1565  ! top and bottom temp/salt interfaces are just the layer
1566  ! average values
1567  tint(:,:,1) = tv%T(:,:,1) ; tint(:,:,nz+1) = tv%T(:,:,nz)
1568  sint(:,:,1) = tv%S(:,:,1) ; sint(:,:,nz+1) = tv%S(:,:,nz)
1569 
1570  ! set the bottom interface depth
1571  zint(:,:,nz+1) = zint(:,:,nz) + h(:,:,nz)
1572 
1573  ! calculate horizontal density derivatives (alpha/beta)
1574  ! between cells in a 5-point stencil, columnwise
1575  do j = g%jsc-1,g%jec+1 ; do i = g%isc-1,g%iec+1
1576  if (g%mask2dT(i,j) < 0.5) then
1577  dzinterface(i,j,:) = 0. ! land point, don't move interfaces, and skip
1578  cycle
1579  endif
1580 
1581  call build_adapt_column(cs%adapt_CS, g, gv, us, tv, i, j, zint, tint, sint, h, znext)
1582 
1583  call filtered_grid_motion(cs, nz, zint(i,j,:), znext, dzinterface(i,j,:))
1584  ! convert from depth to z
1585  do k = 1, nz+1 ; dzinterface(i,j,k) = -dzinterface(i,j,k) ; enddo
1586  call adjust_interface_motion(cs, nz, h(i,j,:), dzinterface(i,j,:))
1587  enddo ; enddo
1588 end subroutine build_grid_adaptive
1589 
1590 !> Builds a grid that tracks density interfaces for water that is denser than
1591 !! the surface density plus an increment of some number of layers, and uses all
1592 !! lighter layers uniformly above this location. Note that this amounts to
1593 !! interpolating to find the depth of an arbitrary (non-integer) interface index
1594 !! which should make the results vary smoothly in space to the extent that the
1595 !! surface density and interior stratification vary smoothly in space. Over
1596 !! shallow topography, this will tend to give a uniform sigma-like coordinate.
1597 !! For sufficiently shallow water, a minimum grid spacing is used to avoid
1598 !! certain instabilities.
1599 subroutine build_grid_slight(G, GV, US, h, tv, dzInterface, CS)
1600  type(ocean_grid_type), intent(in) :: G !< Grid structure
1601  type(verticalgrid_type), intent(in) :: GV !< Ocean vertical grid structure
1602  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
1603  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: h !< Existing model thickness [H ~> m or kg m-2]
1604  type(thermo_var_ptrs), intent(in) :: tv !< Thermodynamics structure
1605  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)+1), intent(inout) :: dzInterface !< Changes in interface position
1606  type(regridding_cs), intent(in) :: CS !< Regridding control structure
1607 
1608  real, dimension(SZK_(GV)+1) :: z_col ! Interface positions relative to the surface [H ~> m or kg m-2]
1609  real, dimension(SZK_(GV)+1) :: z_col_new ! Interface positions relative to the surface [H ~> m or kg m-2]
1610  real, dimension(SZK_(GV)+1) :: dz_col ! The realized change in z_col [H ~> m or kg m-2]
1611  real, dimension(SZK_(GV)) :: p_col ! Layer center pressure [R L2 T-2 ~> Pa]
1612  real :: depth
1613  integer :: i, j, k, nz
1614  real :: h_neglect, h_neglect_edge
1615 
1616  if (.not.cs%remap_answers_2018) then
1617  h_neglect = gv%H_subroundoff ; h_neglect_edge = gv%H_subroundoff
1618  elseif (gv%Boussinesq) then
1619  h_neglect = gv%m_to_H*1.0e-30 ; h_neglect_edge = gv%m_to_H*1.0e-10
1620  else
1621  h_neglect = gv%kg_m2_to_H*1.0e-30 ; h_neglect_edge = gv%kg_m2_to_H*1.0e-10
1622  endif
1623 
1624  nz = gv%ke
1625 
1626  if (.not.cs%target_density_set) call mom_error(fatal, "build_grid_SLight : "//&
1627  "Target densities must be set before build_grid_SLight is called.")
1628 
1629  ! Build grid based on target interface densities
1630  do j = g%jsc-1,g%jec+1 ; do i = g%isc-1,g%iec+1
1631  if (g%mask2dT(i,j)>0.) then
1632 
1633  depth = g%bathyT(i,j) * gv%Z_to_H
1634  z_col(1) = 0. ! Work downward rather than bottom up
1635  do k=1,nz
1636  z_col(k+1) = z_col(k) + h(i,j,k)
1637  p_col(k) = tv%P_Ref + cs%compressibility_fraction * &
1638  ( 0.5 * ( z_col(k) + z_col(k+1) ) * (gv%H_to_RZ*gv%g_Earth) - tv%P_Ref )
1639  enddo
1640 
1641  call build_slight_column(cs%slight_CS, tv%eqn_of_state, gv%H_to_RZ*gv%g_Earth, &
1642  gv%H_subroundoff, nz, depth, h(i, j, :), &
1643  tv%T(i, j, :), tv%S(i, j, :), p_col, z_col, z_col_new, &
1644  h_neglect=h_neglect, h_neglect_edge=h_neglect_edge)
1645 
1646  ! Calculate the final change in grid position after blending new and old grids
1647  call filtered_grid_motion( cs, nz, z_col, z_col_new, dz_col )
1648  do k=1,nz+1 ; dzinterface(i,j,k) = -dz_col(k) ; enddo
1649 #ifdef __DO_SAFETY_CHECKS__
1650  if (dzinterface(i,j,1) /= 0.) stop 'build_grid_SLight: Surface moved?!'
1651  if (dzinterface(i,j,nz+1) /= 0.) stop 'build_grid_SLight: Bottom moved?!'
1652 #endif
1653 
1654  ! This adjusts things robust to round-off errors
1655  call adjust_interface_motion( cs, nz, h(i,j,:), dzinterface(i,j,:) )
1656 
1657  else ! on land
1658  dzinterface(i,j,:) = 0.
1659  endif ! mask2dT
1660  enddo ; enddo ! i,j
1661 
1662 end subroutine build_grid_slight
1663 
1664 !> Adjust dz_Interface to ensure non-negative future thicknesses
1665 subroutine adjust_interface_motion( CS, nk, h_old, dz_int )
1666  type(regridding_cs), intent(in) :: CS !< Regridding control structure
1667  integer, intent(in) :: nk !< Number of layers in h_old
1668  real, dimension(nk), intent(in) :: h_old !< Minium allowed thickness of h [H ~> m or kg m-2]
1669  real, dimension(CS%nk+1), intent(inout) :: dz_int !< Minium allowed thickness of h [H ~> m or kg m-2]
1670  ! Local variables
1671  integer :: k
1672  real :: h_new, eps, h_total, h_err
1673 
1674  eps = 1. ; eps = epsilon(eps)
1675 
1676  h_total = 0. ; h_err = 0.
1677  do k = 1, min(cs%nk,nk)
1678  h_total = h_total + h_old(k)
1679  h_err = h_err + max( h_old(k), abs(dz_int(k)), abs(dz_int(k+1)) )*eps
1680  h_new = h_old(k) + ( dz_int(k) - dz_int(k+1) )
1681  if (h_new < -3.0*h_err) then
1682  write(0,*) 'h<0 at k=',k,'h_old=',h_old(k), &
1683  'wup=',dz_int(k),'wdn=',dz_int(k+1),'dw_dz=',dz_int(k) - dz_int(k+1), &
1684  'h_new=',h_new,'h_err=',h_err
1685  call mom_error( fatal, 'MOM_regridding: adjust_interface_motion() - '//&
1686  'implied h<0 is larger than roundoff!')
1687  endif
1688  enddo
1689  if (cs%nk>nk) then
1690  do k = nk+1, cs%nk
1691  h_err = h_err + max( abs(dz_int(k)), abs(dz_int(k+1)) )*eps
1692  h_new = ( dz_int(k) - dz_int(k+1) )
1693  if (h_new < -3.0*h_err) then
1694  write(0,*) 'h<0 at k=',k,'h_old was empty',&
1695  'wup=',dz_int(k),'wdn=',dz_int(k+1),'dw_dz=',dz_int(k) - dz_int(k+1), &
1696  'h_new=',h_new,'h_err=',h_err
1697  call mom_error( fatal, 'MOM_regridding: adjust_interface_motion() - '//&
1698  'implied h<0 is larger than roundoff!')
1699  endif
1700  enddo
1701  endif
1702  do k = min(cs%nk,nk),2,-1
1703  h_new = h_old(k) + ( dz_int(k) - dz_int(k+1) )
1704  if (h_new<cs%min_thickness) &
1705  dz_int(k) = ( dz_int(k+1) - h_old(k) ) + cs%min_thickness ! Implies next h_new = min_thickness
1706  h_new = h_old(k) + ( dz_int(k) - dz_int(k+1) )
1707  if (h_new<0.) &
1708  dz_int(k) = ( 1. - eps ) * ( dz_int(k+1) - h_old(k) ) ! Backup in case min_thickness==0
1709  h_new = h_old(k) + ( dz_int(k) - dz_int(k+1) )
1710  if (h_new<0.) then
1711  write(0,*) 'h<0 at k=',k,'h_old=',h_old(k), &
1712  'wup=',dz_int(k),'wdn=',dz_int(k+1),'dw_dz=',dz_int(k) - dz_int(k+1), &
1713  'h_new=',h_new
1714  stop 'Still did not work!'
1715  call mom_error( fatal, 'MOM_regridding: adjust_interface_motion() - '//&
1716  'Repeated adjustment for roundoff h<0 failed!')
1717  endif
1718  enddo
1719  !if (dz_int(1)/=0.) stop 'MOM_regridding: adjust_interface_motion() surface moved'
1720 
1721 end subroutine adjust_interface_motion
1722 
1723 !------------------------------------------------------------------------------
1724 ! Build arbitrary grid
1725 !------------------------------------------------------------------------------
1726 subroutine build_grid_arbitrary( G, GV, h, dzInterface, h_new, CS )
1727 !------------------------------------------------------------------------------
1728 ! This routine builds a grid based on arbitrary rules
1729 !------------------------------------------------------------------------------
1730 
1731  ! Arguments
1732  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
1733  type(verticalgrid_type), intent(in) :: GV !< Ocean vertical grid structure
1734  real, dimension(SZI_(G),SZJ_(G), SZK_(GV)), intent(in) :: h !< Original layer thicknesses [H ~> m or kg m-2]
1735  real, dimension(SZI_(G),SZJ_(G), SZK_(GV)+1), intent(inout) :: dzInterface !< The change in interface
1736  !! depth [H ~> m or kg m-2]
1737  real, intent(inout) :: h_new !< New layer thicknesses [H ~> m or kg m-2]
1738  type(regridding_cs), intent(in) :: CS !< Regridding control structure
1739 
1740  ! Local variables
1741  integer :: i, j, k
1742  integer :: nz
1743  real :: z_inter(SZK_(GV)+1)
1744  real :: total_height
1745  real :: delta_h
1746  real :: max_depth
1747  real :: eta ! local elevation
1748  real :: local_depth
1749  real :: x1, y1, x2, y2
1750  real :: x, t
1751 
1752  nz = gv%ke
1753  max_depth = g%max_depth*gv%Z_to_H
1754 
1755  do j = g%jsc-1,g%jec+1
1756  do i = g%isc-1,g%iec+1
1757 
1758  ! Local depth
1759  local_depth = g%bathyT(i,j)*gv%Z_to_H
1760 
1761  ! Determine water column height
1762  total_height = 0.0
1763  do k = 1,nz
1764  total_height = total_height + h(i,j,k)
1765  enddo
1766 
1767  eta = total_height - local_depth
1768 
1769  ! Compute new thicknesses based on stretched water column
1770  delta_h = (max_depth + eta) / nz
1771 
1772  ! Define interfaces
1773  z_inter(1) = eta
1774  do k = 1,nz
1775  z_inter(k+1) = z_inter(k) - delta_h
1776  enddo
1777 
1778  ! Refine grid in the middle
1779  do k = 1,nz+1
1780  x1 = 0.35; y1 = 0.45; x2 = 0.65; y2 = 0.55
1781 
1782  x = - ( z_inter(k) - eta ) / max_depth
1783 
1784  if ( x <= x1 ) then
1785  t = y1*x/x1
1786  elseif ( (x > x1 ) .and. ( x < x2 )) then
1787  t = y1 + (y2-y1) * (x-x1) / (x2-x1)
1788  else
1789  t = y2 + (1.0-y2) * (x-x2) / (1.0-x2)
1790  endif
1791 
1792  z_inter(k) = -t * max_depth + eta
1793 
1794  enddo
1795 
1796  ! Modify interface heights to account for topography
1797  z_inter(nz+1) = - local_depth
1798 
1799  ! Modify interface heights to avoid layers of zero thicknesses
1800  do k = nz,1,-1
1801  if ( z_inter(k) < (z_inter(k+1) + cs%min_thickness) ) then
1802  z_inter(k) = z_inter(k+1) + cs%min_thickness
1803  endif
1804  enddo
1805 
1806  ! Chnage in interface position
1807  x = 0. ! Left boundary at x=0
1808  dzinterface(i,j,1) = 0.
1809  do k = 2,nz
1810  x = x + h(i,j,k)
1811  dzinterface(i,j,k) = z_inter(k) - x
1812  enddo
1813  dzinterface(i,j,nz+1) = 0.
1814 
1815  enddo
1816  enddo
1817 
1818 stop 'OOOOOOPS' ! For some reason the gnu compiler will not let me delete this
1819  ! routine????
1820 
1821 end subroutine build_grid_arbitrary
1822 
1823 
1824 
1825 !------------------------------------------------------------------------------
1826 ! Check grid integrity
1827 !------------------------------------------------------------------------------
1828 subroutine inflate_vanished_layers_old( CS, G, GV, h )
1829 !------------------------------------------------------------------------------
1830 ! This routine is called when initializing the regridding options. The
1831 ! objective is to make sure all layers are at least as thick as the minimum
1832 ! thickness allowed for regridding purposes (this parameter is set in the
1833 ! MOM_input file or defaulted to 1.0e-3). When layers are too thin, they
1834 ! are inflated up to the minmum thickness.
1835 !------------------------------------------------------------------------------
1836 
1837  ! Arguments
1838  type(regridding_cs), intent(in) :: cs !< Regridding control structure
1839  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
1840  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
1841  real, dimension(SZI_(G),SZJ_(G), SZK_(GV)), intent(inout) :: h !< Layer thicknesses [H ~> m or kg m-2]
1842 
1843  ! Local variables
1844  integer :: i, j, k
1845  real :: htmp(gv%ke)
1846 
1847  do i = g%isc-1,g%iec+1
1848  do j = g%jsc-1,g%jec+1
1849 
1850  ! Build grid for current column
1851  do k = 1,gv%ke
1852  htmp(k) = h(i,j,k)
1853  enddo
1854 
1855  call old_inflate_layers_1d( cs%min_thickness, gv%ke, htmp )
1856 
1857  ! Save modified grid
1858  do k = 1,gv%ke
1859  h(i,j,k) = htmp(k)
1860  enddo
1861 
1862  enddo
1863  enddo
1864 
1865 end subroutine inflate_vanished_layers_old
1866 
1867 !------------------------------------------------------------------------------
1868 !> Achieve convective adjustment by swapping layers
1869 subroutine convective_adjustment(G, GV, h, tv)
1870  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1871  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
1872  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
1873  intent(inout) :: h !< Layer thicknesses [H ~> m or kg m-2]
1874  type(thermo_var_ptrs), intent(inout) :: tv !< A structure pointing to various thermodynamic variables
1875 !------------------------------------------------------------------------------
1876 ! Check each water column to see whether it is stratified. If not, sort the
1877 ! layers by successive swappings of water masses (bubble sort algorithm)
1878 !------------------------------------------------------------------------------
1879 
1880  ! Local variables
1881  integer :: i, j, k
1882  real :: T0, T1 ! temperatures
1883  real :: S0, S1 ! salinities
1884  real :: r0, r1 ! densities
1885  real :: h0, h1
1886  logical :: stratified
1887  real, dimension(GV%ke) :: p_col, densities
1888 
1889  p_col(:) = 0.
1890 
1891  ! Loop on columns
1892  do j = g%jsc-1,g%jec+1 ; do i = g%isc-1,g%iec+1
1893 
1894  ! Compute densities within current water column
1895  call calculate_density( tv%T(i,j,:), tv%S(i,j,:), p_col, densities, tv%eqn_of_state)
1896 
1897  ! Repeat restratification until complete
1898  do
1899  stratified = .true.
1900  do k = 1,gv%ke-1
1901  ! Gather information of current and next cells
1902  t0 = tv%T(i,j,k) ; t1 = tv%T(i,j,k+1)
1903  s0 = tv%S(i,j,k) ; s1 = tv%S(i,j,k+1)
1904  r0 = densities(k) ; r1 = densities(k+1)
1905  h0 = h(i,j,k) ; h1 = h(i,j,k+1)
1906  ! If the density of the current cell is larger than the density
1907  ! below it, we swap the cells and recalculate the densitiies
1908  ! within the swapped cells
1909  if ( r0 > r1 ) then
1910  tv%T(i,j,k) = t1 ; tv%T(i,j,k+1) = t0
1911  tv%S(i,j,k) = s1 ; tv%S(i,j,k+1) = s0
1912  h(i,j,k) = h1 ; h(i,j,k+1) = h0
1913  ! Recompute densities at levels k and k+1
1914  call calculate_density( tv%T(i,j,k), tv%S(i,j,k), p_col(k), densities(k), tv%eqn_of_state)
1915  call calculate_density( tv%T(i,j,k+1), tv%S(i,j,k+1), p_col(k+1), &
1916  densities(k+1), tv%eqn_of_state )
1917  stratified = .false.
1918  endif
1919  enddo ! k
1920 
1921  if ( stratified ) exit
1922  enddo
1923 
1924  enddo ; enddo ! i & j
1925 
1926 end subroutine convective_adjustment
1927 
1928 
1929 !------------------------------------------------------------------------------
1930 !> Return a uniform resolution vector in the units of the coordinate
1931 function uniformresolution(nk,coordMode,maxDepth,rhoLight,rhoHeavy)
1932 !------------------------------------------------------------------------------
1933 ! Calculate a vector of uniform resolution in the units of the coordinate
1934 !------------------------------------------------------------------------------
1935  ! Arguments
1936  integer, intent(in) :: nk !< Number of cells in source grid
1937  character(len=*), intent(in) :: coordmode !< A string indicating the coordinate mode.
1938  !! See the documenttion for regrid_consts
1939  !! for the recognized values.
1940  real, intent(in) :: maxdepth !< The range of the grid values in some modes
1941  real, intent(in) :: rholight !< The minimum value of the grid in RHO mode
1942  real, intent(in) :: rhoheavy !< The maximum value of the grid in RHO mode
1943 
1944  real :: uniformresolution(nk) !< The returned uniform resolution grid.
1945 
1946  ! Local variables
1947  integer :: scheme
1948 
1949  scheme = coordinatemode(coordmode)
1950  select case ( scheme )
1951 
1952  case ( regridding_zstar, regridding_hycom1, regridding_slight, regridding_sigma_shelf_zstar, &
1953  regridding_adaptive )
1954  uniformresolution(:) = maxdepth / real(nk)
1955 
1956  case ( regridding_rho )
1957  uniformresolution(:) = (rhoheavy - rholight) / real(nk)
1958 
1959  case ( regridding_sigma )
1960  uniformresolution(:) = 1. / real(nk)
1961 
1962  case default
1963  call mom_error(fatal, "MOM_regridding, uniformResolution: "//&
1964  "Unrecognized choice for coordinate mode ("//trim(coordmode)//").")
1965 
1966  end select ! type of grid
1967 
1968 end function uniformresolution
1969 
1970 !> Initialize the coordinate resolutions by calling the appropriate initialization
1971 !! routine for the specified coordinate mode.
1972 subroutine initcoord(CS, GV, US, coord_mode)
1973  type(regridding_cs), intent(inout) :: CS !< Regridding control structure
1974  character(len=*), intent(in) :: coord_mode !< A string indicating the coordinate mode.
1975  !! See the documentation for regrid_consts
1976  !! for the recognized values.
1977  type(verticalgrid_type), intent(in) :: GV !< Ocean vertical grid structure
1978  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
1979 
1980  select case (coordinatemode(coord_mode))
1981  case (regridding_zstar)
1982  call init_coord_zlike(cs%zlike_CS, cs%nk, cs%coordinateResolution)
1983  case (regridding_sigma_shelf_zstar)
1984  call init_coord_zlike(cs%zlike_CS, cs%nk, cs%coordinateResolution)
1985  case (regridding_sigma)
1986  call init_coord_sigma(cs%sigma_CS, cs%nk, cs%coordinateResolution)
1987  case (regridding_rho)
1988  call init_coord_rho(cs%rho_CS, cs%nk, cs%ref_pressure, cs%target_density, cs%interp_CS)
1989  case (regridding_hycom1)
1990  call init_coord_hycom(cs%hycom_CS, cs%nk, cs%coordinateResolution, cs%target_density, &
1991  cs%interp_CS)
1992  case (regridding_slight)
1993  call init_coord_slight(cs%slight_CS, cs%nk, cs%ref_pressure, cs%target_density, &
1994  cs%interp_CS, gv%m_to_H)
1995  case (regridding_adaptive)
1996  call init_coord_adapt(cs%adapt_CS, cs%nk, cs%coordinateResolution, gv%m_to_H, us%kg_m3_to_R)
1997  end select
1998 end subroutine initcoord
1999 
2000 !------------------------------------------------------------------------------
2001 !> Set the fixed resolution data
2002 subroutine setcoordinateresolution( dz, CS, scale )
2003  real, dimension(:), intent(in) :: dz !< A vector of vertical grid spacings
2004  type(regridding_cs), intent(inout) :: cs !< Regridding control structure
2005  real, optional, intent(in) :: scale !< A scaling factor converting dz to coordRes
2006 
2007  if (size(dz)/=cs%nk) call mom_error( fatal, &
2008  'setCoordinateResolution: inconsistent number of levels' )
2009 
2010  if (present(scale)) then
2011  cs%coordinateResolution(:) = scale*dz(:)
2012  else
2013  cs%coordinateResolution(:) = dz(:)
2014  endif
2015 
2016 end subroutine setcoordinateresolution
2017 
2018 !> Set target densities based on the old Rlay variable
2019 subroutine set_target_densities_from_gv( GV, US, CS )
2020  type(verticalgrid_type), intent(in) :: gv !< Ocean vertical grid structure
2021  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
2022  type(regridding_cs), intent(inout) :: cs !< Regridding control structure
2023  ! Local variables
2024  integer :: k, nz
2025 
2026  nz = cs%nk
2027  if (nz == 1) then ! Set a broad range of bounds. Regridding may not be meaningful in this case.
2028  cs%target_density(1) = 0.0
2029  cs%target_density(2) = 2.0*gv%Rlay(1)
2030  else
2031  cs%target_density(1) = (gv%Rlay(1) + 0.5*(gv%Rlay(1)-gv%Rlay(2)))
2032  cs%target_density(nz+1) = (gv%Rlay(nz) + 0.5*(gv%Rlay(nz)-gv%Rlay(nz-1)))
2033  do k=2,nz
2034  cs%target_density(k) = cs%target_density(k-1) + cs%coordinateResolution(k)
2035  enddo
2036  endif
2037  cs%target_density_set = .true.
2038 
2039 end subroutine set_target_densities_from_gv
2040 
2041 !> Set target densities based on vector of interface values
2042 subroutine set_target_densities( CS, rho_int )
2043  type(regridding_cs), intent(inout) :: cs !< Regridding control structure
2044  real, dimension(CS%nk+1), intent(in) :: rho_int !< Interface densities [R ~> kg m-3]
2045 
2046  if (size(cs%target_density)/=size(rho_int)) then
2047  call mom_error(fatal, "set_target_densities inconsistent args!")
2048  endif
2049 
2050  cs%target_density(:) = rho_int(:)
2051  cs%target_density_set = .true.
2052 
2053 end subroutine set_target_densities
2054 
2055 !> Set maximum interface depths based on a vector of input values.
2056 subroutine set_regrid_max_depths( CS, max_depths, units_to_H )
2057  type(regridding_cs), intent(inout) :: cs !< Regridding control structure
2058  real, dimension(CS%nk+1), intent(in) :: max_depths !< Maximum interface depths, in arbitrary units
2059  real, optional, intent(in) :: units_to_h !< A conversion factor for max_depths into H units
2060  ! Local variables
2061  real :: val_to_h
2062  integer :: k
2063 
2064  if (.not.allocated(cs%max_interface_depths)) allocate(cs%max_interface_depths(1:cs%nk+1))
2065 
2066  val_to_h = 1.0 ; if (present(units_to_h)) val_to_h = units_to_h
2067  if (max_depths(cs%nk+1) < max_depths(1)) val_to_h = -1.0*val_to_h
2068 
2069  ! Check for sign reversals in the depths.
2070  if (max_depths(cs%nk+1) < max_depths(1)) then
2071  do k=1,cs%nk ; if (max_depths(k+1) > max_depths(k)) &
2072  call mom_error(fatal, "Unordered list of maximum depths sent to set_regrid_max_depths!")
2073  enddo
2074  else
2075  do k=1,cs%nk ; if (max_depths(k+1) < max_depths(k)) &
2076  call mom_error(fatal, "Unordered list of maximum depths sent to set_regrid_max_depths.")
2077  enddo
2078  endif
2079 
2080  do k=1,cs%nk+1
2081  cs%max_interface_depths(k) = val_to_h * max_depths(k)
2082  enddo
2083 
2084  ! set max depths for coordinate
2085  select case (cs%regridding_scheme)
2086  case (regridding_hycom1)
2087  call set_hycom_params(cs%hycom_CS, max_interface_depths=cs%max_interface_depths)
2088  case (regridding_slight)
2089  call set_slight_params(cs%slight_CS, max_interface_depths=cs%max_interface_depths)
2090  end select
2091 end subroutine set_regrid_max_depths
2092 
2093 !> Set maximum layer thicknesses based on a vector of input values.
2094 subroutine set_regrid_max_thickness( CS, max_h, units_to_H )
2095  type(regridding_cs), intent(inout) :: cs !< Regridding control structure
2096  real, dimension(CS%nk+1), intent(in) :: max_h !< Maximum interface depths, in arbitrary units
2097  real, optional, intent(in) :: units_to_h !< A conversion factor for max_h into H units
2098  ! Local variables
2099  real :: val_to_h
2100  integer :: k
2101 
2102  if (.not.allocated(cs%max_layer_thickness)) allocate(cs%max_layer_thickness(1:cs%nk))
2103 
2104  val_to_h = 1.0 ; if (present( units_to_h)) val_to_h = units_to_h
2105 
2106  do k=1,cs%nk
2107  cs%max_layer_thickness(k) = val_to_h * max_h(k)
2108  enddo
2109 
2110  ! set max thickness for coordinate
2111  select case (cs%regridding_scheme)
2112  case (regridding_hycom1)
2113  call set_hycom_params(cs%hycom_CS, max_layer_thickness=cs%max_layer_thickness)
2114  case (regridding_slight)
2115  call set_slight_params(cs%slight_CS, max_layer_thickness=cs%max_layer_thickness)
2116  end select
2117 end subroutine set_regrid_max_thickness
2118 
2119 
2120 !------------------------------------------------------------------------------
2121 !> Query the fixed resolution data
2122 function getcoordinateresolution( CS, undo_scaling )
2123  type(regridding_cs), intent(in) :: cs !< Regridding control structure
2124  logical, optional, intent(in) :: undo_scaling !< If present and true, undo any internal
2125  !! rescaling of the resolution data.
2126  real, dimension(CS%nk) :: getcoordinateresolution
2127 
2128  logical :: unscale
2129  unscale = .false. ; if (present(undo_scaling)) unscale = undo_scaling
2130 
2131  if (unscale) then
2132  getcoordinateresolution(:) = cs%coord_scale * cs%coordinateResolution(:)
2133  else
2134  getcoordinateresolution(:) = cs%coordinateResolution(:)
2135  endif
2136 
2137 end function getcoordinateresolution
2138 
2139 !> Query the target coordinate interface positions
2140 function getcoordinateinterfaces( CS, undo_scaling )
2141  type(regridding_cs), intent(in) :: cs !< Regridding control structure
2142  logical, optional, intent(in) :: undo_scaling !< If present and true, undo any internal
2143  !! rescaling of the resolution data.
2144  real, dimension(CS%nk+1) :: getcoordinateinterfaces !< Interface positions in target coordinate
2145 
2146  integer :: k
2147  logical :: unscale
2148  unscale = .false. ; if (present(undo_scaling)) unscale = undo_scaling
2149 
2150  ! When using a coordinate with target densities, we need to get the actual
2151  ! densities, rather than computing the interfaces based on resolution
2152  if (cs%regridding_scheme == regridding_rho) then
2153  if (.not. cs%target_density_set) &
2154  call mom_error(fatal, 'MOM_regridding, getCoordinateInterfaces: '//&
2155  'target densities not set!')
2156 
2157  if (unscale) then
2158  getcoordinateinterfaces(:) = cs%coord_scale * cs%target_density(:)
2159  else
2160  getcoordinateinterfaces(:) = cs%target_density(:)
2161  endif
2162  else
2163  if (unscale) then
2164  getcoordinateinterfaces(1) = 0.
2165  do k = 1, cs%nk
2166  getcoordinateinterfaces(k+1) = getcoordinateinterfaces(k) - &
2167  cs%coord_scale * cs%coordinateResolution(k)
2168  enddo
2169  else
2170  getcoordinateinterfaces(1) = 0.
2171  do k = 1, cs%nk
2172  getcoordinateinterfaces(k+1) = getcoordinateinterfaces(k) - &
2173  cs%coordinateResolution(k)
2174  enddo
2175  endif
2176  ! The following line has an "abs()" to allow ferret users to reference
2177  ! data by index. It is a temporary work around... :( -AJA
2178  getcoordinateinterfaces(:) = abs( getcoordinateinterfaces(:) )
2179  endif
2180 
2181 end function getcoordinateinterfaces
2182 
2183 !------------------------------------------------------------------------------
2184 !> Query the target coordinate units
2185 function getcoordinateunits( CS )
2186  type(regridding_cs), intent(in) :: cs !< Regridding control structure
2187  character(len=20) :: getcoordinateunits
2188 
2189  select case ( cs%regridding_scheme )
2190  case ( regridding_zstar, regridding_hycom1, regridding_slight, regridding_adaptive )
2191  getcoordinateunits = 'meter'
2192  case ( regridding_sigma_shelf_zstar )
2193  getcoordinateunits = 'meter/fraction'
2194  case ( regridding_sigma )
2195  getcoordinateunits = 'fraction'
2196  case ( regridding_rho )
2197  getcoordinateunits = 'kg/m3'
2198  case ( regridding_arbitrary )
2199  getcoordinateunits = 'unknown'
2200  case default
2201  call mom_error(fatal,'MOM_regridding, getCoordinateUnits: '//&
2202  'Unknown regridding scheme selected!')
2203  end select ! type of grid
2204 
2205 end function getcoordinateunits
2206 
2207 !------------------------------------------------------------------------------
2208 !> Query the short name of the coordinate
2209 function getcoordinateshortname( CS )
2210  type(regridding_cs), intent(in) :: cs !< Regridding control structure
2211  character(len=20) :: getcoordinateshortname
2212 
2213  select case ( cs%regridding_scheme )
2214  case ( regridding_zstar )
2215  !getCoordinateShortName = 'z*'
2216  ! The following line is a temporary work around... :( -AJA
2217  getcoordinateshortname = 'pseudo-depth, -z*'
2218  case ( regridding_sigma_shelf_zstar )
2219  getcoordinateshortname = 'pseudo-depth, -z*/sigma'
2220  case ( regridding_sigma )
2221  getcoordinateshortname = 'sigma'
2222  case ( regridding_rho )
2223  getcoordinateshortname = 'rho'
2224  case ( regridding_arbitrary )
2225  getcoordinateshortname = 'coordinate'
2226  case ( regridding_hycom1 )
2227  getcoordinateshortname = 'z-rho'
2228  case ( regridding_slight )
2229  getcoordinateshortname = 's-rho'
2230  case ( regridding_adaptive )
2231  getcoordinateshortname = 'adaptive'
2232  case default
2233  call mom_error(fatal,'MOM_regridding, getCoordinateShortName: '//&
2234  'Unknown regridding scheme selected!')
2235  end select ! type of grid
2236 
2237 end function getcoordinateshortname
2238 
2239 !> Can be used to set any of the parameters for MOM_regridding.
2240 subroutine set_regrid_params( CS, boundary_extrapolation, min_thickness, old_grid_weight, &
2241  interp_scheme, depth_of_time_filter_shallow, depth_of_time_filter_deep, &
2242  compress_fraction, ref_pressure, dz_min_surface, nz_fixed_surface, Rho_ML_avg_depth, &
2243  nlay_ML_to_interior, fix_haloclines, halocline_filt_len, &
2244  halocline_strat_tol, integrate_downward_for_e, remap_answers_2018, &
2245  adaptTimeRatio, adaptZoom, adaptZoomCoeff, adaptBuoyCoeff, adaptAlpha, adaptDoMin, adaptDrho0)
2246  type(regridding_cs), intent(inout) :: cs !< Regridding control structure
2247  logical, optional, intent(in) :: boundary_extrapolation !< Extrapolate in boundary cells
2248  real, optional, intent(in) :: min_thickness !< Minimum thickness allowed when building the
2249  !! new grid [H ~> m or kg m-2]
2250  real, optional, intent(in) :: old_grid_weight !< Weight given to old coordinate when time-filtering grid
2251  character(len=*), optional, intent(in) :: interp_scheme !< Interpolation method for state-dependent coordinates
2252  real, optional, intent(in) :: depth_of_time_filter_shallow !< Depth to start cubic [H ~> m or kg m-2]
2253  real, optional, intent(in) :: depth_of_time_filter_deep !< Depth to end cubic [H ~> m or kg m-2]
2254  real, optional, intent(in) :: compress_fraction !< Fraction of compressibility to add to potential density [nondim]
2255  real, optional, intent(in) :: ref_pressure !< The reference pressure for density-dependent
2256  !! coordinates [R L2 T-2 ~> Pa]
2257  real, optional, intent(in) :: dz_min_surface !< The fixed resolution in the topmost
2258  !! SLight_nkml_min layers [H ~> m or kg m-2]
2259  integer, optional, intent(in) :: nz_fixed_surface !< The number of fixed-thickness layers at the top of the model
2260  real, optional, intent(in) :: rho_ml_avg_depth !< Averaging depth over which to determine mixed layer potential
2261  !! density [H ~> m or kg m-2]
2262  real, optional, intent(in) :: nlay_ml_to_interior !< Number of layers to offset the mixed layer density to find
2263  !! resolved stratification [nondim]
2264  logical, optional, intent(in) :: fix_haloclines !< Detect regions with much weaker stratification in the coordinate
2265  real, optional, intent(in) :: halocline_filt_len !< Length scale over which to filter T & S when looking for
2266  !! spuriously unstable water mass profiles [H ~> m or kg m-2]
2267  real, optional, intent(in) :: halocline_strat_tol !< Value of the stratification ratio that defines a problematic
2268  !! halocline region.
2269  logical, optional, intent(in) :: integrate_downward_for_e !< If true, integrate for interface positions downward
2270  !! from the top.
2271  logical, optional, intent(in) :: remap_answers_2018 !< If true, use the order of arithmetic and expressions
2272  !! that recover the remapping answers from 2018. Otherwise
2273  !! use more robust but mathematically equivalent expressions.
2274  real, optional, intent(in) :: adapttimeratio !< Ratio of the ALE timestep to the grid timescale [nondim].
2275  real, optional, intent(in) :: adaptzoom !< Depth of near-surface zooming region [H ~> m or kg m-2].
2276  real, optional, intent(in) :: adaptzoomcoeff !< Coefficient of near-surface zooming diffusivity [nondim].
2277  real, optional, intent(in) :: adaptbuoycoeff !< Coefficient of buoyancy diffusivity [nondim].
2278  real, optional, intent(in) :: adaptalpha !< Scaling factor on optimization tendency [nondim].
2279  logical, optional, intent(in) :: adaptdomin !< If true, make a HyCOM-like mixed layer by
2280  !! preventing interfaces from being shallower than
2281  !! the depths specified by the regridding coordinate.
2282  real, optional, intent(in) :: adaptdrho0 !< Reference density difference for stratification-dependent
2283  !! diffusion. [R ~> kg m-3]
2284 
2285  if (present(interp_scheme)) call set_interp_scheme(cs%interp_CS, interp_scheme)
2286  if (present(boundary_extrapolation)) call set_interp_extrap(cs%interp_CS, boundary_extrapolation)
2287 
2288  if (present(old_grid_weight)) then
2289  if (old_grid_weight<0. .or. old_grid_weight>1.) &
2290  call mom_error(fatal,'MOM_regridding, set_regrid_params: Weight is out side the range 0..1!')
2291  cs%old_grid_weight = old_grid_weight
2292  endif
2293  if (present(depth_of_time_filter_shallow)) cs%depth_of_time_filter_shallow = depth_of_time_filter_shallow
2294  if (present(depth_of_time_filter_deep)) cs%depth_of_time_filter_deep = depth_of_time_filter_deep
2295  if (present(depth_of_time_filter_shallow) .or. present(depth_of_time_filter_deep)) then
2296  if (cs%depth_of_time_filter_deep<cs%depth_of_time_filter_shallow) call mom_error(fatal,'MOM_regridding, '//&
2297  'set_regrid_params: depth_of_time_filter_deep<depth_of_time_filter_shallow!')
2298  endif
2299 
2300  if (present(min_thickness)) cs%min_thickness = min_thickness
2301  if (present(compress_fraction)) cs%compressibility_fraction = compress_fraction
2302  if (present(ref_pressure)) cs%ref_pressure = ref_pressure
2303  if (present(integrate_downward_for_e)) cs%integrate_downward_for_e = integrate_downward_for_e
2304  if (present(remap_answers_2018)) cs%remap_answers_2018 = remap_answers_2018
2305 
2306  select case (cs%regridding_scheme)
2307  case (regridding_zstar)
2308  if (present(min_thickness)) call set_zlike_params(cs%zlike_CS, min_thickness=min_thickness)
2309  case (regridding_sigma_shelf_zstar)
2310  if (present(min_thickness)) call set_zlike_params(cs%zlike_CS, min_thickness=min_thickness)
2311  case (regridding_sigma)
2312  if (present(min_thickness)) call set_sigma_params(cs%sigma_CS, min_thickness=min_thickness)
2313  case (regridding_rho)
2314  if (present(min_thickness)) call set_rho_params(cs%rho_CS, min_thickness=min_thickness)
2315  if (present(integrate_downward_for_e)) &
2316  call set_rho_params(cs%rho_CS, integrate_downward_for_e=integrate_downward_for_e)
2317  if (associated(cs%rho_CS) .and. (present(interp_scheme) .or. present(boundary_extrapolation))) &
2318  call set_rho_params(cs%rho_CS, interp_cs=cs%interp_CS)
2319  case (regridding_hycom1)
2320  if (associated(cs%hycom_CS) .and. (present(interp_scheme) .or. present(boundary_extrapolation))) &
2321  call set_hycom_params(cs%hycom_CS, interp_cs=cs%interp_CS)
2322  case (regridding_slight)
2323  if (present(min_thickness)) call set_slight_params(cs%slight_CS, min_thickness=min_thickness)
2324  if (present(dz_min_surface)) call set_slight_params(cs%slight_CS, dz_ml_min=dz_min_surface)
2325  if (present(nz_fixed_surface)) call set_slight_params(cs%slight_CS, nz_fixed_surface=nz_fixed_surface)
2326  if (present(rho_ml_avg_depth)) call set_slight_params(cs%slight_CS, rho_ml_avg_depth=rho_ml_avg_depth)
2327  if (present(nlay_ml_to_interior)) call set_slight_params(cs%slight_CS, nlay_ml_offset=nlay_ml_to_interior)
2328  if (present(fix_haloclines)) call set_slight_params(cs%slight_CS, fix_haloclines=fix_haloclines)
2329  if (present(halocline_filt_len)) call set_slight_params(cs%slight_CS, halocline_filter_length=halocline_filt_len)
2330  if (present(halocline_strat_tol)) call set_slight_params(cs%slight_CS, halocline_strat_tol=halocline_strat_tol)
2331  if (present(compress_fraction)) call set_slight_params(cs%slight_CS, compressibility_fraction=compress_fraction)
2332  if (associated(cs%slight_CS) .and. (present(interp_scheme) .or. present(boundary_extrapolation))) &
2333  call set_slight_params(cs%slight_CS, interp_cs=cs%interp_CS)
2334  case (regridding_adaptive)
2335  if (present(adapttimeratio)) call set_adapt_params(cs%adapt_CS, adapttimeratio=adapttimeratio)
2336  if (present(adaptzoom)) call set_adapt_params(cs%adapt_CS, adaptzoom=adaptzoom)
2337  if (present(adaptzoomcoeff)) call set_adapt_params(cs%adapt_CS, adaptzoomcoeff=adaptzoomcoeff)
2338  if (present(adaptbuoycoeff)) call set_adapt_params(cs%adapt_CS, adaptbuoycoeff=adaptbuoycoeff)
2339  if (present(adaptalpha)) call set_adapt_params(cs%adapt_CS, adaptalpha=adaptalpha)
2340  if (present(adaptdomin)) call set_adapt_params(cs%adapt_CS, adaptdomin=adaptdomin)
2341  if (present(adaptdrho0)) call set_adapt_params(cs%adapt_CS, adaptdrho0=adaptdrho0)
2342  end select
2343 
2344 end subroutine set_regrid_params
2345 
2346 !> Returns the number of levels/layers in the regridding control structure
2347 integer function get_regrid_size(CS)
2348  type(regridding_cs), intent(inout) :: cs !< Regridding control structure
2349 
2350  get_regrid_size = cs%nk
2351 
2352 end function get_regrid_size
2353 
2354 !> This returns a copy of the zlike_CS stored in the regridding control structure.
2355 function get_zlike_cs(CS)
2356  type(regridding_cs), intent(in) :: cs !< Regridding control structure
2357  type(zlike_cs) :: get_zlike_cs
2358 
2359  get_zlike_cs = cs%zlike_CS
2360 end function get_zlike_cs
2361 
2362 !> This returns a copy of the sigma_CS stored in the regridding control structure.
2363 function get_sigma_cs(CS)
2364  type(regridding_cs), intent(in) :: cs !< Regridding control structure
2365  type(sigma_cs) :: get_sigma_cs
2366 
2367  get_sigma_cs = cs%sigma_CS
2368 end function get_sigma_cs
2369 
2370 !> This returns a copy of the rho_CS stored in the regridding control structure.
2371 function get_rho_cs(CS)
2372  type(regridding_cs), intent(in) :: cs !< Regridding control structure
2373  type(rho_cs) :: get_rho_cs
2374 
2375  get_rho_cs = cs%rho_CS
2376 end function get_rho_cs
2377 
2378 !------------------------------------------------------------------------------
2379 !> Return coordinate-derived thicknesses for fixed coordinate systems
2380 function getstaticthickness( CS, SSH, depth )
2381  type(regridding_cs), intent(in) :: cs !< Regridding control structure
2382  real, intent(in) :: ssh !< The sea surface height, in the same units as depth
2383  real, intent(in) :: depth !< The maximum depth of the grid, often [Z ~> m]
2384  real, dimension(CS%nk) :: getstaticthickness !< The returned thicknesses in the units of depth
2385  ! Local
2386  integer :: k
2387  real :: z, dz
2388 
2389  select case ( cs%regridding_scheme )
2390  case ( regridding_zstar, regridding_sigma_shelf_zstar, regridding_hycom1, regridding_slight, regridding_adaptive )
2391  if (depth>0.) then
2392  z = ssh
2393  do k = 1, cs%nk
2394  dz = cs%coordinateResolution(k) * ( 1. + ssh/depth ) ! Nominal dz*
2395  dz = max(dz, 0.) ! Avoid negative incase ssh=-depth
2396  dz = min(dz, depth - z) ! Clip if below topography
2397  z = z + dz ! Bottom of layer
2398  getstaticthickness(k) = dz
2399  enddo
2400  else
2401  getstaticthickness(:) = 0. ! On land ...
2402  endif
2403  case ( regridding_sigma )
2404  getstaticthickness(:) = cs%coordinateResolution(:) * ( depth + ssh )
2405  case ( regridding_rho )
2406  getstaticthickness(:) = 0. ! Not applicable
2407  case ( regridding_arbitrary )
2408  getstaticthickness(:) = 0. ! Not applicable
2409  case default
2410  call mom_error(fatal,'MOM_regridding, getStaticThickness: '//&
2411  'Unknown regridding scheme selected!')
2412  end select ! type of grid
2413 
2414 end function getstaticthickness
2415 
2416 !> Parses a string and generates a dz(:) profile that goes like k**power.
2417 subroutine dz_function1( string, dz )
2418  character(len=*), intent(in) :: string !< String with list of parameters in form
2419  !! dz_min, H_total, power, precision
2420  real, dimension(:), intent(inout) :: dz !< Profile of nominal thicknesses
2421  ! Local variables
2422  integer :: nk, k
2423  real :: dz_min, power, prec, H_total
2424 
2425  nk = size(dz) ! Number of cells
2426  prec = -1024.
2427  read( string, *) dz_min, h_total, power, prec
2428  if (prec == -1024.) call mom_error(fatal,"dz_function1: "// &
2429  "Problem reading FNC1: string ="//trim(string))
2430  ! Create profile of ( dz - dz_min )
2431  do k = 1, nk
2432  dz(k) = (real(k-1)/real(nk-1))**power
2433  enddo
2434  dz(:) = ( h_total - real(nk) * dz_min ) * ( dz(:) / sum(dz) ) ! Rescale to so total is H_total
2435  dz(:) = anint( dz(:) / prec ) * prec ! Rounds to precision prec
2436  dz(:) = ( h_total - real(nk) * dz_min ) * ( dz(:) / sum(dz) ) ! Rescale to so total is H_total
2437  dz(:) = anint( dz(:) / prec ) * prec ! Rounds to precision prec
2438  dz(nk) = dz(nk) + ( h_total - sum( dz(:) + dz_min ) ) ! Adjust bottommost layer
2439  dz(:) = anint( dz(:) / prec ) * prec ! Rounds to precision prec
2440  dz(:) = dz(:) + dz_min ! Finally add in the constant dz_min
2441 
2442 end subroutine dz_function1
2443 
2444 !> Parses a string and generates a rho_target(:) profile with refined resolution downward
2445 !! and returns the number of levels
2446 integer function rho_function1( string, rho_target )
2447  character(len=*), intent(in) :: string !< String with list of parameters in form
2448  !! dz_min, H_total, power, precision
2449  real, dimension(:), allocatable, intent(inout) :: rho_target !< Profile of interface densities [kg m-3]
2450  ! Local variables
2451  integer :: nki, k, nk
2452  real :: ddx, dx, rho_1, rho_2, rho_3, drho, rho_4, drho_min
2453 
2454  read( string, *) nk, rho_1, rho_2, rho_3, drho, rho_4, drho_min
2455  allocate(rho_target(nk+1))
2456  nki = nk + 1 - 4 ! Number of interfaces minus 4 specified values
2457  rho_target(1) = rho_1
2458  rho_target(2) = rho_2
2459  dx = 0.
2460  do k = 0, nki
2461  ddx = max( drho_min, real(nki-k)/real(nki*nki) )
2462  dx = dx + ddx
2463  rho_target(3+k) = rho_3 + (2. * drho) * dx
2464  enddo
2465  rho_target(nki+4) = rho_4
2466 
2467  rho_function1 = nk
2468 
2469 end function rho_function1
2470 
2471 !> \namespace mom_regridding
2472 !!
2473 !! A vertical grid is defined solely by the cell thicknesses, \f$h\f$.
2474 !! Most calculations in this module start with the coordinate at the bottom
2475 !! of the column set to -depth, and use a increasing value of coordinate with
2476 !! decreasing k. This is consistent with the rest of MOM6 that uses position,
2477 !! \f$z\f$ which is a negative quantity for most of the ocean.
2478 !!
2479 !! A change in grid is define through a change in position of the interfaces:
2480 !! \f[
2481 !! z^n_{k+1/2} = z^{n-1}_{k+1/2} + \Delta z_{k+1/2}
2482 !! \f]
2483 !! with the positive upward coordinate convention
2484 !! \f[
2485 !! z_{k-1/2} = z_{k+1/2} + h_k
2486 !! \f]
2487 !! so that
2488 !! \f[
2489 !! h^n_k = h^{n-1}_k + ( \Delta z_{k-1/2} - \Delta z_{k+1/2} )
2490 !! \f]
2491 !!
2492 !! Original date of creation: 2008.06.09 by L. White
2493 
2494 end module mom_regridding
coord_slight::slight_cs
Control structure containing required parameters for the SLight coordinate.
Definition: coord_slight.F90:15
coord_rho
Regrid columns for the continuous isopycnal (rho) coordinate.
Definition: coord_rho.F90:2
mom_regridding::regridding_cs
Regridding control structure.
Definition: MOM_regridding.F90:45
mom_verticalgrid
Provides a transparent vertical ocean grid type and supporting routines.
Definition: MOM_verticalGrid.F90:2
regrid_consts::state_dependent
Returns true if the coordinate is dependent on the state density, returns false otherwise.
Definition: regrid_consts.F90:44
mom_variables::thermo_var_ptrs
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
Definition: MOM_variables.F90:80
coord_zlike::zlike_cs
Control structure containing required parameters for a z-like coordinate.
Definition: coord_zlike.F90:11
mom_string_functions
Handy functions for manipulating strings.
Definition: MOM_string_functions.F90:2
regrid_consts::coordinateunits
Returns a string with the coordinate units associated with the coordinate mode.
Definition: regrid_consts.F90:38
mom_eos
Provides subroutines for quantities specific to the equation of state.
Definition: MOM_EOS.F90:2
coord_rho::rho_cs
Control structure containing required parameters for the rho coordinate.
Definition: coord_rho.F90:14
mom_remapping::remapping_cs
Container for remapping parameters.
Definition: MOM_remapping.F90:24
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_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
regrid_consts
Contains constants for interpreting input parameters that control regridding.
Definition: regrid_consts.F90:2
mom_remapping
Provides column-wise vertical remapping functions.
Definition: MOM_remapping.F90:2
coord_hycom
Regrid columns for the HyCOM coordinate.
Definition: coord_hycom.F90:2
regrid_interp
Vertical interpolation for regridding.
Definition: regrid_interp.F90:2
mom_eos::eos_type
A control structure for the equation of state.
Definition: MOM_EOS.F90:108
mom_verticalgrid::verticalgrid_type
Describes the vertical ocean grid, including unit conversion factors.
Definition: MOM_verticalGrid.F90:24
coord_zlike
Regrid columns for a z-like coordinate (z-star, z-level)
Definition: coord_zlike.F90:2
mom_variables
Provides transparent structures with groups of MOM6 variables and supporting routines.
Definition: MOM_variables.F90:2
mom_io::mom_read_data
Read a data field from a file.
Definition: MOM_io.F90:74
coord_adapt::adapt_cs
Control structure for adaptive coordinates (coord_adapt).
Definition: coord_adapt.F90:17
mom_file_parser
The MOM6 facility to parse input files for runtime parameters.
Definition: MOM_file_parser.F90:2
coord_sigma::sigma_cs
Control structure containing required parameters for the sigma coordinate.
Definition: coord_sigma.F90:11
mom_unit_scaling
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Definition: MOM_unit_scaling.F90:2
regrid_interp::interp_cs_type
Control structure for regrid_interp module.
Definition: regrid_interp.F90:23
mom_io::file_exists
Indicate whether a file exists, perhaps with domain decomposition.
Definition: MOM_io.F90:68
mom_regridding
Generates vertical grids as part of the ALE algorithm.
Definition: MOM_regridding.F90:2
mom_file_parser::log_param
An overloaded interface to log the values of various types of parameters.
Definition: MOM_file_parser.F90:96
coord_sigma
Regrid columns for the sigma coordinate.
Definition: coord_sigma.F90:2
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2
coord_hycom::hycom_cs
Control structure containing required parameters for the HyCOM coordinate.
Definition: coord_hycom.F90:13
coord_slight
Regrid columns for the SLight coordinate.
Definition: coord_slight.F90:2
coord_adapt
Regrid columns for the adaptive coordinate.
Definition: coord_adapt.F90:2
mom_eos::calculate_density
Calculates density of sea water from T, S and P.
Definition: MOM_EOS.F90:68