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