mom_hybgen_regrid module reference

This module contains the hybgen regridding routines from HYCOM, with minor modifications to follow the MOM6 coding conventions.

More…

Data Types

hybgen_regrid_cs

Control structure containing required parameters for the hybgen coordinate generator.

Functions/Subroutines

init_hybgen_regrid()

Initialise a hybgen_regrid_CS control structure and store its parameters.

write_hybgen_coord_file()

Writes out a file containing any available data related to the vertical grid used by the MOM ocean model.

end_hybgen_regrid()

This subroutine deallocates memory in the control structure for the hybgen module.

get_hybgen_regrid_params()

This subroutine can be used to retrieve the parameters for the hybgen regrid module.

hybgen_regrid()

Modify the input grid to give a new vertical grid based on the HYCOM hybgen code.

hybgen_column_init()

Initialize some of the variables that are used for regridding or unmixing, including the stretched contraits on where the new interfaces can be.

cushn()

The cushion function from Bleck & Benjamin, 1992, which returns a smoothly varying but limited value that goes between dp0 and delp.

hybgen_column_regrid()

Create a new grid for a column of water using the Hybgen algorithm.

Detailed Description

This module contains the hybgen regridding routines from HYCOM, with minor modifications to follow the MOM6 coding conventions.

Type Documentation

type mom_hybgen_regrid/hybgen_regrid_cs

Control structure containing required parameters for the hybgen coordinate generator.

Type fields:
  • % min_thickness [real] :: Minimum thickness allowed for layers [H ~> m or kg m-2].

  • % nk [integer] :: Number of layers on the target grid.

  • % ref_pressure [real] :: Reference pressure for density calculations [R L2 T-2 ~> Pa].

  • % hybiso [real] :: Hybgen uses PCM if layer is within hybiso of target density [R ~> kg m-3].

  • % nsigma [integer] :: Number of sigma levels used by HYBGEN.

  • % dp00i [real] :: Deep isopycnal spacing minimum thickness [H ~> m or kg m-2].

  • % qhybrlx [real] :: Fractional relaxation within a regridding step [nondim].

  • % dp0k [real(:),allocatable] :: minimum deep z-layer separation [H ~> m or kg m-2]

  • % ds0k [real(:),allocatable] :: minimum shallow z-layer separation [H ~> m or kg m-2]

  • % coord_scale [real] :: A scaling factor to restores the depth coordinates to values in m [m H-1 ~> 1 or m3 kg-1].

  • % rho_coord_scale [real] :: A scaling factor to restores the denesity coordinates to values in kg m-3 [kg m-3 R-1 ~> 1].

  • % dpns [real] :: depth to start terrain following [H ~> m or kg m-2]

  • % dsns [real] :: depth to stop terrain following [H ~> m or kg m-2]

  • % min_dilate [real] :: The minimum amount of dilation that is permitted when converting target coordinates from z to z* [nondim]. This limit applies when wetting occurs.

  • % max_dilate [real] :: The maximum amount of dilation that is permitted when converting target coordinates from z to z* [nondim]. This limit applies when drying occurs.

  • % thkbot [real] :: Thickness of a bottom boundary layer, within which hybgen does something different. [H ~> m or kg m-2].

  • % topiso_const [real] :: Shallowest depth for isopycnal layers [H ~> m or kg m-2].

  • % target_density [real(:),allocatable] :: Nominal density of interfaces [R ~> kg m-3].

  • % dp_far_from_sfc [real] :: A distance that determines when an interface is suffiently far from the surface that certain adjustments can be made in the Hybgen regridding code [H ~> m or kg m-2]. In Hycom, this is set to tenm (nominally 10 m).

  • % dp_far_from_bot [real] :: A distance that determines when an interface is suffiently far from the bottom that certain adjustments can be made in the Hybgen regridding code [H ~> m or kg m-2]. In Hycom, this is set to onem (nominally 1 m).

  • % h_thin [real] :: A layer thickness below which a layer is considered to be too thin for certain adjustments to be made in the Hybgen regridding code [H ~> m or kg m-2]. In Hycom, this is set to onemm (nominally 0.001 m).

  • % rho_eps [real] :: A small nonzero density that is used to prevent division by zero in several expressions in the Hybgen regridding code [R ~> kg m-3].

  • % onem [real] :: Nominally one m in thickness units [H ~> m or kg m-2], used only in certain debugging tests.

Function/Subroutine Documentation

subroutine mom_hybgen_regrid/init_hybgen_regrid(CS, GV, US, param_file)

Initialise a hybgen_regrid_CS control structure and store its parameters.

Parameters:
  • cs :: Unassociated pointer to hold the control structure

  • gv :: [in] Ocean vertical grid structure

  • us :: [in] A dimensional unit scaling type

  • param_file :: [in] Parameter file

Call to:

mom_error_handler::mom_error mom_string_functions::slasher

subroutine mom_hybgen_regrid/write_hybgen_coord_file(GV, CS, filepath)

Writes out a file containing any available data related to the vertical grid used by the MOM ocean model.

Parameters:
  • gv :: [in] The ocean’s vertical grid structure

  • cs :: [in] Control structure for this module

  • filepath :: [in] The full path to the file to write

Call to:

mom_io::var_desc

Called from:

mom_regridding::write_regrid_file

subroutine mom_hybgen_regrid/end_hybgen_regrid(CS)

This subroutine deallocates memory in the control structure for the hybgen module.

Parameters:

cs :: Coordinate control structure

subroutine mom_hybgen_regrid/get_hybgen_regrid_params(CS, nk, ref_pressure, hybiso, nsigma, dp00i, qhybrlx, dp0k, ds0k, dpns, dsns, min_dilate, max_dilate, thkbot, topiso_const, target_density)

This subroutine can be used to retrieve the parameters for the hybgen regrid module.

Parameters:
  • cs :: Coordinate regridding control structure

  • nk :: [out] Number of layers on the target grid

  • ref_pressure :: [out] Reference pressure for density calculations [R L2 T-2 ~> Pa]

  • hybiso :: [out] Hybgen uses PCM if layer is within hybiso of target density [R ~> kg m-3]

  • nsigma :: [out] Number of sigma levels used by HYBGEN

  • dp00i :: [out] Deep isopycnal spacing minimum thickness [H ~> m or kg m-2]

  • qhybrlx :: [out] Fractional relaxation amount per timestep, 0 < qyhbrlx <= 1 [nondim]

  • dp0k :: [out] minimum deep z-layer separation [H ~> m or kg m-2]

  • ds0k :: [out] minimum shallow z-layer separation [H ~> m or kg m-2]

  • dpns :: [out] depth to start terrain following [H ~> m or kg m-2]

  • dsns :: [out] depth to stop terrain following [H ~> m or kg m-2]

  • min_dilate :: [out] The minimum amount of dilation that is permitted when converting target coordinates from z to z* [nondim]. This limit applies when wetting occurs.

  • max_dilate :: [out] The maximum amount of dilation that is permitted when converting target coordinates from z to z* [nondim]. This limit applies when drying occurs.

  • thkbot :: [out] Thickness of a bottom boundary layer, within which hybgen does something different. [H ~> m or kg m-2]

  • topiso_const :: [out] Shallowest depth for isopycnal layers [H ~> m or kg m-2]

  • target_density :: [out] Nominal density of interfaces [R ~> kg m-3]

Call to:

mom_error_handler::mom_error

Called from:

mom_hybgen_unmix::init_hybgen_unmix

subroutine mom_hybgen_regrid/hybgen_regrid(G, GV, US, dp, nom_depth_H, tv, CS, dzInterface, PCM_cell)

Modify the input grid to give a new vertical grid based on the HYCOM hybgen code.

Parameters:
  • g :: [in] Ocean grid structure

  • gv :: [in] Ocean vertical grid structure

  • us :: [in] A dimensional unit scaling type

  • dp :: [in] Source grid layer thicknesses [H ~> m or kg m-2]

  • nom_depth_h :: [in] The bathymetric depth of this column

  • tv :: [in] Thermodynamics structure

  • cs :: [in] hybgen control structure

  • dzinterface :: [inout] The change in height of each interface,

  • pcm_cell :: [inout] If true, PCM remapping should be used in a cell.

Call to:

hybgen_column_init hybgen_column_regrid

subroutine mom_hybgen_regrid/hybgen_column_init(nk, nsigma, dp0k, ds0k, dp00i, topiso_i_j, qhybrlx, dpns, dsns, h_tot, dilate, h_col, fixlay, qhrlx, dp0ij, dp0cum)

Initialize some of the variables that are used for regridding or unmixing, including the stretched contraits on where the new interfaces can be.

Parameters:
  • nk :: [in] The number of layers in the new grid

  • nsigma :: [in] The number of sigma levels

  • dp0k :: [in] Layer deep z-level spacing minimum thicknesses [H ~> m or kg m-2]

  • ds0k :: [in] Layer shallow z-level spacing minimum thicknesses [H ~> m or kg m-2]

  • dp00i :: [in] Deep isopycnal spacing minimum thickness [H ~> m or kg m-2]

  • topiso_i_j :: [in] Shallowest depth for isopycnal layers [H ~> m or kg m-2]

  • qhybrlx :: [in] Fractional relaxation amount per timestep, 0 < qyhbrlx <= 1 [nondim]

  • h_tot :: [in] The sum of the initial layer thicknesses [H ~> m or kg m-2]

  • dilate :: [in] A factor by which to dilate the target positions from z to z* [nondim]

  • h_col :: [in] Initial layer thicknesses [H ~> m or kg m-2]

  • dpns :: [in] Vertical sum of dp0k [H ~> m or kg m-2]

  • dsns :: [in] Vertical sum of ds0k [H ~> m or kg m-2]

  • fixlay :: [out] Deepest fixed coordinate layer

  • qhrlx :: [out] Fractional relaxation within a timestep (between 0 and 1) [nondim]

  • dp0ij :: [out] minimum layer thickness [H ~> m or kg m-2]

  • dp0cum :: [out] minimum interface depth [H ~> m or kg m-2]

Called from:

hybgen_regrid

function mom_hybgen_regrid/cushn(delp, dp0) [real]

The cushion function from Bleck & Benjamin, 1992, which returns a smoothly varying but limited value that goes between dp0 and delp.

Called from:

hybgen_column_regrid

subroutine mom_hybgen_regrid/hybgen_column_regrid(CS, nk, thkbot, Rcv_tgt, fixlay, qhrlx, dp0ij, dp0cum, Rcv, h_in, dp_int)

Create a new grid for a column of water using the Hybgen algorithm.

Parameters:
  • cs :: [in] hybgen regridding control structure

  • nk :: [in] number of layers

  • thkbot :: [in] thickness of bottom boundary layer [H ~> m or kg m-2]

  • rcv_tgt :: [in] Target potential density [R ~> kg m-3]

  • fixlay :: [in] deepest fixed coordinate layer

  • qhrlx :: [in] relaxation coefficient per timestep [nondim]

  • dp0ij :: [in] minimum layer thickness [H ~> m or kg m-2]

  • dp0cum :: [in] minimum interface depth [H ~> m or kg m-2]

  • rcv :: [in] Coordinate potential density [R ~> kg m-3]

  • h_in :: [in] Layer thicknesses [H ~> m or kg m-2]

  • dp_int :: [out] The change in interface positions [H ~> m or kg m-2]

Call to:

cushn mom_error_handler::mom_error

Called from:

hybgen_regrid