MOM6
|
Generates vertical grids as part of the ALE algorithm.
A vertical grid is defined solely by the cell thicknesses, \(h\). Most calculations in this module start with the coordinate at the bottom of the column set to -depth, and use a increasing value of coordinate with decreasing k. This is consistent with the rest of MOM6 that uses position, \(z\) which is a negative quantity for most of the ocean.
A change in grid is define through a change in position of the interfaces:
\[ z^n_{k+1/2} = z^{n-1}_{k+1/2} + \Delta z_{k+1/2} \]
with the positive upward coordinate convention
\[ z_{k-1/2} = z_{k+1/2} + h_k \]
so that
\[ h^n_k = h^{n-1}_k + ( \Delta z_{k-1/2} - \Delta z_{k+1/2} ) \]
Original date of creation: 2008.06.09 by L. White
Data Types | |
type | regridding_cs |
Regridding control structure. More... | |
Functions/Subroutines | |
subroutine, public | initialize_regridding (CS, GV, max_depth, param_file, mod, coord_mode, param_prefix, param_suffix) |
Initialization and configures a regridding control structure based on customizable run-time parameters. More... | |
subroutine | check_grid_def (filename, varname, expected_units, msg, ierr) |
Do some basic checks on the vertical grid definition file, variable. More... | |
subroutine, public | end_regridding (CS) |
Deallocation of regridding memory. More... | |
subroutine, public | regridding_main (remapCS, CS, G, GV, h, tv, h_new, dzInterface, frac_shelf_h, conv_adjust) |
subroutine | calc_h_new_by_dz (G, GV, h, dzInterface, h_new) |
Calculates h_new from h + delta_k dzInterface. More... | |
subroutine, public | check_remapping_grid (G, GV, h, dzInterface, msg) |
Check that the total thickness of two grids match. More... | |
subroutine, public | check_grid_column (nk, depth, h, dzInterface, msg) |
Check that the total thickness of new and old grids are consistent. More... | |
subroutine | filtered_grid_motion (CS, nk, z_old, z_new, dz_g) |
Returns the change in interface position motion after filtering and assuming the top and bottom interfaces do not move. The filtering is a function of depth, and is applied as the integrated average filtering over the trajectory of the interface. By design, this code can not give tangled interfaces provided that z_old and z_new are not already tangled. More... | |
subroutine | build_zstar_grid (CS, G, GV, h, dzInterface, frac_shelf_h) |
Builds a z*-ccordinate grid with partial steps (Adcroft and Campin, 2004). z* is defined as z* = (z-eta)/(H+eta)*H s.t. z*=0 when z=eta and z*=-H when z=-H . More... | |
subroutine | build_sigma_grid (CS, G, GV, h, dzInterface) |
subroutine | build_rho_grid (G, GV, h, tv, dzInterface, remapCS, CS) |
subroutine | build_grid_hycom1 (G, GV, h, tv, dzInterface, CS) |
Builds a simple HyCOM-like grid with the deepest location of potential density interpolated from the column profile and a clipping of depth for each interface to a fixed z* or p* grid. This should probably be (optionally?) changed to find the nearest location of the target density. More... | |
subroutine | build_grid_adaptive (G, GV, h, tv, dzInterface, remapCS, CS) |
subroutine | build_grid_slight (G, GV, h, tv, dzInterface, CS) |
Builds a grid that tracks density interfaces for water that is denser than the surface density plus an increment of some number of layers, and uses all lighter layers uniformly above this location. Note that this amounts to interpolating to find the depth of an arbitrary (non-integer) interface index which should make the results vary smoothly in space to the extent that the surface density and interior stratification vary smoothly in space. Over shallow topography, this will tend to give a uniform sigma-like coordinate. For sufficiently shallow water, a minimum grid spacing is used to avoid certain instabilities. More... | |
subroutine, public | adjust_interface_motion (nk, min_thickness, h_old, dz_int) |
Adjust dz_Interface to ensure non-negative future thicknesses. More... | |
subroutine | build_grid_arbitrary (G, GV, h, dzInterface, h_new, CS) |
subroutine, public | inflate_vanished_layers_old (CS, G, GV, h) |
subroutine | convective_adjustment (G, GV, h, tv) |
real function, dimension(nk), public | uniformresolution (nk, coordMode, maxDepth, rhoLight, rhoHeavy) |
subroutine | initcoord (CS, coord_mode) |
subroutine, public | setcoordinateresolution (dz, CS) |
subroutine, public | set_target_densities_from_gv (GV, CS) |
Set target densities based on the old Rlay variable. More... | |
subroutine, public | set_target_densities (CS, rho_int) |
Set target densities based on vector of interface values. More... | |
subroutine, public | set_regrid_max_depths (CS, max_depths, units_to_H) |
Set maximum interface depths based on a vector of input values. More... | |
subroutine, public | set_regrid_max_thickness (CS, max_h, units_to_H) |
Set maximum layer thicknesses based on a vector of input values. More... | |
real function, dimension(cs%nk), public | getcoordinateresolution (CS) |
real function, dimension(cs%nk+1), public | getcoordinateinterfaces (CS) |
Query the target coordinate interface positions. More... | |
character(len=20) function, public | getcoordinateunits (CS) |
character(len=20) function, public | getcoordinateshortname (CS) |
subroutine, public | set_regrid_params (CS, boundary_extrapolation, min_thickness, old_grid_weight, interp_scheme, depth_of_time_filter_shallow, depth_of_time_filter_deep, compress_fraction, dz_min_surface, nz_fixed_surface, Rho_ML_avg_depth, nlay_ML_to_interior, fix_haloclines, halocline_filt_len, halocline_strat_tol, integrate_downward_for_e, adaptTimeRatio, adaptZoom, adaptZoomCoeff, adaptBuoyCoeff, adaptAlpha, adaptDoMin) |
Can be used to set any of the parameters for MOM_regridding. More... | |
integer function, public | get_regrid_size (CS) |
Returns the number of levels/layers in the regridding control structure. More... | |
type(zlike_cs) function, public | get_zlike_cs (CS) |
type(sigma_cs) function, public | get_sigma_cs (CS) |
type(rho_cs) function, public | get_rho_cs (CS) |
real function, dimension(cs%nk), public | getstaticthickness (CS, SSH, depth) |
subroutine | dz_function1 (string, dz) |
Parses a string and generates a dz(:) profile that goes like k**power. More... | |
Variables | |
character(len= *), parameter, public | regriddingcoordinatemodedoc = " LAYER - Isopycnal or stacked shallow water layers\n"// " ZSTAR, Z* - stetched geopotential z*\n"// " SIGMA_SHELF_ZSTAR - stetched geopotential z* ignoring shelf\n"// " SIGMA - terrain following coordinates\n"// " RHO - continuous isopycnal\n"// " HYCOM1 - HyCOM-like hybrid coordinate\n"// " SLIGHT - stretched coordinates above continuous isopycnal\n"// " ADAPTIVE - optimize for smooth neutral density surfaces" |
Documentation for coordinate options. More... | |
character(len= *), parameter, public | regriddinginterpschemedoc = " P1M_H2 (2nd-order accurate)\n"// " P1M_H4 (2nd-order accurate)\n"// " P1M_IH4 (2nd-order accurate)\n"// " PLM (2nd-order accurate)\n"// " PPM_H4 (3rd-order accurate)\n"// " PPM_IH4 (3rd-order accurate)\n"// " P3M_IH4IH3 (4th-order accurate)\n"// " P3M_IH6IH5 (4th-order accurate)\n"// " PQM_IH4IH3 (4th-order accurate)\n"// " PQM_IH6IH5 (5th-order accurate)" |
character(len= *), parameter, public | regriddingdefaultinterpscheme = "P1M_H2" |
logical, parameter, public | regriddingdefaultboundaryextrapolation = .false. |
real, parameter, public | regriddingdefaultminthickness = 1.e-3 |
subroutine, public mom_regridding::adjust_interface_motion | ( | integer, intent(in) | nk, |
real, intent(in) | min_thickness, | ||
real, dimension(nk), intent(in) | h_old, | ||
real, dimension(nk+1), intent(inout) | dz_int | ||
) |
Adjust dz_Interface to ensure non-negative future thicknesses.
[in] | nk | Number of layers |
[in] | min_thickness | Minium allowed thickness of h (H units) |
[in] | h_old | Minium allowed thickness of h (H units) |
[in,out] | dz_int | Minium allowed thickness of h (H units) |
Definition at line 1550 of file MOM_regridding.F90.
Referenced by build_grid_adaptive(), build_grid_hycom1(), build_grid_slight(), and build_zstar_grid().
|
private |
[in] | g | The ocean's grid structure |
[in] | gv | The ocean's vertical grid structure |
[in] | h | Layer thicknesses, in H (usually m or kg m-2) |
[in] | tv | A structure pointing to various thermodynamic variables |
Definition at line 1432 of file MOM_regridding.F90.
References adjust_interface_motion(), and filtered_grid_motion().
Referenced by regridding_main().
|
private |
[in] | g | Ocean grid structure |
[in] | gv | Ocean vertical grid structure |
[in] | h | Original ayer thicknesses, in H |
[in,out] | dzinterface | The change in interface depth in H |
[in,out] | h_new | New layer thicknesses, in H |
[in] | cs | Regridding control structure |
Definition at line 1596 of file MOM_regridding.F90.
Referenced by regridding_main().
|
private |
Builds a simple HyCOM-like grid with the deepest location of potential density interpolated from the column profile and a clipping of depth for each interface to a fixed z* or p* grid. This should probably be (optionally?) changed to find the nearest location of the target density.
[in] | g | Grid structure |
[in] | gv | Ocean vertical grid structure |
[in] | h | Existing model thickness, in H units |
[in] | tv | Thermodynamics structure |
[in,out] | dzinterface | Changes in interface position |
[in] | cs | Regridding control structure |
Definition at line 1382 of file MOM_regridding.F90.
References adjust_interface_motion(), coord_hycom::build_hycom1_column(), filtered_grid_motion(), and mom_error_handler::mom_error().
Referenced by regridding_main().
|
private |
Builds a grid that tracks density interfaces for water that is denser than the surface density plus an increment of some number of layers, and uses all lighter layers uniformly above this location. Note that this amounts to interpolating to find the depth of an arbitrary (non-integer) interface index which should make the results vary smoothly in space to the extent that the surface density and interior stratification vary smoothly in space. Over shallow topography, this will tend to give a uniform sigma-like coordinate. For sufficiently shallow water, a minimum grid spacing is used to avoid certain instabilities.
[in] | g | Grid structure |
[in] | gv | Ocean vertical grid structure |
[in] | h | Existing model thickness, in H units |
[in] | tv | Thermodynamics structure |
[in,out] | dzinterface | Changes in interface position |
[in] | cs | Regridding control structure |
Definition at line 1496 of file MOM_regridding.F90.
References adjust_interface_motion(), and filtered_grid_motion().
Referenced by regridding_main().
|
private |
[in] | g | Ocean grid structure |
[in] | gv | Ocean vertical grid structure |
[in] | h | Layer thicknesses, in H |
[in] | tv | Thermodynamics structure |
[in,out] | dzinterface | The change in interface depth in H |
[in] | remapcs | The remapping control structure |
[in] | cs | Regridding control structure |
Definition at line 1259 of file MOM_regridding.F90.
References filtered_grid_motion(), and mom_error_handler::mom_error().
Referenced by regridding_main().
|
private |
[in] | cs | Regridding control structure |
[in] | g | Ocean grid structure |
[in] | gv | ocean vertical grid structure |
[in] | h | Layer thicknesses, in H |
[in,out] | dzinterface | The change in interface depth in H. |
Definition at line 1182 of file MOM_regridding.F90.
References coord_sigma::build_sigma_column(), filtered_grid_motion(), and mom_error_handler::mom_error().
Referenced by regridding_main().
|
private |
Builds a z*-ccordinate grid with partial steps (Adcroft and Campin, 2004). z* is defined as z* = (z-eta)/(H+eta)*H s.t. z*=0 when z=eta and z*=-H when z=-H .
[in] | cs | Regridding control structure |
[in] | g | Ocean grid structure |
[in] | gv | ocean vertical grid structure |
[in] | h | Layer thicknesses, in H |
[in,out] | dzinterface | The change in interface depth in H. |
frac_shelf_h | Fractional ice shelf coverage. |
Definition at line 1090 of file MOM_regridding.F90.
References adjust_interface_motion(), coord_zlike::build_zstar_column(), filtered_grid_motion(), and mom_error_handler::mom_error().
Referenced by regridding_main().
|
private |
Calculates h_new from h + delta_k dzInterface.
[in] | g | Grid structure |
[in] | gv | Ocean vertical grid structure |
[in] | h | Old layer thicknesses (m) |
[in] | dzinterface | Change in interface positions (m) |
[in,out] | h_new | New layer thicknesses (m) |
Definition at line 855 of file MOM_regridding.F90.
Referenced by regridding_main().
subroutine, public mom_regridding::check_grid_column | ( | integer, intent(in) | nk, |
real, intent(in) | depth, | ||
real, dimension(nk), intent(in) | h, | ||
real, dimension(nk+1), intent(in) | dzInterface, | ||
character(len=*), intent(in) | msg | ||
) |
Check that the total thickness of new and old grids are consistent.
[in] | nk | Number of cells |
[in] | depth | Depth of bottom (m) |
[in] | h | Cell thicknesses (m) |
[in] | dzinterface | Change in interface positions (m) |
[in] | msg | Message to append to errors |
Definition at line 899 of file MOM_regridding.F90.
References mom_error_handler::mom_error().
Referenced by check_remapping_grid().
|
private |
Do some basic checks on the vertical grid definition file, variable.
[in] | filename | File name |
[in] | varname | Variable name |
[in] | expected_units | Expected units of variable |
[in,out] | msg | Message to use for errors |
[out] | ierr | True if an error occurs |
Definition at line 685 of file MOM_regridding.F90.
Referenced by initialize_regridding().
subroutine, public mom_regridding::check_remapping_grid | ( | type(ocean_grid_type), intent(in) | G, |
type(verticalgrid_type), intent(in) | GV, | ||
real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(in) | h, | ||
real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke+1), intent(in) | dzInterface, | ||
character(len=*), intent(in) | msg | ||
) |
Check that the total thickness of two grids match.
[in] | g | Grid structure |
[in] | gv | Ocean vertical grid structure |
[in] | h | Layer thicknesses (m) |
[in] | dzinterface | Change in interface positions (m) |
[in] | msg | Message to append to errors |
Definition at line 880 of file MOM_regridding.F90.
References check_grid_column().
Referenced by regridding_main().
|
private |
[in] | g | The ocean's grid structure |
[in] | gv | The ocean's vertical grid structure |
[in,out] | h | Layer thicknesses, in H (usually m or kg m-2) |
[in,out] | tv | A structure pointing to various thermodynamic variables |
Definition at line 1742 of file MOM_regridding.F90.
Referenced by regridding_main().
|
private |
Parses a string and generates a dz(:) profile that goes like k**power.
[in] | string | String with list of parameters in form dz_min, H_total, power, precision |
[in,out] | dz | Profile of nominal thicknesses |
Definition at line 2219 of file MOM_regridding.F90.
Referenced by initialize_regridding().
subroutine, public mom_regridding::end_regridding | ( | type(regridding_cs), intent(inout) | CS | ) |
Deallocation of regridding memory.
[in,out] | cs | Regridding control structure |
Definition at line 741 of file MOM_regridding.F90.
References coord_adapt::end_coord_adapt(), coord_hycom::end_coord_hycom(), coord_sigma::end_coord_sigma(), coord_slight::end_coord_slight(), and coord_zlike::end_coord_zlike().
|
private |
Returns the change in interface position motion after filtering and assuming the top and bottom interfaces do not move. The filtering is a function of depth, and is applied as the integrated average filtering over the trajectory of the interface. By design, this code can not give tangled interfaces provided that z_old and z_new are not already tangled.
[in] | cs | Regridding control structure |
[in] | nk | Number of cells |
[in] | z_old | Old grid position (m) |
[in] | z_new | New grid position (m) |
[in,out] | dz_g | Change in interface positions (m) |
Definition at line 960 of file MOM_regridding.F90.
References mom_error_handler::mom_error().
Referenced by build_grid_adaptive(), build_grid_hycom1(), build_grid_slight(), build_rho_grid(), build_sigma_grid(), and build_zstar_grid().
integer function, public mom_regridding::get_regrid_size | ( | type(regridding_cs), intent(inout) | CS | ) |
Returns the number of levels/layers in the regridding control structure.
[in,out] | cs | Regridding control structure |
Definition at line 2151 of file MOM_regridding.F90.
type(rho_cs) function, public mom_regridding::get_rho_cs | ( | type(regridding_cs), intent(in) | CS | ) |
Definition at line 2172 of file MOM_regridding.F90.
Referenced by mom_diag_remap::diag_remap_update().
type(sigma_cs) function, public mom_regridding::get_sigma_cs | ( | type(regridding_cs), intent(in) | CS | ) |
Definition at line 2165 of file MOM_regridding.F90.
Referenced by mom_diag_remap::diag_remap_update().
type(zlike_cs) function, public mom_regridding::get_zlike_cs | ( | type(regridding_cs), intent(in) | CS | ) |
Definition at line 2158 of file MOM_regridding.F90.
Referenced by mom_diag_remap::diag_remap_update().
real function, dimension(cs%nk+1), public mom_regridding::getcoordinateinterfaces | ( | type(regridding_cs), intent(in) | CS | ) |
Query the target coordinate interface positions.
[in] | cs | Regridding control structure |
Definition at line 1982 of file MOM_regridding.F90.
real function, dimension(cs%nk), public mom_regridding::getcoordinateresolution | ( | type(regridding_cs), intent(in) | CS | ) |
Definition at line 1973 of file MOM_regridding.F90.
Referenced by mom_state_initialization::mom_temp_salt_initialize_from_z().
character(len=20) function, public mom_regridding::getcoordinateshortname | ( | type(regridding_cs), intent(in) | CS | ) |
Definition at line 2037 of file MOM_regridding.F90.
character(len=20) function, public mom_regridding::getcoordinateunits | ( | type(regridding_cs), intent(in) | CS | ) |
Definition at line 2012 of file MOM_regridding.F90.
real function, dimension(cs%nk), public mom_regridding::getstaticthickness | ( | type(regridding_cs), intent(in) | CS, |
real, intent(in) | SSH, | ||
real, intent(in) | depth | ||
) |
Definition at line 2182 of file MOM_regridding.F90.
Referenced by mom_ale::ale_initthicknesstocoord().
subroutine, public mom_regridding::inflate_vanished_layers_old | ( | type(regridding_cs), intent(in) | CS, |
type(ocean_grid_type), intent(in) | G, | ||
type(verticalgrid_type), intent(in) | GV, | ||
real, dimension(szi_(g),szj_(g), szk_(gv)), intent(inout) | h | ||
) |
[in] | g | The ocean's grid structure |
[in] | gv | The ocean's vertical grid structure |
[in,out] | h | Layer thicknesses, in H (usually m or kg m-2) |
Definition at line 1700 of file MOM_regridding.F90.
|
private |
Definition at line 1842 of file MOM_regridding.F90.
Referenced by initialize_regridding().
subroutine, public mom_regridding::initialize_regridding | ( | type(regridding_cs), intent(inout) | CS, |
type(verticalgrid_type), intent(in) | GV, | ||
real, intent(in) | max_depth, | ||
type(param_file_type), intent(in) | param_file, | ||
character(len=*), intent(in) | mod, | ||
character(len=*), intent(in) | coord_mode, | ||
character(len=*), intent(in) | param_prefix, | ||
character(len=*), intent(in) | param_suffix | ||
) |
Initialization and configures a regridding control structure based on customizable run-time parameters.
[in,out] | cs | Regridding control structure |
[in] | gv | Ocean vertical grid structure |
[in] | max_depth | The maximum depth of the ocean, in m. |
[in] | param_file | Parameter file |
[in] | mod | Name of calling module. |
[in] | coord_mode | Coordinate mode |
[in] | param_prefix | String to prefix to parameter names. If empty, causes main model parameters to be used. |
[in] | param_suffix | String to append to parameter names. |
Definition at line 164 of file MOM_regridding.F90.
References check_grid_def(), dz_function1(), mom_string_functions::extract_integer(), mom_string_functions::extract_real(), mom_string_functions::extractword(), initcoord(), mom_error_handler::mom_error(), regrid_consts::regridding_adaptive, regrid_consts::regridding_hycom1, regrid_consts::regridding_slight, regriddingdefaultboundaryextrapolation, regriddingdefaultinterpscheme, regriddingdefaultminthickness, regriddinginterpschemedoc, set_regrid_max_depths(), set_regrid_max_thickness(), set_regrid_params(), set_target_densities(), set_target_densities_from_gv(), setcoordinateresolution(), and uniformresolution().
subroutine, public mom_regridding::regridding_main | ( | type(remapping_cs), intent(in) | remapCS, |
type(regridding_cs), intent(in) | CS, | ||
type(ocean_grid_type), intent(in) | G, | ||
type(verticalgrid_type), intent(in) | GV, | ||
real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(inout) | h, | ||
type(thermo_var_ptrs), intent(inout) | tv, | ||
real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(inout) | h_new, | ||
real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke+1), intent(inout) | dzInterface, | ||
real, dimension(:,:), optional, pointer | frac_shelf_h, | ||
logical, intent(in), optional | conv_adjust | ||
) |
[in] | remapcs | Remapping parameters and options |
[in] | cs | Regridding control structure |
[in] | g | Ocean grid structure |
[in] | gv | Ocean vertical grid structure |
[in,out] | h | Current 3D grid obtained after the last time step |
[in,out] | tv | Thermodynamical variables (T, S, ...) |
[in,out] | h_new | New 3D grid consistent with target coordinate |
[in,out] | dzinterface | The change in position of each interface |
frac_shelf_h | Fractional ice shelf coverage |
Definition at line 761 of file MOM_regridding.F90.
References build_grid_adaptive(), build_grid_arbitrary(), build_grid_hycom1(), build_grid_slight(), build_rho_grid(), build_sigma_grid(), build_zstar_grid(), calc_h_new_by_dz(), check_remapping_grid(), convective_adjustment(), mom_error_handler::mom_error(), regrid_consts::regridding_adaptive, regrid_consts::regridding_hycom1, and regrid_consts::regridding_slight.
subroutine, public mom_regridding::set_regrid_max_depths | ( | type(regridding_cs), intent(inout) | CS, |
real, dimension(cs%nk+1), intent(in) | max_depths, | ||
real, intent(in), optional | units_to_H | ||
) |
Set maximum interface depths based on a vector of input values.
[in,out] | cs | Regridding control structure |
[in] | max_depths | Maximum interface depths, in arbitrary units |
[in] | units_to_h | A conversion factor for max_depths into H units |
Definition at line 1906 of file MOM_regridding.F90.
Referenced by initialize_regridding().
subroutine, public mom_regridding::set_regrid_max_thickness | ( | type(regridding_cs), intent(inout) | CS, |
real, dimension(cs%nk+1), intent(in) | max_h, | ||
real, intent(in), optional | units_to_H | ||
) |
Set maximum layer thicknesses based on a vector of input values.
[in,out] | cs | Regridding control structure |
[in] | max_h | Maximum interface depths, in arbitrary units |
[in] | units_to_h | A conversion factor for max_h into H units |
Definition at line 1944 of file MOM_regridding.F90.
Referenced by initialize_regridding().
subroutine, public mom_regridding::set_regrid_params | ( | type(regridding_cs), intent(inout) | CS, |
logical, intent(in), optional | boundary_extrapolation, | ||
real, intent(in), optional | min_thickness, | ||
real, intent(in), optional | old_grid_weight, | ||
character(len=*), intent(in), optional | interp_scheme, | ||
real, intent(in), optional | depth_of_time_filter_shallow, | ||
real, intent(in), optional | depth_of_time_filter_deep, | ||
real, intent(in), optional | compress_fraction, | ||
real, intent(in), optional | dz_min_surface, | ||
integer, intent(in), optional | nz_fixed_surface, | ||
real, intent(in), optional | Rho_ML_avg_depth, | ||
real, intent(in), optional | nlay_ML_to_interior, | ||
logical, intent(in), optional | fix_haloclines, | ||
real, intent(in), optional | halocline_filt_len, | ||
real, intent(in), optional | halocline_strat_tol, | ||
logical, intent(in), optional | integrate_downward_for_e, | ||
real, intent(in), optional | adaptTimeRatio, | ||
real, intent(in), optional | adaptZoom, | ||
real, intent(in), optional | adaptZoomCoeff, | ||
real, intent(in), optional | adaptBuoyCoeff, | ||
real, intent(in), optional | adaptAlpha, | ||
logical, intent(in), optional | adaptDoMin | ||
) |
Can be used to set any of the parameters for MOM_regridding.
[in,out] | cs | Regridding control structure |
[in] | boundary_extrapolation | Extrapolate in boundary cells |
[in] | min_thickness | Minimum thickness allowed when building the new grid (m) |
[in] | old_grid_weight | Weight given to old coordinate when time-filtering grid |
[in] | interp_scheme | Interpolation method for state-dependent coordinates |
[in] | depth_of_time_filter_shallow | Depth to start cubic (H units) |
[in] | depth_of_time_filter_deep | Depth to end cubic (H units) |
[in] | compress_fraction | Fraction of compressibility to add to potential density |
[in] | dz_min_surface | The fixed resolution in the topmost SLight_nkml_min layers (m) |
[in] | nz_fixed_surface | The number of fixed-thickess layers at the top of the model |
[in] | rho_ml_avg_depth | Averaging depth over which to determine mixed layer potential density (m) |
[in] | nlay_ml_to_interior | Number of layers to offset the mixed layer density to find resolved stratification (nondim) |
[in] | fix_haloclines | Detect regions with much weaker stratification in the coordinate |
[in] | halocline_filt_len | Length scale over which to filter T & S when looking for spuriously unstable water mass profiles (m) |
[in] | halocline_strat_tol | Value of the stratification ratio that defines a problematic halocline region. |
[in] | integrate_downward_for_e | If true, integrate for interface positions downward from the top. |
Definition at line 2073 of file MOM_regridding.F90.
Referenced by initialize_regridding(), and mom_state_initialization::mom_temp_salt_initialize_from_z().
subroutine, public mom_regridding::set_target_densities | ( | type(regridding_cs), intent(inout) | CS, |
real, dimension(cs%nk+1), intent(in) | rho_int | ||
) |
Set target densities based on vector of interface values.
[in,out] | cs | Regridding control structure |
[in] | rho_int | Interface densities |
Definition at line 1896 of file MOM_regridding.F90.
Referenced by initialize_regridding().
subroutine, public mom_regridding::set_target_densities_from_gv | ( | type(verticalgrid_type), intent(in) | GV, |
type(regridding_cs), intent(inout) | CS | ||
) |
Set target densities based on the old Rlay variable.
[in] | gv | Ocean vertical grid structure |
[in,out] | cs | Regridding control structure |
Definition at line 1879 of file MOM_regridding.F90.
Referenced by initialize_regridding().
subroutine, public mom_regridding::setcoordinateresolution | ( | real, dimension(:), intent(in) | dz, |
type(regridding_cs), intent(inout) | CS | ||
) |
Definition at line 1867 of file MOM_regridding.F90.
Referenced by initialize_regridding().
real function, dimension(nk), public mom_regridding::uniformresolution | ( | integer, intent(in) | nk, |
character(len=*), intent(in) | coordMode, | ||
real, intent(in) | maxDepth, | ||
real, intent(in) | rhoLight, | ||
real, intent(in) | rhoHeavy | ||
) |
Definition at line 1808 of file MOM_regridding.F90.
Referenced by initialize_regridding().
character(len=*), parameter, public mom_regridding::regriddingcoordinatemodedoc = " LAYER - Isopycnal or stacked shallow water layers\n"// " ZSTAR, Z* - stetched geopotential z*\n"// " SIGMA_SHELF_ZSTAR - stetched geopotential z* ignoring shelf\n"// " SIGMA - terrain following coordinates\n"// " RHO - continuous isopycnal\n"// " HYCOM1 - HyCOM-like hybrid coordinate\n"// " SLIGHT - stretched coordinates above continuous isopycnal\n"// " ADAPTIVE - optimize for smooth neutral density surfaces" |
logical, parameter, public mom_regridding::regriddingdefaultboundaryextrapolation = .false. |
Definition at line 155 of file MOM_regridding.F90.
Referenced by initialize_regridding().
character(len=*), parameter, public mom_regridding::regriddingdefaultinterpscheme = "P1M_H2" |
Definition at line 154 of file MOM_regridding.F90.
Referenced by initialize_regridding().
real, parameter, public mom_regridding::regriddingdefaultminthickness = 1.e-3 |
Definition at line 156 of file MOM_regridding.F90.
Referenced by initialize_regridding().
character(len=*), parameter, public mom_regridding::regriddinginterpschemedoc = " P1M_H2 (2nd-order accurate)\n"// " P1M_H4 (2nd-order accurate)\n"// " P1M_IH4 (2nd-order accurate)\n"// " PLM (2nd-order accurate)\n"// " PPM_H4 (3rd-order accurate)\n"// " PPM_IH4 (3rd-order accurate)\n"// " P3M_IH4IH3 (4th-order accurate)\n"// " P3M_IH6IH5 (4th-order accurate)\n"// " PQM_IH4IH3 (4th-order accurate)\n"// " PQM_IH6IH5 (5th-order accurate)" |