18 use regrid_consts,
only : coordinatemode, default_coordinate_mode
21 use regrid_consts,
only : regridding_arbitrary, regridding_sigma_shelf_zstar
22 use regrid_consts,
only : regridding_hycom1, regridding_slight, regridding_adaptive
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
35 implicit none ;
private
37 #include <MOM_memory.h>
51 real,
dimension(:),
allocatable :: coordinateresolution
55 real :: coord_scale = 1.0
62 real,
dimension(:),
allocatable :: target_density
65 logical :: target_density_set = .false.
69 real,
dimension(:),
allocatable :: max_interface_depths
73 real,
dimension(:),
allocatable :: max_layer_thickness
79 integer :: regridding_scheme
88 real :: ref_pressure = 2.e7
94 real :: old_grid_weight = 0.
97 real :: depth_of_time_filter_shallow = 0.
100 real :: depth_of_time_filter_deep = 0.
104 real :: compressibility_fraction = 0.
108 logical :: set_maximum_depths = .false.
112 real :: max_depth_index_scale = 2.0
116 logical :: integrate_downward_for_e = .true.
120 logical :: remap_answers_2018 = .true.
132 public initialize_regridding, end_regridding, regridding_main
133 public inflate_vanished_layers_old, check_remapping_grid, check_grid_column
134 public set_regrid_params, get_regrid_size
135 public uniformresolution, setcoordinateresolution
136 public build_rho_column
137 public set_target_densities_from_gv, set_target_densities
138 public set_regrid_max_depths, set_regrid_max_thickness
139 public getcoordinateresolution, getcoordinateinterfaces
140 public getcoordinateunits, getcoordinateshortname, getstaticthickness
141 public default_coordinate_mode
142 public get_zlike_cs, get_sigma_cs, get_rho_cs
145 character(len=*),
parameter,
public :: regriddingcoordinatemodedoc = &
146 " LAYER - Isopycnal or stacked shallow water layers\n"//&
147 " ZSTAR, Z* - stretched geopotential z*\n"//&
148 " SIGMA_SHELF_ZSTAR - stretched geopotential z* ignoring shelf\n"//&
149 " SIGMA - terrain following coordinates\n"//&
150 " RHO - continuous isopycnal\n"//&
151 " HYCOM1 - HyCOM-like hybrid coordinate\n"//&
152 " SLIGHT - stretched coordinates above continuous isopycnal\n"//&
153 " ADAPTIVE - optimize for smooth neutral density surfaces"
156 character(len=*),
parameter,
public :: regriddinginterpschemedoc = &
157 " P1M_H2 (2nd-order accurate)\n"//&
158 " P1M_H4 (2nd-order accurate)\n"//&
159 " P1M_IH4 (2nd-order accurate)\n"//&
160 " PLM (2nd-order accurate)\n"//&
161 " PPM_H4 (3rd-order accurate)\n"//&
162 " PPM_IH4 (3rd-order accurate)\n"//&
163 " P3M_IH4IH3 (4th-order accurate)\n"//&
164 " P3M_IH6IH5 (4th-order accurate)\n"//&
165 " PQM_IH4IH3 (4th-order accurate)\n"//&
166 " PQM_IH6IH5 (5th-order accurate)"
169 character(len=*),
parameter,
public :: regriddingdefaultinterpscheme =
"P1M_H2"
171 logical,
parameter,
public :: regriddingdefaultboundaryextrapolation = .false.
173 real,
parameter,
public :: regriddingdefaultminthickness = 1.e-3
175 #undef __DO_SAFETY_CHECKS__
180 subroutine initialize_regridding(CS, GV, US, max_depth, param_file, mdl, coord_mode, param_prefix, param_suffix)
184 real,
intent(in) :: max_depth
186 character(len=*),
intent(in) :: mdl
187 character(len=*),
intent(in) :: coord_mode
188 character(len=*),
intent(in) :: param_prefix
190 character(len=*),
intent(in) :: param_suffix
194 character(len=80) :: string, string2, varname
195 character(len=40) :: coord_units, param_name, coord_res_param
196 character(len=200) :: inputdir, filename
197 character(len=320) :: message
198 character(len=12) :: expected_units
199 logical :: tmplogical, fix_haloclines, set_max, do_sum, main_parameters
200 logical :: coord_is_state_dependent, ierr
201 logical :: default_2018_answers, remap_answers_2018
202 real :: filt_len, strat_tol, index_scale, tmpreal, p_ref
203 real :: maximum_depth
204 real :: dz_fixed_sfc, rho_avg_depth, nlay_sfc_int
205 real :: adapttimeratio, adaptzoom, adaptzoomcoeff, adaptbuoycoeff, adaptalpha
207 integer :: nz_fixed_sfc, k, nzf(4)
208 real,
dimension(:),
allocatable :: dz
210 real,
dimension(:),
allocatable :: h_max
211 real,
dimension(:),
allocatable :: z_max
213 real,
dimension(:),
allocatable :: dz_max
215 real,
dimension(:),
allocatable :: rho_target
217 real,
dimension(40) :: woa09_dz = (/ 5., 10., 10., 15., 22.5, 25., 25., 25., &
218 37.5, 50., 50., 75., 100., 100., 100., 100., &
219 100., 100., 100., 100., 100., 100., 100., 175., &
220 250., 375., 500., 500., 500., 500., 500., 500., &
221 500., 500., 500., 500., 500., 500., 500., 500. /)
223 call get_param(param_file, mdl,
"INPUTDIR", inputdir, default=
".")
224 inputdir = slasher(inputdir)
226 main_parameters=.false.
227 if (len_trim(param_prefix)==0) main_parameters=.true.
228 if (main_parameters .and. len_trim(param_suffix)>0)
call mom_error(fatal,trim(mdl)//
', initialize_regridding: '// &
229 'Suffix provided without prefix for parameter names!')
232 cs%regridding_scheme = coordinatemode(coord_mode)
234 maximum_depth = us%Z_to_m*max_depth
236 if (main_parameters)
then
238 call get_param(param_file, mdl,
"REGRIDDING_COORDINATE_UNITS", coord_units, &
239 "Units of the regridding coordinate.", default=
coordinateunits(coord_mode))
244 if (coord_is_state_dependent)
then
245 if (main_parameters)
then
246 param_name =
"INTERPOLATION_SCHEME"
247 string2 = regriddingdefaultinterpscheme
249 param_name = trim(param_prefix)//
"_INTERP_SCHEME_"//trim(param_suffix)
252 call get_param(param_file, mdl,
"INTERPOLATION_SCHEME", string, &
253 "This sets the interpolation scheme to use to "//&
254 "determine the new grid. These parameters are "//&
255 "only relevant when REGRIDDING_COORDINATE_MODE is "//&
256 "set to a function of state. Otherwise, it is not "//&
257 "used. It can be one of the following schemes: "//&
258 trim(regriddinginterpschemedoc), default=trim(string2))
259 call set_regrid_params(cs, interp_scheme=string)
261 call get_param(param_file, mdl,
"DEFAULT_2018_ANSWERS", default_2018_answers, &
262 "This sets the default value for the various _2018_ANSWERS parameters.", &
264 call get_param(param_file, mdl,
"REMAPPING_2018_ANSWERS", remap_answers_2018, &
265 "If true, use the order of arithmetic and expressions that recover the "//&
266 "answers from the end of 2018. Otherwise, use updated and more robust "//&
267 "forms of the same expressions.", default=default_2018_answers)
268 call set_regrid_params(cs, remap_answers_2018=remap_answers_2018)
271 if (main_parameters .and. coord_is_state_dependent)
then
272 call get_param(param_file, mdl,
"BOUNDARY_EXTRAPOLATION", tmplogical, &
273 "When defined, a proper high-order reconstruction "//&
274 "scheme is used within boundary cells rather "//&
275 "than PCM. E.g., if PPM is used for remapping, a "//&
276 "PPM reconstruction will also be used within "//&
277 "boundary cells.", default=regriddingdefaultboundaryextrapolation)
278 call set_regrid_params(cs, boundary_extrapolation=tmplogical)
280 call set_regrid_params(cs, boundary_extrapolation=.false.)
284 if (main_parameters)
then
285 param_name =
"ALE_COORDINATE_CONFIG"
286 coord_res_param =
"ALE_RESOLUTION"
289 param_name = trim(param_prefix)//
"_DEF_"//trim(param_suffix)
290 coord_res_param = trim(param_prefix)//
"_RES_"//trim(param_suffix)
292 if (maximum_depth>3000.) string2=
'WOA09'
294 call get_param(param_file, mdl, param_name, string, &
295 "Determines how to specify the coordinate "//&
296 "resolution. Valid options are:\n"//&
297 " PARAM - use the vector-parameter "//trim(coord_res_param)//
"\n"//&
298 " UNIFORM[:N] - uniformly distributed\n"//&
299 " FILE:string - read from a file. The string specifies\n"//&
300 " the filename and variable name, separated\n"//&
301 " by a comma or space, e.g. FILE:lev.nc,dz\n"//&
302 " or FILE:lev.nc,interfaces=zw\n"//&
303 " WOA09[:N] - the WOA09 vertical grid (approximately)\n"//&
304 " FNC1:string - FNC1:dz_min,H_total,power,precision\n"//&
305 " HYBRID:string - read from a file. The string specifies\n"//&
306 " the filename and two variable names, separated\n"//&
307 " by a comma or space, for sigma-2 and dz. e.g.\n"//&
308 " HYBRID:vgrid.nc,sigma2,dz",&
309 default=trim(string2))
310 message =
"The distribution of vertical resolution for the target\n"//&
311 "grid used for Eulerian-like coordinates. For example,\n"//&
312 "in z-coordinate mode, the parameter is a list of level\n"//&
313 "thicknesses (in m). In sigma-coordinate mode, the list\n"//&
314 "is of non-dimensional fractions of the water column."
315 if (index(trim(string),
'UNIFORM')==1)
then
316 if (len_trim(string)==7)
then
318 tmpreal = maximum_depth
319 elseif (index(trim(string),
'UNIFORM:')==1 .and. len_trim(string)>8)
then
321 ke = extract_integer(string(9:len_trim(string)),
'',1)
322 tmpreal = extract_real(string(9:len_trim(string)),
',',2,missing_value=maximum_depth)
324 call mom_error(fatal,trim(mdl)//
', initialize_regridding: '// &
325 'Unable to interpret "'//trim(string)//
'".')
328 dz(:) = uniformresolution(ke, coord_mode, tmpreal, &
329 us%R_to_kg_m3*(gv%Rlay(1) + 0.5*(gv%Rlay(1)-gv%Rlay(min(2,ke)))), &
330 us%R_to_kg_m3*(gv%Rlay(ke) + 0.5*(gv%Rlay(ke)-gv%Rlay(max(ke-1,1)))) )
331 if (main_parameters)
call log_param(param_file, mdl,
"!"//coord_res_param, dz, &
332 trim(message), units=trim(coord_units))
333 elseif (trim(string)==
'PARAM')
then
337 call get_param(param_file, mdl, coord_res_param, dz, &
338 trim(message), units=trim(coord_units), fail_if_missing=.true.)
339 elseif (index(trim(string),
'FILE:')==1)
then
342 if (string(6:6)==
'.' .or. string(6:6)==
'/')
then
344 filename = trim( extractword(trim(string(6:80)), 1) )
347 filename = trim(inputdir) // trim( extractword(trim(string(6:80)), 1) )
349 if (.not.
file_exists(filename))
call mom_error(fatal,trim(mdl)//
", initialize_regridding: "// &
350 "Specified file not found: Looking for '"//trim(filename)//
"' ("//trim(string)//
")")
352 varname = trim( extractword(trim(string(6:)), 2) )
353 if (len_trim(varname)==0)
then
354 if (field_exists(filename,
'dz')) then; varname =
'dz'
355 elseif (field_exists(filename,
'dsigma')) then; varname =
'dsigma'
356 elseif (field_exists(filename,
'ztest')) then; varname =
'ztest'
357 else ;
call mom_error(fatal,trim(mdl)//
", initialize_regridding: "// &
358 "Coordinate variable not specified and none could be guessed.")
364 if (cs%regridding_scheme == regridding_sigma)
then
365 expected_units =
'nondim'
366 elseif (cs%regridding_scheme == regridding_rho)
then
367 expected_units =
'kg m-3'
369 expected_units =
'meters'
371 if (index(trim(varname),
'interfaces=')==1)
then
372 varname=trim(varname(12:))
373 call check_grid_def(filename, varname, expected_units, message, ierr)
374 if (ierr)
call mom_error(fatal,trim(mdl)//
", initialize_regridding: "//&
375 "Unsupported format in grid definition '"//trim(filename)//
"'. Error message "//trim(message))
376 call field_size(trim(filename), trim(varname), nzf)
378 if (ke < 1)
call mom_error(fatal, trim(mdl)//
" initialize_regridding via Var "//trim(varname)//&
379 "in FILE "//trim(filename)//
" requires at least 2 target interface values.")
380 if (cs%regridding_scheme == regridding_rho)
then
381 allocate(rho_target(ke+1))
382 call mom_read_data(trim(filename), trim(varname), rho_target)
385 allocate(z_max(ke+1))
387 dz(:) = abs(z_max(1:ke) - z_max(2:ke+1))
392 call field_size(trim(filename), trim(varname), nzf)
397 if (main_parameters .and. (ke/=gv%ke))
then
398 call mom_error(fatal,trim(mdl)//
', initialize_regridding: '// &
399 'Mismatch in number of model levels and "'//trim(string)//
'".')
401 if (main_parameters)
call log_param(param_file, mdl,
"!"//coord_res_param, dz, &
403 elseif (index(trim(string),
'FNC1:')==1)
then
404 ke = gv%ke;
allocate(dz(ke))
405 call dz_function1( trim(string(6:)), dz )
406 if (main_parameters)
call log_param(param_file, mdl,
"!"//coord_res_param, dz, &
408 elseif (index(trim(string),
'RFNC1:')==1)
then
410 ke = rho_function1( trim(string(7:)), rho_target )
411 elseif (index(trim(string),
'HYBRID:')==1)
then
412 ke = gv%ke;
allocate(dz(ke))
414 allocate(rho_target(ke+1))
415 filename = trim( extractword(trim(string(8:)), 1) )
416 if (filename(1:1)/=
'.' .and. filename(1:1)/=
'/') filename = trim(inputdir) // trim( filename )
417 if (.not.
file_exists(filename))
call mom_error(fatal,trim(mdl)//
", initialize_regridding: HYBRID "// &
418 "Specified file not found: Looking for '"//trim(filename)//
"' ("//trim(string)//
")")
419 varname = trim( extractword(trim(string(8:)), 2) )
420 if (.not. field_exists(filename,varname))
call mom_error(fatal,trim(mdl)//
", initialize_regridding: HYBRID "// &
421 "Specified field not found: Looking for '"//trim(varname)//
"' ("//trim(string)//
")")
422 call mom_read_data(trim(filename), trim(varname), rho_target)
423 varname = trim( extractword(trim(string(8:)), 3) )
424 if (varname(1:5) ==
'FNC1:')
then
425 call dz_function1( trim(string((index(trim(string),
'FNC1:')+5):)), dz )
427 if (.not. field_exists(filename,varname))
call mom_error(fatal,trim(mdl)//
", initialize_regridding: HYBRID "// &
428 "Specified field not found: Looking for '"//trim(varname)//
"' ("//trim(string)//
")")
431 if (main_parameters)
then
432 call log_param(param_file, mdl,
"!"//coord_res_param, dz, &
434 call log_param(param_file, mdl,
"!TARGET_DENSITIES", rho_target, &
435 'HYBRID target densities for interfaces', units=
coordinateunits(coord_mode))
437 elseif (index(trim(string),
'WOA09')==1)
then
438 if (len_trim(string)==5)
then
439 tmpreal = 0. ; ke = 0
440 do while (tmpreal<maximum_depth)
442 tmpreal = tmpreal + woa09_dz(ke)
444 elseif (index(trim(string),
'WOA09:')==1)
then
445 if (len_trim(string)==6)
call mom_error(fatal,trim(mdl)//
', initialize_regridding: '// &
446 'Expected string of form "WOA09:N" but got "'//trim(string)//
'".')
447 ke = extract_integer(string(7:len_trim(string)),
'',1)
449 if (ke>40 .or. ke<1)
call mom_error(fatal,trim(mdl)//
', initialize_regridding: '// &
450 'For "WOA05:N" N must 0<N<41 but got "'//trim(string)//
'".')
452 dz(1:ke) = woa09_dz(1:ke)
454 call mom_error(fatal,trim(mdl)//
", initialize_regridding: "// &
455 "Unrecognized coordinate configuration"//trim(string))
458 if (main_parameters)
then
460 if (coordinatemode(coord_mode) == regridding_zstar .or. &
461 coordinatemode(coord_mode) == regridding_hycom1 .or. &
462 coordinatemode(coord_mode) == regridding_slight .or. &
463 coordinatemode(coord_mode) == regridding_adaptive)
then
465 tmpreal = sum( dz(:) )
466 if (tmpreal < maximum_depth)
then
467 dz(ke) = dz(ke) + ( maximum_depth - tmpreal )
468 elseif (tmpreal > maximum_depth)
then
469 if ( dz(ke) + ( maximum_depth - tmpreal ) > 0. )
then
470 dz(ke) = dz(ke) + ( maximum_depth - tmpreal )
472 call mom_error(fatal,trim(mdl)//
", initialize_regridding: "// &
473 "MAXIMUM_DEPTH was too shallow to adjust bottom layer of DZ!"//trim(string))
482 allocate( cs%coordinateResolution(cs%nk) ); cs%coordinateResolution(:) = -1.e30
485 allocate( cs%target_density(cs%nk+1) ); cs%target_density(:) = -1.e30*us%kg_m3_to_R
488 if (
allocated(dz))
then
489 if (coordinatemode(coord_mode) == regridding_sigma)
then
490 call setcoordinateresolution(dz, cs, scale=1.0)
491 elseif (coordinatemode(coord_mode) == regridding_rho)
then
492 call setcoordinateresolution(dz, cs, scale=us%kg_m3_to_R)
493 cs%coord_scale = us%R_to_kg_m3
494 elseif (coordinatemode(coord_mode) == regridding_adaptive)
then
495 call setcoordinateresolution(dz, cs, scale=gv%m_to_H)
496 cs%coord_scale = gv%H_to_m
498 call setcoordinateresolution(dz, cs, scale=us%m_to_Z)
499 cs%coord_scale = us%Z_to_m
503 if (
allocated(rho_target))
then
504 call set_target_densities(cs, us%kg_m3_to_R*rho_target)
505 deallocate(rho_target)
508 elseif (coordinatemode(coord_mode) == regridding_rho)
then
509 call set_target_densities_from_gv(gv, us, cs)
510 call log_param(param_file, mdl,
"!TARGET_DENSITIES", us%R_to_kg_m3*cs%target_density(:), &
511 'RHO target densities for interfaces', units=
coordinateunits(coord_mode))
515 call initcoord(cs, gv, us, coord_mode)
517 if (main_parameters .and. coord_is_state_dependent)
then
518 call get_param(param_file, mdl,
"P_REF", p_ref, &
519 "The pressure that is used for calculating the coordinate "//&
520 "density. (1 Pa = 1e4 dbar, so 2e7 is commonly used.) "//&
521 "This is only used if USE_EOS and ENABLE_THERMODYNAMICS are true.", &
522 units=
"Pa", default=2.0e7, scale=us%kg_m3_to_R*us%m_s_to_L_T**2)
523 call get_param(param_file, mdl,
"REGRID_COMPRESSIBILITY_FRACTION", tmpreal, &
524 "When interpolating potential density profiles we can add "//&
525 "some artificial compressibility solely to make homogeneous "//&
526 "regions appear stratified.", units=
"nondim", default=0.)
527 call set_regrid_params(cs, compress_fraction=tmpreal, ref_pressure=p_ref)
530 if (main_parameters)
then
531 call get_param(param_file, mdl,
"MIN_THICKNESS", tmpreal, &
532 "When regridding, this is the minimum layer "//&
533 "thickness allowed.", units=
"m", scale=gv%m_to_H, &
534 default=regriddingdefaultminthickness )
535 call set_regrid_params(cs, min_thickness=tmpreal)
537 call set_regrid_params(cs, min_thickness=0.)
540 if (coordinatemode(coord_mode) == regridding_slight)
then
542 call get_param(param_file, mdl,
"SLIGHT_DZ_SURFACE", dz_fixed_sfc, &
543 "The nominal thickness of fixed thickness near-surface "//&
544 "layers with the SLight coordinate.", units=
"m", default=1.0, scale=gv%m_to_H)
545 call get_param(param_file, mdl,
"SLIGHT_NZ_SURFACE_FIXED", nz_fixed_sfc, &
546 "The number of fixed-depth surface layers with the SLight "//&
547 "coordinate.", units=
"nondimensional", default=2)
548 call get_param(param_file, mdl,
"SLIGHT_SURFACE_AVG_DEPTH", rho_avg_depth, &
549 "The thickness of the surface region over which to average "//&
550 "when calculating the density to use to define the interior "//&
551 "with the SLight coordinate.", units=
"m", default=1.0, scale=gv%m_to_H)
552 call get_param(param_file, mdl,
"SLIGHT_NLAY_TO_INTERIOR", nlay_sfc_int, &
553 "The number of layers to offset the surface density when "//&
554 "defining where the interior ocean starts with SLight.", &
555 units=
"nondimensional", default=2.0)
556 call get_param(param_file, mdl,
"SLIGHT_FIX_HALOCLINES", fix_haloclines, &
557 "If true, identify regions above the reference pressure "//&
558 "where the reference pressure systematically underestimates "//&
559 "the stratification and use this in the definition of the "//&
560 "interior with the SLight coordinate.", default=.false.)
562 call set_regrid_params(cs, dz_min_surface=dz_fixed_sfc, &
563 nz_fixed_surface=nz_fixed_sfc, rho_ml_avg_depth=rho_avg_depth, &
564 nlay_ml_to_interior=nlay_sfc_int, fix_haloclines=fix_haloclines)
565 if (fix_haloclines)
then
567 call get_param(param_file, mdl,
"HALOCLINE_FILTER_LENGTH", filt_len, &
568 "A length scale over which to smooth the temperature and "//&
569 "salinity before identifying erroneously unstable haloclines.", &
570 units=
"m", default=2.0, scale=gv%m_to_H)
571 call get_param(param_file, mdl,
"HALOCLINE_STRAT_TOL", strat_tol, &
572 "A tolerance for the ratio of the stratification of the "//&
573 "apparent coordinate stratification to the actual value "//&
574 "that is used to identify erroneously unstable haloclines. "//&
575 "This ratio is 1 when they are equal, and sensible values "//&
576 "are between 0 and 0.5.", units=
"nondimensional", default=0.2)
577 call set_regrid_params(cs, halocline_filt_len=filt_len, &
578 halocline_strat_tol=strat_tol)
583 if (coordinatemode(coord_mode) == regridding_adaptive)
then
584 call get_param(param_file, mdl,
"ADAPT_TIME_RATIO", adapttimeratio, &
585 "Ratio of ALE timestep to grid timescale.", units=
"nondim", default=1.0e-1)
586 call get_param(param_file, mdl,
"ADAPT_ZOOM_DEPTH", adaptzoom, &
587 "Depth of near-surface zooming region.", units=
"m", default=200.0, scale=gv%m_to_H)
588 call get_param(param_file, mdl,
"ADAPT_ZOOM_COEFF", adaptzoomcoeff, &
589 "Coefficient of near-surface zooming diffusivity.", units=
"nondim", default=0.2)
590 call get_param(param_file, mdl,
"ADAPT_BUOY_COEFF", adaptbuoycoeff, &
591 "Coefficient of buoyancy diffusivity.", units=
"nondim", default=0.8)
592 call get_param(param_file, mdl,
"ADAPT_ALPHA", adaptalpha, &
593 "Scaling on optimization tendency.", units=
"nondim", default=1.0)
594 call get_param(param_file, mdl,
"ADAPT_DO_MIN_DEPTH", tmplogical, &
595 "If true, make a HyCOM-like mixed layer by preventing interfaces "//&
596 "from being shallower than the depths specified by the regridding coordinate.", &
598 call get_param(param_file, mdl,
"ADAPT_DRHO0", adaptdrho0, &
599 "Reference density difference for stratification-dependent diffusion.", &
600 units=
"kg m-3", default=0.5, scale=us%kg_m3_to_R)
602 call set_regrid_params(cs, adapttimeratio=adapttimeratio, adaptzoom=adaptzoom, &
603 adaptzoomcoeff=adaptzoomcoeff, adaptbuoycoeff=adaptbuoycoeff, adaptalpha=adaptalpha, &
604 adaptdomin=tmplogical, adaptdrho0=adaptdrho0)
607 if (main_parameters .and. coord_is_state_dependent)
then
608 call get_param(param_file, mdl,
"MAXIMUM_INT_DEPTH_CONFIG", string, &
609 "Determines how to specify the maximum interface depths.\n"//&
610 "Valid options are:\n"//&
611 " NONE - there are no maximum interface depths\n"//&
612 " PARAM - use the vector-parameter MAXIMUM_INTERFACE_DEPTHS\n"//&
613 " FILE:string - read from a file. The string specifies\n"//&
614 " the filename and variable name, separated\n"//&
615 " by a comma or space, e.g. FILE:lev.nc,Z\n"//&
616 " FNC1:string - FNC1:dz_min,H_total,power,precision",&
618 message =
"The list of maximum depths for each interface."
619 allocate(z_max(ke+1))
621 if ( trim(string) ==
"NONE")
then
623 elseif ( trim(string) ==
"PARAM")
then
624 call get_param(param_file, mdl,
"MAXIMUM_INTERFACE_DEPTHS", z_max, &
625 trim(message), units=
"m", scale=gv%m_to_H, fail_if_missing=.true.)
626 call set_regrid_max_depths(cs, z_max)
627 elseif (index(trim(string),
'FILE:')==1)
then
628 if (string(6:6)==
'.' .or. string(6:6)==
'/')
then
630 filename = trim( extractword(trim(string(6:80)), 1) )
633 filename = trim(inputdir) // trim( extractword(trim(string(6:80)), 1) )
635 if (.not.
file_exists(filename))
call mom_error(fatal,trim(mdl)//
", initialize_regridding: "// &
636 "Specified file not found: Looking for '"//trim(filename)//
"' ("//trim(string)//
")")
639 varname = trim( extractword(trim(string(6:)), 2) )
640 if (.not. field_exists(filename,varname))
call mom_error(fatal,trim(mdl)//
", initialize_regridding: "// &
641 "Specified field not found: Looking for '"//trim(varname)//
"' ("//trim(string)//
")")
642 if (len_trim(varname)==0)
then
643 if (field_exists(filename,
'z_max')) then; varname =
'z_max'
644 elseif (field_exists(filename,
'dz')) then; varname =
'dz' ; do_sum = .true.
645 elseif (field_exists(filename,
'dz_max')) then; varname =
'dz_max' ; do_sum = .true.
646 else ;
call mom_error(fatal,trim(mdl)//
", initialize_regridding: "// &
647 "MAXIMUM_INT_DEPTHS variable not specified and none could be guessed.")
652 z_max(1) = 0.0 ;
do k=1,ke ; z_max(k+1) = z_max(k) + dz_max(k) ;
enddo
656 call log_param(param_file, mdl,
"!MAXIMUM_INT_DEPTHS", z_max, &
658 call set_regrid_max_depths(cs, z_max, gv%m_to_H)
659 elseif (index(trim(string),
'FNC1:')==1)
then
660 call dz_function1( trim(string(6:)), dz_max )
661 if ((coordinatemode(coord_mode) == regridding_slight) .and. &
662 (dz_fixed_sfc > 0.0))
then
663 do k=1,nz_fixed_sfc ; dz_max(k) = dz_fixed_sfc ;
enddo
665 z_max(1) = 0.0 ;
do k=1,ke ; z_max(k+1) = z_max(k) + dz_max(k) ;
enddo
666 call log_param(param_file, mdl,
"!MAXIMUM_INT_DEPTHS", z_max, &
668 call set_regrid_max_depths(cs, z_max, gv%m_to_H)
670 call mom_error(fatal,trim(mdl)//
", initialize_regridding: "// &
671 "Unrecognized MAXIMUM_INT_DEPTH_CONFIG "//trim(string))
678 call get_param(param_file, mdl,
"MAX_LAYER_THICKNESS_CONFIG", string, &
679 "Determines how to specify the maximum layer thicknesses.\n"//&
680 "Valid options are:\n"//&
681 " NONE - there are no maximum layer thicknesses\n"//&
682 " PARAM - use the vector-parameter MAX_LAYER_THICKNESS\n"//&
683 " FILE:string - read from a file. The string specifies\n"//&
684 " the filename and variable name, separated\n"//&
685 " by a comma or space, e.g. FILE:lev.nc,Z\n"//&
686 " FNC1:string - FNC1:dz_min,H_total,power,precision",&
688 message =
"The list of maximum thickness for each layer."
690 if ( trim(string) ==
"NONE")
then
692 elseif ( trim(string) ==
"PARAM")
then
693 call get_param(param_file, mdl,
"MAX_LAYER_THICKNESS", h_max, &
694 trim(message), units=
"m", fail_if_missing=.true., scale=gv%m_to_H)
695 call set_regrid_max_thickness(cs, h_max)
696 elseif (index(trim(string),
'FILE:')==1)
then
697 if (string(6:6)==
'.' .or. string(6:6)==
'/')
then
699 filename = trim( extractword(trim(string(6:80)), 1) )
702 filename = trim(inputdir) // trim( extractword(trim(string(6:80)), 1) )
704 if (.not.
file_exists(filename))
call mom_error(fatal,trim(mdl)//
", initialize_regridding: "// &
705 "Specified file not found: Looking for '"//trim(filename)//
"' ("//trim(string)//
")")
707 varname = trim( extractword(trim(string(6:)), 2) )
708 if (.not. field_exists(filename,varname))
call mom_error(fatal,trim(mdl)//
", initialize_regridding: "// &
709 "Specified field not found: Looking for '"//trim(varname)//
"' ("//trim(string)//
")")
710 if (len_trim(varname)==0)
then
711 if (field_exists(filename,
'h_max')) then; varname =
'h_max'
712 elseif (field_exists(filename,
'dz_max')) then; varname =
'dz_max'
713 else ;
call mom_error(fatal,trim(mdl)//
", initialize_regridding: "// &
714 "MAXIMUM_INT_DEPTHS variable not specified and none could be guessed.")
718 call log_param(param_file, mdl,
"!MAX_LAYER_THICKNESS", h_max, &
720 call set_regrid_max_thickness(cs, h_max, gv%m_to_H)
721 elseif (index(trim(string),
'FNC1:')==1)
then
722 call dz_function1( trim(string(6:)), h_max )
723 call log_param(param_file, mdl,
"!MAX_LAYER_THICKNESS", h_max, &
725 call set_regrid_max_thickness(cs, h_max, gv%m_to_H)
727 call mom_error(fatal,trim(mdl)//
", initialize_regridding: "// &
728 "Unrecognized MAX_LAYER_THICKNESS_CONFIG "//trim(string))
733 if (
allocated(dz))
deallocate(dz)
734 end subroutine initialize_regridding
737 subroutine check_grid_def(filename, varname, expected_units, msg, ierr)
738 character(len=*),
intent(in) :: filename
739 character(len=*),
intent(in) :: varname
740 character(len=*),
intent(in) :: expected_units
741 character(len=*),
intent(inout) :: msg
742 logical,
intent(out) :: ierr
744 character (len=200) :: units, long_name
745 integer :: ncid, status, intid, vid
749 status = nf90_open(trim(filename), nf90_nowrite, ncid)
750 if (status /= nf90_noerr)
then
752 msg =
'File not found: '//trim(filename)
756 status = nf90_inq_varid(ncid, trim(varname), vid)
757 if (status /= nf90_noerr)
then
759 msg =
'Var not found: '//trim(varname)
763 status = nf90_get_att(ncid, vid,
"units", units)
764 if (status /= nf90_noerr)
then
766 msg =
'Attribute not found: units'
772 do i=1,len_trim(units)
773 if (units(i:i) == char(0)) units(i:i) =
" "
776 if (trim(units) /= trim(expected_units))
then
777 if (trim(expected_units) ==
"meters")
then
778 if (trim(units) /=
"m")
then
787 msg =
'Units incorrect: '//trim(units)//
' /= '//trim(expected_units)
790 end subroutine check_grid_def
793 subroutine end_regridding(CS)
796 if (
associated(cs%zlike_CS))
call end_coord_zlike(cs%zlike_CS)
797 if (
associated(cs%sigma_CS))
call end_coord_sigma(cs%sigma_CS)
798 if (
associated(cs%rho_CS))
call end_coord_rho(cs%rho_CS)
799 if (
associated(cs%hycom_CS))
call end_coord_hycom(cs%hycom_CS)
800 if (
associated(cs%slight_CS))
call end_coord_slight(cs%slight_CS)
801 if (
associated(cs%adapt_CS))
call end_coord_adapt(cs%adapt_CS)
803 deallocate( cs%coordinateResolution )
804 if (
allocated(cs%target_density))
deallocate( cs%target_density )
805 if (
allocated(cs%max_interface_depths) )
deallocate( cs%max_interface_depths )
806 if (
allocated(cs%max_layer_thickness) )
deallocate( cs%max_layer_thickness )
808 end subroutine end_regridding
812 subroutine regridding_main( remapCS, CS, G, GV, h, tv, h_new, dzInterface, frac_shelf_h, conv_adjust)
833 type(ocean_grid_type),
intent(in) :: g
835 real,
dimension(SZI_(G),SZJ_(G), SZK_(GV)),
intent(inout) :: h
838 real,
dimension(SZI_(G),SZJ_(G), CS%nk),
intent(inout) :: h_new
839 real,
dimension(SZI_(G),SZJ_(G), CS%nk+1),
intent(inout) :: dzinterface
840 real,
dimension(:,:),
optional,
pointer :: frac_shelf_h
841 logical,
optional,
intent(in ) :: conv_adjust
843 real :: trickgnucompiler
844 logical :: use_ice_shelf
845 logical :: do_convective_adjustment
847 do_convective_adjustment = .true.
848 if (
present(conv_adjust)) do_convective_adjustment = conv_adjust
850 use_ice_shelf = .false.
851 if (
present(frac_shelf_h))
then
852 if (
associated(frac_shelf_h)) use_ice_shelf = .true.
855 select case ( cs%regridding_scheme )
857 case ( regridding_zstar )
858 if (use_ice_shelf)
then
859 call build_zstar_grid( cs, g, gv, h, dzinterface, frac_shelf_h )
861 call build_zstar_grid( cs, g, gv, h, dzinterface )
863 call calc_h_new_by_dz(cs, g, gv, h, dzinterface, h_new)
865 case ( regridding_sigma_shelf_zstar)
866 call build_zstar_grid( cs, g, gv, h, dzinterface )
867 call calc_h_new_by_dz(cs, g, gv, h, dzinterface, h_new)
869 case ( regridding_sigma )
870 call build_sigma_grid( cs, g, gv, h, dzinterface )
871 call calc_h_new_by_dz(cs, g, gv, h, dzinterface, h_new)
873 case ( regridding_rho )
874 if (do_convective_adjustment)
call convective_adjustment(g, gv, h, tv)
875 call build_rho_grid( g, gv, g%US, h, tv, dzinterface, remapcs, cs )
876 call calc_h_new_by_dz(cs, g, gv, h, dzinterface, h_new)
878 case ( regridding_arbitrary )
879 call build_grid_arbitrary( g, gv, h, dzinterface, trickgnucompiler, cs )
880 call calc_h_new_by_dz(cs, g, gv, h, dzinterface, h_new)
882 case ( regridding_hycom1 )
883 call build_grid_hycom1( g, gv, g%US, h, tv, h_new, dzinterface, cs )
885 case ( regridding_slight )
886 call build_grid_slight( g, gv, g%US, h, tv, dzinterface, cs )
887 call calc_h_new_by_dz(cs, g, gv, h, dzinterface, h_new)
889 case ( regridding_adaptive )
890 call build_grid_adaptive(g, gv, g%US, h, tv, dzinterface, remapcs, cs)
891 call calc_h_new_by_dz(cs, g, gv, h, dzinterface, h_new)
894 call mom_error(fatal,
'MOM_regridding, regridding_main: '//&
895 'Unknown regridding scheme selected!')
899 #ifdef __DO_SAFETY_CHECKS__
900 call check_remapping_grid(g, gv, h, dzinterface,
'in regridding_main')
903 end subroutine regridding_main
906 subroutine calc_h_new_by_dz(CS, G, GV, h, dzInterface, h_new)
908 type(ocean_grid_type),
intent(in) :: G
910 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)),
intent(in) :: h
911 real,
dimension(SZI_(G),SZJ_(G),CS%nk+1),
intent(in) :: dzInterface
912 real,
dimension(SZI_(G),SZJ_(G),CS%nk),
intent(inout) :: h_new
914 integer :: i, j, k, nki
916 nki = min(cs%nk, gv%ke)
919 do j = g%jsc-1,g%jec+1
920 do i = g%isc-1,g%iec+1
921 if (g%mask2dT(i,j)>0.)
then
923 h_new(i,j,k) = max( 0., h(i,j,k) + ( dzinterface(i,j,k) - dzinterface(i,j,k+1) ) )
925 if (cs%nk > gv%ke)
then
927 h_new(i,j,k) = max( 0., dzinterface(i,j,k) - dzinterface(i,j,k+1) )
931 h_new(i,j,1:nki) = h(i,j,1:nki)
932 if (cs%nk > gv%ke) h_new(i,j,nki+1:cs%nk) = 0.
938 end subroutine calc_h_new_by_dz
941 subroutine check_remapping_grid( G, GV, h, dzInterface, msg )
942 type(ocean_grid_type),
intent(in) :: g
944 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)),
intent(in) :: h
945 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)+1),
intent(in) :: dzinterface
947 character(len=*),
intent(in) :: msg
952 do j = g%jsc-1,g%jec+1 ;
do i = g%isc-1,g%iec+1
953 if (g%mask2dT(i,j)>0.)
call check_grid_column( gv%ke, gv%Z_to_H*g%bathyT(i,j), h(i,j,:), dzinterface(i,j,:), msg )
956 end subroutine check_remapping_grid
959 subroutine check_grid_column( nk, depth, h, dzInterface, msg )
960 integer,
intent(in) :: nk
961 real,
intent(in) :: depth
962 real,
dimension(nk),
intent(in) :: h
963 real,
dimension(nk+1),
intent(in) :: dzinterface
964 character(len=*),
intent(in) :: msg
967 real :: eps, total_h_old, total_h_new, h_new, z_old, z_new
969 eps =1. ; eps = epsilon(eps)
974 total_h_old = total_h_old + h(k)
979 if (depth == 0.) z_old = - total_h_old
983 z_new = z_old + dzinterface(k)
984 h_new = h(k) + ( dzinterface(k) - dzinterface(k+1) )
986 write(0,*)
'k,h,hnew=',k,h(k),h_new
987 write(0,*)
'dzI(k+1),dzI(k)=',dzinterface(k+1),dzinterface(k)
988 call mom_error( fatal,
'MOM_regridding, check_grid_column: '//&
989 'Negative layer thickness implied by re-gridding, '//trim(msg))
991 total_h_new = total_h_new + h_new
996 if (abs(total_h_new-total_h_old)>real(nk-1)*0.5*(total_h_old+total_h_new)*eps)
then
999 write(0,*)
'k,h,hnew=',k,h(k),h(k)+(dzinterface(k)-dzinterface(k+1))
1001 write(0,*)
'Hold,Hnew,Hnew-Hold=',total_h_old,total_h_new,total_h_new-total_h_old
1002 write(0,*)
'eps,(n)/2*eps*H=',eps,real(nk-1)*0.5*(total_h_old+total_h_new)*eps
1003 call mom_error( fatal,
'MOM_regridding, check_grid_column: '//&
1004 'Re-gridding did NOT conserve total thickness to within roundoff '//trim(msg))
1008 if (dzinterface(1) /= 0.)
call mom_error( fatal, &
1009 'MOM_regridding, check_grid_column: Non-zero dzInterface at surface! '//trim(msg))
1010 if (dzinterface(nk+1) /= 0.)
call mom_error( fatal, &
1011 'MOM_regridding, check_grid_column: Non-zero dzInterface at bottom! '//trim(msg))
1013 end subroutine check_grid_column
1020 subroutine filtered_grid_motion( CS, nk, z_old, z_new, dz_g )
1022 integer,
intent(in) :: nk
1023 real,
dimension(nk+1),
intent(in) :: z_old
1024 real,
dimension(CS%nk+1),
intent(in) :: z_new
1025 real,
dimension(CS%nk+1),
intent(inout) :: dz_g
1028 real :: dz_tgt, zr1, z_old_k
1029 real :: Aq, Bq, dz0, z0, F0
1030 real :: zs, zd, dzwt, Idzwt
1032 real :: Int_zs, Int_zd, dInt_zs_zd
1034 real,
dimension(nk+1) :: z_act
1036 logical :: debug = .false.
1039 if ((z_old(nk+1) - z_old(1)) * (z_new(cs%nk+1) - z_new(1)) < 0.0)
then
1040 call mom_error(fatal,
"filtered_grid_motion: z_old and z_new use different sign conventions.")
1041 elseif ((z_old(nk+1) - z_old(1)) * (z_new(cs%nk+1) - z_new(1)) == 0.0)
then
1043 do k=1,cs%nk+1 ; dz_g(k) = 0.0 ;
enddo ;
return
1044 elseif ((z_old(nk+1) - z_old(1)) + (z_new(cs%nk+1) - z_new(1)) > 0.0)
then
1052 if (sgn*(z_new(k)-z_new(k-1)) < -5e-16*(abs(z_new(k))+abs(z_new(k-1))) ) &
1053 call mom_error(fatal,
"filtered_grid_motion: z_new is tangled.")
1056 if (sgn*(z_old(k)-z_old(k-1)) < -5e-16*(abs(z_old(k))+abs(z_old(k-1))) ) &
1057 call mom_error(fatal,
"filtered_grid_motion: z_old is tangled.")
1062 zs = cs%depth_of_time_filter_shallow
1063 zd = cs%depth_of_time_filter_deep
1064 wtd = 1.0 - cs%old_grid_weight
1068 idzwt = 0.0 ;
if (abs(zd - zs) > 0.0) idzwt = 1.0 / (zd - zs)
1069 dint_zs_zd = 0.5*(1.0 + iwtd) * (zd - zs)
1070 aq = 0.5*(iwtd - 1.0)
1075 if (k<=nk+1) z_old_k = z_old(k)
1077 dz_tgt = sgn*(z_new(k) - z_old_k)
1078 zr1 = sgn*(z_old_k - z_old(1))
1082 if ((zr1 > zd) .and. (zr1 + wtd * dz_tgt > zd))
then
1083 dz_g(k) = sgn * wtd * dz_tgt
1084 elseif ((zr1 < zs) .and. (zr1 + dz_tgt < zs))
then
1085 dz_g(k) = sgn * dz_tgt
1094 int_zd = iwtd*(zd - zr1)
1095 int_zs = int_zd - dint_zs_zd
1096 elseif (zr1 <= zs)
then
1098 int_zd = dint_zs_zd + (zs - zr1)
1101 int_zd = (zd - zr1) * (iwtd*(0.5*(zd+zr1) - zs) + 0.5*(zd - zr1)) * idzwt
1102 int_zs = (zs - zr1) * (0.5*iwtd * ((zr1 - zs)) + (zd - 0.5*(zr1+zs))) * idzwt
1106 if (dz_tgt >= int_zd)
then
1107 dz_g(k) = sgn * ((zd-zr1) + wtd*(dz_tgt - int_zd))
1108 elseif (dz_tgt <= int_zs)
then
1109 dz_g(k) = sgn * ((zs-zr1) + (dz_tgt - int_zs))
1115 if (zr1 <= zs)
then ; dz0 = zs-zr1 ; z0 = zs ; f0 = dz_tgt - int_zs
1116 elseif (zr1 >= zd)
then ; dz0 = zd-zr1 ; z0 = zd ; f0 = dz_tgt - int_zd
1117 else ; dz0 = 0.0 ; z0 = zr1 ; f0 = dz_tgt ;
endif
1119 bq = (dzwt + 2.0*aq*(z0-zs))
1122 dz_g(k) = sgn * (dz0 + 2.0*f0*dzwt / (bq + sqrt(bq**2 + 4.0*aq*f0*dzwt) ))
1144 if (k<=nk+1) z_old_k = z_old(k)
1145 z_act(k) = z_old_k + dz_g(k)
1148 if (sgn*((z_act(k))-z_act(k-1)) < -1e-15*(abs(z_act(k))+abs(z_act(k-1))) ) &
1149 call mom_error(fatal,
"filtered_grid_motion: z_output is tangled.")
1153 end subroutine filtered_grid_motion
1158 subroutine build_zstar_grid( CS, G, GV, h, dzInterface, frac_shelf_h)
1162 type(ocean_grid_type),
intent(in) :: G
1164 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)),
intent(in) :: h
1165 real,
dimension(SZI_(G),SZJ_(G), CS%nk+1),
intent(inout) :: dzInterface
1167 real,
dimension(:,:),
optional,
pointer :: frac_shelf_h
1169 real :: nominalDepth, totalThickness, dh
1170 real,
dimension(SZK_(GV)+1) :: zOld, zNew
1171 integer :: i, j, k, nz
1172 logical :: ice_shelf
1176 if (
present(frac_shelf_h))
then
1177 if (
associated(frac_shelf_h)) ice_shelf = .true.
1182 do j = g%jsc-1,g%jec+1
1183 do i = g%isc-1,g%iec+1
1185 if (g%mask2dT(i,j)==0.)
then
1186 dzinterface(i,j,:) = 0.
1191 nominaldepth = g%bathyT(i,j)*gv%Z_to_H
1194 totalthickness = 0.0
1196 totalthickness = totalthickness + h(i,j,k)
1199 zold(nz+1) = - nominaldepth
1201 zold(k) = zold(k+1) + h(i,j,k)
1205 if (frac_shelf_h(i,j) > 0.)
then
1206 call build_zstar_column(cs%zlike_CS, nominaldepth, totalthickness, znew, &
1207 z_rigid_top = totalthickness-nominaldepth, &
1208 eta_orig=zold(1), zscale=gv%Z_to_H)
1210 call build_zstar_column(cs%zlike_CS, nominaldepth, totalthickness, &
1211 znew, zscale=gv%Z_to_H)
1214 call build_zstar_column(cs%zlike_CS, nominaldepth, totalthickness, &
1215 znew, zscale=gv%Z_to_H)
1219 call filtered_grid_motion( cs, nz, zold, znew, dzinterface(i,j,:) )
1221 #ifdef __DO_SAFETY_CHECKS__
1222 dh=max(nominaldepth,totalthickness)
1223 if (abs(znew(1)-zold(1))>(nz-1)*0.5*epsilon(dh)*dh)
then
1224 write(0,*)
'min_thickness=',cs%min_thickness
1225 write(0,*)
'nominalDepth=',nominaldepth,
'totalThickness=',totalthickness
1226 write(0,*)
'dzInterface(1) = ',dzinterface(i,j,1),epsilon(dh),nz
1228 write(0,*) k,zold(k),znew(k)
1231 write(0,*) k,h(i,j,k),znew(k)-znew(k+1),cs%coordinateResolution(k)
1233 call mom_error( fatal, &
1234 'MOM_regridding, build_zstar_grid(): top surface has moved!!!' )
1238 call adjust_interface_motion( cs, nz, h(i,j,:), dzinterface(i,j,:) )
1243 end subroutine build_zstar_grid
1248 subroutine build_sigma_grid( CS, G, GV, h, dzInterface )
1258 type(ocean_grid_type),
intent(in) :: G
1260 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)),
intent(in) :: h
1261 real,
dimension(SZI_(G),SZJ_(G), CS%nk+1),
intent(inout) :: dzInterface
1267 real :: nominalDepth, totalThickness, dh
1268 real,
dimension(SZK_(GV)+1) :: zOld, zNew
1272 do i = g%isc-1,g%iec+1
1273 do j = g%jsc-1,g%jec+1
1275 if (g%mask2dT(i,j)==0.)
then
1276 dzinterface(i,j,:) = 0.
1281 nominaldepth = g%bathyT(i,j)*gv%Z_to_H
1284 totalthickness = 0.0
1286 totalthickness = totalthickness + h(i,j,k)
1289 call build_sigma_column(cs%sigma_CS, nominaldepth, totalthickness, znew)
1292 zold(nz+1) = -nominaldepth
1294 zold(k) = zold(k+1) + h(i, j, k)
1297 call filtered_grid_motion( cs, nz, zold, znew, dzinterface(i,j,:) )
1299 #ifdef __DO_SAFETY_CHECKS__
1300 dh=max(nominaldepth,totalthickness)
1301 if (abs(znew(1)-zold(1))>(cs%nk-1)*0.5*epsilon(dh)*dh)
then
1302 write(0,*)
'min_thickness=',cs%min_thickness
1303 write(0,*)
'nominalDepth=',nominaldepth,
'totalThickness=',totalthickness
1304 write(0,*)
'dzInterface(1) = ',dzinterface(i,j,1),epsilon(dh),nz,cs%nk
1306 write(0,*) k,zold(k),znew(k)
1309 write(0,*) k,h(i,j,k),znew(k)-znew(k+1),totalthickness*cs%coordinateResolution(k),cs%coordinateResolution(k)
1311 call mom_error( fatal, &
1312 'MOM_regridding, build_sigma_grid: top surface has moved!!!' )
1314 dzinterface(i,j,1) = 0.
1315 dzinterface(i,j,cs%nk+1) = 0.
1321 end subroutine build_sigma_grid
1327 subroutine build_rho_grid( G, GV, US, h, tv, dzInterface, remapCS, CS )
1344 type(ocean_grid_type),
intent(in) :: G
1347 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)),
intent(in) :: h
1349 real,
dimension(SZI_(G),SZJ_(G), SZK_(GV)+1),
intent(inout) :: dzInterface
1357 real :: nominalDepth
1358 real,
dimension(SZK_(GV)+1) :: zOld, zNew
1359 real :: h_neglect, h_neglect_edge
1360 #ifdef __DO_SAFETY_CHECKS__
1361 real :: totalThickness
1365 if (.not.cs%remap_answers_2018)
then
1366 h_neglect = gv%H_subroundoff ; h_neglect_edge = gv%H_subroundoff
1367 elseif (gv%Boussinesq)
then
1368 h_neglect = gv%m_to_H*1.0e-30 ; h_neglect_edge = gv%m_to_H*1.0e-10
1370 h_neglect = gv%kg_m2_to_H*1.0e-30 ; h_neglect_edge = gv%kg_m2_to_H*1.0e-10
1375 if (.not.cs%target_density_set)
call mom_error(fatal,
"build_rho_grid: "//&
1376 "Target densities must be set before build_rho_grid is called.")
1379 do j = g%jsc-1,g%jec+1
1380 do i = g%isc-1,g%iec+1
1382 if (g%mask2dT(i,j)==0.)
then
1383 dzinterface(i,j,:) = 0.
1389 nominaldepth = g%bathyT(i,j)*gv%Z_to_H
1391 call build_rho_column(cs%rho_CS, nz, nominaldepth, h(i, j, :), &
1392 tv%T(i, j, :), tv%S(i, j, :), tv%eqn_of_state, znew, &
1393 h_neglect=h_neglect, h_neglect_edge=h_neglect_edge)
1395 if (cs%integrate_downward_for_e)
then
1398 zold(k+1) = zold(k) - h(i,j,k)
1402 zold(nz+1) = - nominaldepth
1404 zold(k) = zold(k+1) + h(i,j,k)
1409 call filtered_grid_motion( cs, nz, zold, znew, dzinterface(i,j,:) )
1411 #ifdef __DO_SAFETY_CHECKS__
1413 if (znew(k) > zold(1))
then
1414 write(0,*)
'zOld=',zold
1415 write(0,*)
'zNew=',znew
1416 call mom_error( fatal,
'MOM_regridding, build_rho_grid: '//&
1417 'interior interface above surface!' )
1419 if (znew(k) > znew(k-1))
then
1420 write(0,*)
'zOld=',zold
1421 write(0,*)
'zNew=',znew
1422 call mom_error( fatal,
'MOM_regridding, build_rho_grid: '//&
1423 'interior interfaces cross!' )
1427 totalthickness = 0.0
1429 totalthickness = totalthickness + h(i,j,k)
1432 dh=max(nominaldepth,totalthickness)
1433 if (abs(znew(1)-zold(1))>(nz-1)*0.5*epsilon(dh)*dh)
then
1434 write(0,*)
'min_thickness=',cs%min_thickness
1435 write(0,*)
'nominalDepth=',nominaldepth,
'totalThickness=',totalthickness
1436 write(0,*)
'zNew(1)-zOld(1) = ',znew(1)-zold(1),epsilon(dh),nz
1438 write(0,*) k,zold(k),znew(k)
1441 write(0,*) k,h(i,j,k),znew(k)-znew(k+1)
1443 call mom_error( fatal, &
1444 'MOM_regridding, build_rho_grid: top surface has moved!!!' )
1451 end subroutine build_rho_grid
1460 subroutine build_grid_hycom1( G, GV, US, h, tv, h_new, dzInterface, CS )
1461 type(ocean_grid_type),
intent(in) :: G
1464 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)),
intent(in) :: h
1467 real,
dimension(SZI_(G),SZJ_(G),CS%nk),
intent(inout) :: h_new
1468 real,
dimension(SZI_(G),SZJ_(G),CS%nk+1),
intent(inout) :: dzInterface
1471 real,
dimension(SZK_(GV)+1) :: z_col
1472 real,
dimension(CS%nk+1) :: z_col_new
1473 real,
dimension(SZK_(GV)+1) :: dz_col
1474 real,
dimension(SZK_(GV)) :: p_col
1476 integer :: i, j, k, nki
1478 real :: h_neglect, h_neglect_edge
1480 if (.not.cs%remap_answers_2018)
then
1481 h_neglect = gv%H_subroundoff ; h_neglect_edge = gv%H_subroundoff
1482 elseif (gv%Boussinesq)
then
1483 h_neglect = gv%m_to_H*1.0e-30 ; h_neglect_edge = gv%m_to_H*1.0e-10
1485 h_neglect = gv%kg_m2_to_H*1.0e-30 ; h_neglect_edge = gv%kg_m2_to_H*1.0e-10
1488 if (.not.cs%target_density_set)
call mom_error(fatal,
"build_grid_HyCOM1 : "//&
1489 "Target densities must be set before build_grid_HyCOM1 is called.")
1491 nki = min(gv%ke, cs%nk)
1494 do j = g%jsc-1,g%jec+1 ;
do i = g%isc-1,g%iec+1
1495 if (g%mask2dT(i,j)>0.)
then
1497 depth = g%bathyT(i,j) * gv%Z_to_H
1501 z_col(k+1) = z_col(k) + h(i,j,k)
1502 p_col(k) = tv%P_Ref + cs%compressibility_fraction * &
1503 ( 0.5 * ( z_col(k) + z_col(k+1) ) * (gv%H_to_RZ*gv%g_Earth) - tv%P_Ref )
1506 call build_hycom1_column(cs%hycom_CS, tv%eqn_of_state, gv%ke, depth, &
1507 h(i,j,:), tv%T(i,j,:), tv%S(i,j,:), p_col, &
1508 z_col, z_col_new, zscale=gv%Z_to_H, &
1509 h_neglect=h_neglect, h_neglect_edge=h_neglect_edge)
1512 call filtered_grid_motion( cs, gv%ke, z_col, z_col_new, dz_col )
1515 dz_col(:) = -dz_col(:)
1516 call adjust_interface_motion( cs, gv%ke, h(i,j,:), dz_col(:) )
1518 dzinterface(i,j,1:nki+1) = dz_col(1:nki+1)
1519 if (nki<cs%nk) dzinterface(i,j,nki+2:cs%nk+1) = 0.
1522 dzinterface(i,j,:) = 0.
1526 call calc_h_new_by_dz(cs, g, gv, h, dzinterface, h_new)
1528 end subroutine build_grid_hycom1
1532 subroutine build_grid_adaptive(G, GV, US, h, tv, dzInterface, remapCS, CS)
1533 type(ocean_grid_type),
intent(in) :: G
1536 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)),
intent(in) :: h
1539 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)+1),
intent(inout) :: dzInterface
1545 integer :: i, j, k, nz
1547 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)+1) :: tInt, sInt
1550 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)+1) :: zInt
1551 real,
dimension(SZK_(GV)+1) :: zNext
1559 do k = 2, nz ;
do j = g%jsc-2,g%jec+2 ;
do i = g%isc-2,g%iec+2
1560 tint(i,j,k) = 0.5 * (tv%T(i,j,k-1) + tv%T(i,j,k))
1561 sint(i,j,k) = 0.5 * (tv%S(i,j,k-1) + tv%S(i,j,k))
1562 zint(i,j,k) = zint(i,j,k-1) + h(i,j,k-1)
1563 enddo ;
enddo ;
enddo
1567 tint(:,:,1) = tv%T(:,:,1) ; tint(:,:,nz+1) = tv%T(:,:,nz)
1568 sint(:,:,1) = tv%S(:,:,1) ; sint(:,:,nz+1) = tv%S(:,:,nz)
1571 zint(:,:,nz+1) = zint(:,:,nz) + h(:,:,nz)
1575 do j = g%jsc-1,g%jec+1 ;
do i = g%isc-1,g%iec+1
1576 if (g%mask2dT(i,j) < 0.5)
then
1577 dzinterface(i,j,:) = 0.
1581 call build_adapt_column(cs%adapt_CS, g, gv, us, tv, i, j, zint, tint, sint, h, znext)
1583 call filtered_grid_motion(cs, nz, zint(i,j,:), znext, dzinterface(i,j,:))
1585 do k = 1, nz+1 ; dzinterface(i,j,k) = -dzinterface(i,j,k) ;
enddo
1586 call adjust_interface_motion(cs, nz, h(i,j,:), dzinterface(i,j,:))
1588 end subroutine build_grid_adaptive
1599 subroutine build_grid_slight(G, GV, US, h, tv, dzInterface, CS)
1600 type(ocean_grid_type),
intent(in) :: G
1603 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)),
intent(in) :: h
1605 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)+1),
intent(inout) :: dzInterface
1608 real,
dimension(SZK_(GV)+1) :: z_col
1609 real,
dimension(SZK_(GV)+1) :: z_col_new
1610 real,
dimension(SZK_(GV)+1) :: dz_col
1611 real,
dimension(SZK_(GV)) :: p_col
1613 integer :: i, j, k, nz
1614 real :: h_neglect, h_neglect_edge
1616 if (.not.cs%remap_answers_2018)
then
1617 h_neglect = gv%H_subroundoff ; h_neglect_edge = gv%H_subroundoff
1618 elseif (gv%Boussinesq)
then
1619 h_neglect = gv%m_to_H*1.0e-30 ; h_neglect_edge = gv%m_to_H*1.0e-10
1621 h_neglect = gv%kg_m2_to_H*1.0e-30 ; h_neglect_edge = gv%kg_m2_to_H*1.0e-10
1626 if (.not.cs%target_density_set)
call mom_error(fatal,
"build_grid_SLight : "//&
1627 "Target densities must be set before build_grid_SLight is called.")
1630 do j = g%jsc-1,g%jec+1 ;
do i = g%isc-1,g%iec+1
1631 if (g%mask2dT(i,j)>0.)
then
1633 depth = g%bathyT(i,j) * gv%Z_to_H
1636 z_col(k+1) = z_col(k) + h(i,j,k)
1637 p_col(k) = tv%P_Ref + cs%compressibility_fraction * &
1638 ( 0.5 * ( z_col(k) + z_col(k+1) ) * (gv%H_to_RZ*gv%g_Earth) - tv%P_Ref )
1641 call build_slight_column(cs%slight_CS, tv%eqn_of_state, gv%H_to_RZ*gv%g_Earth, &
1642 gv%H_subroundoff, nz, depth, h(i, j, :), &
1643 tv%T(i, j, :), tv%S(i, j, :), p_col, z_col, z_col_new, &
1644 h_neglect=h_neglect, h_neglect_edge=h_neglect_edge)
1647 call filtered_grid_motion( cs, nz, z_col, z_col_new, dz_col )
1648 do k=1,nz+1 ; dzinterface(i,j,k) = -dz_col(k) ;
enddo
1649 #ifdef __DO_SAFETY_CHECKS__
1650 if (dzinterface(i,j,1) /= 0.) stop
'build_grid_SLight: Surface moved?!'
1651 if (dzinterface(i,j,nz+1) /= 0.) stop
'build_grid_SLight: Bottom moved?!'
1655 call adjust_interface_motion( cs, nz, h(i,j,:), dzinterface(i,j,:) )
1658 dzinterface(i,j,:) = 0.
1662 end subroutine build_grid_slight
1665 subroutine adjust_interface_motion( CS, nk, h_old, dz_int )
1667 integer,
intent(in) :: nk
1668 real,
dimension(nk),
intent(in) :: h_old
1669 real,
dimension(CS%nk+1),
intent(inout) :: dz_int
1672 real :: h_new, eps, h_total, h_err
1674 eps = 1. ; eps = epsilon(eps)
1676 h_total = 0. ; h_err = 0.
1677 do k = 1, min(cs%nk,nk)
1678 h_total = h_total + h_old(k)
1679 h_err = h_err + max( h_old(k), abs(dz_int(k)), abs(dz_int(k+1)) )*eps
1680 h_new = h_old(k) + ( dz_int(k) - dz_int(k+1) )
1681 if (h_new < -3.0*h_err)
then
1682 write(0,*)
'h<0 at k=',k,
'h_old=',h_old(k), &
1683 'wup=',dz_int(k),
'wdn=',dz_int(k+1),
'dw_dz=',dz_int(k) - dz_int(k+1), &
1684 'h_new=',h_new,
'h_err=',h_err
1685 call mom_error( fatal,
'MOM_regridding: adjust_interface_motion() - '//&
1686 'implied h<0 is larger than roundoff!')
1691 h_err = h_err + max( abs(dz_int(k)), abs(dz_int(k+1)) )*eps
1692 h_new = ( dz_int(k) - dz_int(k+1) )
1693 if (h_new < -3.0*h_err)
then
1694 write(0,*)
'h<0 at k=',k,
'h_old was empty',&
1695 'wup=',dz_int(k),
'wdn=',dz_int(k+1),
'dw_dz=',dz_int(k) - dz_int(k+1), &
1696 'h_new=',h_new,
'h_err=',h_err
1697 call mom_error( fatal,
'MOM_regridding: adjust_interface_motion() - '//&
1698 'implied h<0 is larger than roundoff!')
1702 do k = min(cs%nk,nk),2,-1
1703 h_new = h_old(k) + ( dz_int(k) - dz_int(k+1) )
1704 if (h_new<cs%min_thickness) &
1705 dz_int(k) = ( dz_int(k+1) - h_old(k) ) + cs%min_thickness
1706 h_new = h_old(k) + ( dz_int(k) - dz_int(k+1) )
1708 dz_int(k) = ( 1. - eps ) * ( dz_int(k+1) - h_old(k) )
1709 h_new = h_old(k) + ( dz_int(k) - dz_int(k+1) )
1711 write(0,*)
'h<0 at k=',k,
'h_old=',h_old(k), &
1712 'wup=',dz_int(k),
'wdn=',dz_int(k+1),
'dw_dz=',dz_int(k) - dz_int(k+1), &
1714 stop
'Still did not work!'
1715 call mom_error( fatal,
'MOM_regridding: adjust_interface_motion() - '//&
1716 'Repeated adjustment for roundoff h<0 failed!')
1721 end subroutine adjust_interface_motion
1726 subroutine build_grid_arbitrary( G, GV, h, dzInterface, h_new, CS )
1732 type(ocean_grid_type),
intent(in) :: G
1734 real,
dimension(SZI_(G),SZJ_(G), SZK_(GV)),
intent(in) :: h
1735 real,
dimension(SZI_(G),SZJ_(G), SZK_(GV)+1),
intent(inout) :: dzInterface
1737 real,
intent(inout) :: h_new
1743 real :: z_inter(SZK_(GV)+1)
1744 real :: total_height
1749 real :: x1, y1, x2, y2
1753 max_depth = g%max_depth*gv%Z_to_H
1755 do j = g%jsc-1,g%jec+1
1756 do i = g%isc-1,g%iec+1
1759 local_depth = g%bathyT(i,j)*gv%Z_to_H
1764 total_height = total_height + h(i,j,k)
1767 eta = total_height - local_depth
1770 delta_h = (max_depth + eta) / nz
1775 z_inter(k+1) = z_inter(k) - delta_h
1780 x1 = 0.35; y1 = 0.45; x2 = 0.65; y2 = 0.55
1782 x = - ( z_inter(k) - eta ) / max_depth
1786 elseif ( (x > x1 ) .and. ( x < x2 ))
then
1787 t = y1 + (y2-y1) * (x-x1) / (x2-x1)
1789 t = y2 + (1.0-y2) * (x-x2) / (1.0-x2)
1792 z_inter(k) = -t * max_depth + eta
1797 z_inter(nz+1) = - local_depth
1801 if ( z_inter(k) < (z_inter(k+1) + cs%min_thickness) )
then
1802 z_inter(k) = z_inter(k+1) + cs%min_thickness
1808 dzinterface(i,j,1) = 0.
1811 dzinterface(i,j,k) = z_inter(k) - x
1813 dzinterface(i,j,nz+1) = 0.
1821 end subroutine build_grid_arbitrary
1828 subroutine inflate_vanished_layers_old( CS, G, GV, h )
1839 type(ocean_grid_type),
intent(in) :: g
1841 real,
dimension(SZI_(G),SZJ_(G), SZK_(GV)),
intent(inout) :: h
1847 do i = g%isc-1,g%iec+1
1848 do j = g%jsc-1,g%jec+1
1855 call old_inflate_layers_1d( cs%min_thickness, gv%ke, htmp )
1865 end subroutine inflate_vanished_layers_old
1869 subroutine convective_adjustment(G, GV, h, tv)
1870 type(ocean_grid_type),
intent(in) :: G
1872 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
1886 logical :: stratified
1887 real,
dimension(GV%ke) :: p_col, densities
1892 do j = g%jsc-1,g%jec+1 ;
do i = g%isc-1,g%iec+1
1895 call calculate_density( tv%T(i,j,:), tv%S(i,j,:), p_col, densities, tv%eqn_of_state)
1902 t0 = tv%T(i,j,k) ; t1 = tv%T(i,j,k+1)
1903 s0 = tv%S(i,j,k) ; s1 = tv%S(i,j,k+1)
1904 r0 = densities(k) ; r1 = densities(k+1)
1905 h0 = h(i,j,k) ; h1 = h(i,j,k+1)
1910 tv%T(i,j,k) = t1 ; tv%T(i,j,k+1) = t0
1911 tv%S(i,j,k) = s1 ; tv%S(i,j,k+1) = s0
1912 h(i,j,k) = h1 ; h(i,j,k+1) = h0
1914 call calculate_density( tv%T(i,j,k), tv%S(i,j,k), p_col(k), densities(k), tv%eqn_of_state)
1916 densities(k+1), tv%eqn_of_state )
1917 stratified = .false.
1921 if ( stratified )
exit
1926 end subroutine convective_adjustment
1931 function uniformresolution(nk,coordMode,maxDepth,rhoLight,rhoHeavy)
1936 integer,
intent(in) :: nk
1937 character(len=*),
intent(in) :: coordmode
1940 real,
intent(in) :: maxdepth
1941 real,
intent(in) :: rholight
1942 real,
intent(in) :: rhoheavy
1944 real :: uniformresolution(nk)
1949 scheme = coordinatemode(coordmode)
1950 select case ( scheme )
1952 case ( regridding_zstar, regridding_hycom1, regridding_slight, regridding_sigma_shelf_zstar, &
1953 regridding_adaptive )
1954 uniformresolution(:) = maxdepth / real(nk)
1956 case ( regridding_rho )
1957 uniformresolution(:) = (rhoheavy - rholight) / real(nk)
1959 case ( regridding_sigma )
1960 uniformresolution(:) = 1. / real(nk)
1963 call mom_error(fatal,
"MOM_regridding, uniformResolution: "//&
1964 "Unrecognized choice for coordinate mode ("//trim(coordmode)//
").")
1968 end function uniformresolution
1972 subroutine initcoord(CS, GV, US, coord_mode)
1974 character(len=*),
intent(in) :: coord_mode
1980 select case (coordinatemode(coord_mode))
1981 case (regridding_zstar)
1982 call init_coord_zlike(cs%zlike_CS, cs%nk, cs%coordinateResolution)
1983 case (regridding_sigma_shelf_zstar)
1984 call init_coord_zlike(cs%zlike_CS, cs%nk, cs%coordinateResolution)
1985 case (regridding_sigma)
1986 call init_coord_sigma(cs%sigma_CS, cs%nk, cs%coordinateResolution)
1987 case (regridding_rho)
1988 call init_coord_rho(cs%rho_CS, cs%nk, cs%ref_pressure, cs%target_density, cs%interp_CS)
1989 case (regridding_hycom1)
1990 call init_coord_hycom(cs%hycom_CS, cs%nk, cs%coordinateResolution, cs%target_density, &
1992 case (regridding_slight)
1993 call init_coord_slight(cs%slight_CS, cs%nk, cs%ref_pressure, cs%target_density, &
1994 cs%interp_CS, gv%m_to_H)
1995 case (regridding_adaptive)
1996 call init_coord_adapt(cs%adapt_CS, cs%nk, cs%coordinateResolution, gv%m_to_H, us%kg_m3_to_R)
1998 end subroutine initcoord
2002 subroutine setcoordinateresolution( dz, CS, scale )
2003 real,
dimension(:),
intent(in) :: dz
2005 real,
optional,
intent(in) :: scale
2007 if (
size(dz)/=cs%nk)
call mom_error( fatal, &
2008 'setCoordinateResolution: inconsistent number of levels' )
2010 if (
present(scale))
then
2011 cs%coordinateResolution(:) = scale*dz(:)
2013 cs%coordinateResolution(:) = dz(:)
2016 end subroutine setcoordinateresolution
2019 subroutine set_target_densities_from_gv( GV, US, CS )
2028 cs%target_density(1) = 0.0
2029 cs%target_density(2) = 2.0*gv%Rlay(1)
2031 cs%target_density(1) = (gv%Rlay(1) + 0.5*(gv%Rlay(1)-gv%Rlay(2)))
2032 cs%target_density(nz+1) = (gv%Rlay(nz) + 0.5*(gv%Rlay(nz)-gv%Rlay(nz-1)))
2034 cs%target_density(k) = cs%target_density(k-1) + cs%coordinateResolution(k)
2037 cs%target_density_set = .true.
2039 end subroutine set_target_densities_from_gv
2042 subroutine set_target_densities( CS, rho_int )
2044 real,
dimension(CS%nk+1),
intent(in) :: rho_int
2046 if (
size(cs%target_density)/=
size(rho_int))
then
2047 call mom_error(fatal,
"set_target_densities inconsistent args!")
2050 cs%target_density(:) = rho_int(:)
2051 cs%target_density_set = .true.
2053 end subroutine set_target_densities
2056 subroutine set_regrid_max_depths( CS, max_depths, units_to_H )
2058 real,
dimension(CS%nk+1),
intent(in) :: max_depths
2059 real,
optional,
intent(in) :: units_to_h
2064 if (.not.
allocated(cs%max_interface_depths))
allocate(cs%max_interface_depths(1:cs%nk+1))
2066 val_to_h = 1.0 ;
if (
present(units_to_h)) val_to_h = units_to_h
2067 if (max_depths(cs%nk+1) < max_depths(1)) val_to_h = -1.0*val_to_h
2070 if (max_depths(cs%nk+1) < max_depths(1))
then
2071 do k=1,cs%nk ;
if (max_depths(k+1) > max_depths(k)) &
2072 call mom_error(fatal,
"Unordered list of maximum depths sent to set_regrid_max_depths!")
2075 do k=1,cs%nk ;
if (max_depths(k+1) < max_depths(k)) &
2076 call mom_error(fatal,
"Unordered list of maximum depths sent to set_regrid_max_depths.")
2081 cs%max_interface_depths(k) = val_to_h * max_depths(k)
2085 select case (cs%regridding_scheme)
2086 case (regridding_hycom1)
2087 call set_hycom_params(cs%hycom_CS, max_interface_depths=cs%max_interface_depths)
2088 case (regridding_slight)
2089 call set_slight_params(cs%slight_CS, max_interface_depths=cs%max_interface_depths)
2091 end subroutine set_regrid_max_depths
2094 subroutine set_regrid_max_thickness( CS, max_h, units_to_H )
2096 real,
dimension(CS%nk+1),
intent(in) :: max_h
2097 real,
optional,
intent(in) :: units_to_h
2102 if (.not.
allocated(cs%max_layer_thickness))
allocate(cs%max_layer_thickness(1:cs%nk))
2104 val_to_h = 1.0 ;
if (
present( units_to_h)) val_to_h = units_to_h
2107 cs%max_layer_thickness(k) = val_to_h * max_h(k)
2111 select case (cs%regridding_scheme)
2112 case (regridding_hycom1)
2113 call set_hycom_params(cs%hycom_CS, max_layer_thickness=cs%max_layer_thickness)
2114 case (regridding_slight)
2115 call set_slight_params(cs%slight_CS, max_layer_thickness=cs%max_layer_thickness)
2117 end subroutine set_regrid_max_thickness
2122 function getcoordinateresolution( CS, undo_scaling )
2124 logical,
optional,
intent(in) :: undo_scaling
2126 real,
dimension(CS%nk) :: getcoordinateresolution
2129 unscale = .false. ;
if (
present(undo_scaling)) unscale = undo_scaling
2132 getcoordinateresolution(:) = cs%coord_scale * cs%coordinateResolution(:)
2134 getcoordinateresolution(:) = cs%coordinateResolution(:)
2137 end function getcoordinateresolution
2140 function getcoordinateinterfaces( CS, undo_scaling )
2142 logical,
optional,
intent(in) :: undo_scaling
2144 real,
dimension(CS%nk+1) :: getcoordinateinterfaces
2148 unscale = .false. ;
if (
present(undo_scaling)) unscale = undo_scaling
2152 if (cs%regridding_scheme == regridding_rho)
then
2153 if (.not. cs%target_density_set) &
2154 call mom_error(fatal,
'MOM_regridding, getCoordinateInterfaces: '//&
2155 'target densities not set!')
2158 getcoordinateinterfaces(:) = cs%coord_scale * cs%target_density(:)
2160 getcoordinateinterfaces(:) = cs%target_density(:)
2164 getcoordinateinterfaces(1) = 0.
2166 getcoordinateinterfaces(k+1) = getcoordinateinterfaces(k) - &
2167 cs%coord_scale * cs%coordinateResolution(k)
2170 getcoordinateinterfaces(1) = 0.
2172 getcoordinateinterfaces(k+1) = getcoordinateinterfaces(k) - &
2173 cs%coordinateResolution(k)
2178 getcoordinateinterfaces(:) = abs( getcoordinateinterfaces(:) )
2181 end function getcoordinateinterfaces
2185 function getcoordinateunits( CS )
2187 character(len=20) :: getcoordinateunits
2189 select case ( cs%regridding_scheme )
2190 case ( regridding_zstar, regridding_hycom1, regridding_slight, regridding_adaptive )
2191 getcoordinateunits =
'meter'
2192 case ( regridding_sigma_shelf_zstar )
2193 getcoordinateunits =
'meter/fraction'
2194 case ( regridding_sigma )
2195 getcoordinateunits =
'fraction'
2196 case ( regridding_rho )
2197 getcoordinateunits =
'kg/m3'
2198 case ( regridding_arbitrary )
2199 getcoordinateunits =
'unknown'
2201 call mom_error(fatal,
'MOM_regridding, getCoordinateUnits: '//&
2202 'Unknown regridding scheme selected!')
2205 end function getcoordinateunits
2209 function getcoordinateshortname( CS )
2211 character(len=20) :: getcoordinateshortname
2213 select case ( cs%regridding_scheme )
2214 case ( regridding_zstar )
2217 getcoordinateshortname =
'pseudo-depth, -z*'
2218 case ( regridding_sigma_shelf_zstar )
2219 getcoordinateshortname =
'pseudo-depth, -z*/sigma'
2220 case ( regridding_sigma )
2221 getcoordinateshortname =
'sigma'
2222 case ( regridding_rho )
2223 getcoordinateshortname =
'rho'
2224 case ( regridding_arbitrary )
2225 getcoordinateshortname =
'coordinate'
2226 case ( regridding_hycom1 )
2227 getcoordinateshortname =
'z-rho'
2228 case ( regridding_slight )
2229 getcoordinateshortname =
's-rho'
2230 case ( regridding_adaptive )
2231 getcoordinateshortname =
'adaptive'
2233 call mom_error(fatal,
'MOM_regridding, getCoordinateShortName: '//&
2234 'Unknown regridding scheme selected!')
2237 end function getcoordinateshortname
2240 subroutine set_regrid_params( CS, boundary_extrapolation, min_thickness, old_grid_weight, &
2241 interp_scheme, depth_of_time_filter_shallow, depth_of_time_filter_deep, &
2242 compress_fraction, ref_pressure, dz_min_surface, nz_fixed_surface, Rho_ML_avg_depth, &
2243 nlay_ML_to_interior, fix_haloclines, halocline_filt_len, &
2244 halocline_strat_tol, integrate_downward_for_e, remap_answers_2018, &
2245 adaptTimeRatio, adaptZoom, adaptZoomCoeff, adaptBuoyCoeff, adaptAlpha, adaptDoMin, adaptDrho0)
2247 logical,
optional,
intent(in) :: boundary_extrapolation
2248 real,
optional,
intent(in) :: min_thickness
2250 real,
optional,
intent(in) :: old_grid_weight
2251 character(len=*),
optional,
intent(in) :: interp_scheme
2252 real,
optional,
intent(in) :: depth_of_time_filter_shallow
2253 real,
optional,
intent(in) :: depth_of_time_filter_deep
2254 real,
optional,
intent(in) :: compress_fraction
2255 real,
optional,
intent(in) :: ref_pressure
2257 real,
optional,
intent(in) :: dz_min_surface
2259 integer,
optional,
intent(in) :: nz_fixed_surface
2260 real,
optional,
intent(in) :: rho_ml_avg_depth
2262 real,
optional,
intent(in) :: nlay_ml_to_interior
2264 logical,
optional,
intent(in) :: fix_haloclines
2265 real,
optional,
intent(in) :: halocline_filt_len
2267 real,
optional,
intent(in) :: halocline_strat_tol
2269 logical,
optional,
intent(in) :: integrate_downward_for_e
2271 logical,
optional,
intent(in) :: remap_answers_2018
2274 real,
optional,
intent(in) :: adapttimeratio
2275 real,
optional,
intent(in) :: adaptzoom
2276 real,
optional,
intent(in) :: adaptzoomcoeff
2277 real,
optional,
intent(in) :: adaptbuoycoeff
2278 real,
optional,
intent(in) :: adaptalpha
2279 logical,
optional,
intent(in) :: adaptdomin
2282 real,
optional,
intent(in) :: adaptdrho0
2285 if (
present(interp_scheme))
call set_interp_scheme(cs%interp_CS, interp_scheme)
2286 if (
present(boundary_extrapolation))
call set_interp_extrap(cs%interp_CS, boundary_extrapolation)
2288 if (
present(old_grid_weight))
then
2289 if (old_grid_weight<0. .or. old_grid_weight>1.) &
2290 call mom_error(fatal,
'MOM_regridding, set_regrid_params: Weight is out side the range 0..1!')
2291 cs%old_grid_weight = old_grid_weight
2293 if (
present(depth_of_time_filter_shallow)) cs%depth_of_time_filter_shallow = depth_of_time_filter_shallow
2294 if (
present(depth_of_time_filter_deep)) cs%depth_of_time_filter_deep = depth_of_time_filter_deep
2295 if (
present(depth_of_time_filter_shallow) .or.
present(depth_of_time_filter_deep))
then
2296 if (cs%depth_of_time_filter_deep<cs%depth_of_time_filter_shallow)
call mom_error(fatal,
'MOM_regridding, '//&
2297 'set_regrid_params: depth_of_time_filter_deep<depth_of_time_filter_shallow!')
2300 if (
present(min_thickness)) cs%min_thickness = min_thickness
2301 if (
present(compress_fraction)) cs%compressibility_fraction = compress_fraction
2302 if (
present(ref_pressure)) cs%ref_pressure = ref_pressure
2303 if (
present(integrate_downward_for_e)) cs%integrate_downward_for_e = integrate_downward_for_e
2304 if (
present(remap_answers_2018)) cs%remap_answers_2018 = remap_answers_2018
2306 select case (cs%regridding_scheme)
2307 case (regridding_zstar)
2308 if (
present(min_thickness))
call set_zlike_params(cs%zlike_CS, min_thickness=min_thickness)
2309 case (regridding_sigma_shelf_zstar)
2310 if (
present(min_thickness))
call set_zlike_params(cs%zlike_CS, min_thickness=min_thickness)
2311 case (regridding_sigma)
2312 if (
present(min_thickness))
call set_sigma_params(cs%sigma_CS, min_thickness=min_thickness)
2313 case (regridding_rho)
2314 if (
present(min_thickness))
call set_rho_params(cs%rho_CS, min_thickness=min_thickness)
2315 if (
present(integrate_downward_for_e)) &
2316 call set_rho_params(cs%rho_CS, integrate_downward_for_e=integrate_downward_for_e)
2317 if (
associated(cs%rho_CS) .and. (
present(interp_scheme) .or.
present(boundary_extrapolation))) &
2318 call set_rho_params(cs%rho_CS, interp_cs=cs%interp_CS)
2319 case (regridding_hycom1)
2320 if (
associated(cs%hycom_CS) .and. (
present(interp_scheme) .or.
present(boundary_extrapolation))) &
2321 call set_hycom_params(cs%hycom_CS, interp_cs=cs%interp_CS)
2322 case (regridding_slight)
2323 if (
present(min_thickness))
call set_slight_params(cs%slight_CS, min_thickness=min_thickness)
2324 if (
present(dz_min_surface))
call set_slight_params(cs%slight_CS, dz_ml_min=dz_min_surface)
2325 if (
present(nz_fixed_surface))
call set_slight_params(cs%slight_CS, nz_fixed_surface=nz_fixed_surface)
2326 if (
present(rho_ml_avg_depth))
call set_slight_params(cs%slight_CS, rho_ml_avg_depth=rho_ml_avg_depth)
2327 if (
present(nlay_ml_to_interior))
call set_slight_params(cs%slight_CS, nlay_ml_offset=nlay_ml_to_interior)
2328 if (
present(fix_haloclines))
call set_slight_params(cs%slight_CS, fix_haloclines=fix_haloclines)
2329 if (
present(halocline_filt_len))
call set_slight_params(cs%slight_CS, halocline_filter_length=halocline_filt_len)
2330 if (
present(halocline_strat_tol))
call set_slight_params(cs%slight_CS, halocline_strat_tol=halocline_strat_tol)
2331 if (
present(compress_fraction))
call set_slight_params(cs%slight_CS, compressibility_fraction=compress_fraction)
2332 if (
associated(cs%slight_CS) .and. (
present(interp_scheme) .or.
present(boundary_extrapolation))) &
2333 call set_slight_params(cs%slight_CS, interp_cs=cs%interp_CS)
2334 case (regridding_adaptive)
2335 if (
present(adapttimeratio))
call set_adapt_params(cs%adapt_CS, adapttimeratio=adapttimeratio)
2336 if (
present(adaptzoom))
call set_adapt_params(cs%adapt_CS, adaptzoom=adaptzoom)
2337 if (
present(adaptzoomcoeff))
call set_adapt_params(cs%adapt_CS, adaptzoomcoeff=adaptzoomcoeff)
2338 if (
present(adaptbuoycoeff))
call set_adapt_params(cs%adapt_CS, adaptbuoycoeff=adaptbuoycoeff)
2339 if (
present(adaptalpha))
call set_adapt_params(cs%adapt_CS, adaptalpha=adaptalpha)
2340 if (
present(adaptdomin))
call set_adapt_params(cs%adapt_CS, adaptdomin=adaptdomin)
2341 if (
present(adaptdrho0))
call set_adapt_params(cs%adapt_CS, adaptdrho0=adaptdrho0)
2344 end subroutine set_regrid_params
2347 integer function get_regrid_size(CS)
2350 get_regrid_size = cs%nk
2352 end function get_regrid_size
2355 function get_zlike_cs(CS)
2359 get_zlike_cs = cs%zlike_CS
2360 end function get_zlike_cs
2363 function get_sigma_cs(CS)
2367 get_sigma_cs = cs%sigma_CS
2368 end function get_sigma_cs
2371 function get_rho_cs(CS)
2373 type(
rho_cs) :: get_rho_cs
2375 get_rho_cs = cs%rho_CS
2376 end function get_rho_cs
2380 function getstaticthickness( CS, SSH, depth )
2382 real,
intent(in) :: ssh
2383 real,
intent(in) :: depth
2384 real,
dimension(CS%nk) :: getstaticthickness
2389 select case ( cs%regridding_scheme )
2390 case ( regridding_zstar, regridding_sigma_shelf_zstar, regridding_hycom1, regridding_slight, regridding_adaptive )
2394 dz = cs%coordinateResolution(k) * ( 1. + ssh/depth )
2396 dz = min(dz, depth - z)
2398 getstaticthickness(k) = dz
2401 getstaticthickness(:) = 0.
2403 case ( regridding_sigma )
2404 getstaticthickness(:) = cs%coordinateResolution(:) * ( depth + ssh )
2405 case ( regridding_rho )
2406 getstaticthickness(:) = 0.
2407 case ( regridding_arbitrary )
2408 getstaticthickness(:) = 0.
2410 call mom_error(fatal,
'MOM_regridding, getStaticThickness: '//&
2411 'Unknown regridding scheme selected!')
2414 end function getstaticthickness
2417 subroutine dz_function1( string, dz )
2418 character(len=*),
intent(in) :: string
2420 real,
dimension(:),
intent(inout) :: dz
2423 real :: dz_min, power, prec, H_total
2427 read( string, *) dz_min, h_total, power, prec
2428 if (prec == -1024.)
call mom_error(fatal,
"dz_function1: "// &
2429 "Problem reading FNC1: string ="//trim(string))
2432 dz(k) = (real(k-1)/real(nk-1))**power
2434 dz(:) = ( h_total - real(nk) * dz_min ) * ( dz(:) / sum(dz) )
2435 dz(:) = anint( dz(:) / prec ) * prec
2436 dz(:) = ( h_total - real(nk) * dz_min ) * ( dz(:) / sum(dz) )
2437 dz(:) = anint( dz(:) / prec ) * prec
2438 dz(nk) = dz(nk) + ( h_total - sum( dz(:) + dz_min ) )
2439 dz(:) = anint( dz(:) / prec ) * prec
2440 dz(:) = dz(:) + dz_min
2442 end subroutine dz_function1
2446 integer function rho_function1( string, rho_target )
2447 character(len=*),
intent(in) :: string
2449 real,
dimension(:),
allocatable,
intent(inout) :: rho_target
2451 integer :: nki, k, nk
2452 real :: ddx, dx, rho_1, rho_2, rho_3, drho, rho_4, drho_min
2454 read( string, *) nk, rho_1, rho_2, rho_3, drho, rho_4, drho_min
2455 allocate(rho_target(nk+1))
2457 rho_target(1) = rho_1
2458 rho_target(2) = rho_2
2461 ddx = max( drho_min, real(nki-k)/real(nki*nki) )
2463 rho_target(3+k) = rho_3 + (2. * drho) * dx
2465 rho_target(nki+4) = rho_4
2469 end function rho_function1