mom_hybgen_unmix module reference
This module contains the hybgen unmixing routines from HYCOM, with modifications to follow the MOM6 coding conventions and several bugs fixed.
Data Types
Control structure containing required parameters for the hybgen coordinate generator. |
Functions/Subroutines
Initialise a hybgen_unmix_CS control structure and store its parameters. |
|
This subroutine deallocates memory in the control structure for the hybgen unmixing module. |
|
This subroutine can be used to set the parameters for the hybgen module. |
|
Unmix the properties in the lowest layer with mass if it is too light, and make any other changes to the water column to prepare for regridding. |
|
Unmix the properties in the lowest layer if it is too light. |
Detailed Description
This module contains the hybgen unmixing routines from HYCOM, with modifications to follow the MOM6 coding conventions and several bugs fixed.
Type Documentation
- type mom_hybgen_unmix/hybgen_unmix_cs
Control structure containing required parameters for the hybgen coordinate generator.
- Type fields:
% nsigma ::
integerNumber of sigma levels used by HYBGEN.% hybiso ::
realHybgen uses PCM if layer is within hybiso of target density [R ~> kg m-3].% dp00i ::
realDeep isopycnal spacing minimum thickness [H ~> m or kg m-2].% qhybrlx ::
realHybgen relaxation amount per thermodynamic time steps [nondim].% dp0k ::
real, dimension(:), allocatableminimum deep z-layer separation [H ~> m or kg m-2]% ds0k ::
real, dimension(:), allocatableminimum shallow z-layer separation [H ~> m or kg m-2]% dpns ::
realdepth to start terrain following [H ~> m or kg m-2]% dsns ::
realdepth to stop terrain following [H ~> m or kg m-2]% min_dilate ::
realThe minimum amount of dilation that is permitted when converting target coordinates from z to z* [nondim]. This limit applies when wetting occurs.% max_dilate ::
realThe maximum amount of dilation that is permitted when converting target coordinates from z to z* [nondim]. This limit applies when drying occurs.% topiso_const ::
realShallowest depth for isopycnal layers [H ~> m or kg m-2].% ref_pressure ::
realReference pressure for density calculations [R L2 T-2 ~> Pa].% target_density ::
real, dimension(:), allocatableNominal density of interfaces [R ~> kg m-3].
Function/Subroutine Documentation
- subroutine mom_hybgen_unmix/init_hybgen_unmix(CS, GV, US, param_file, hybgen_regridCS)
Initialise a hybgen_unmix_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 ::
param_file[in] Parameter filehybgen_regridcs :: Control structure for hybgen regridding for sharing parameters.
- Call to:
mom_hybgen_regrid::get_hybgen_regrid_paramsmom_error_handler::mom_error- Called from:
- subroutine mom_hybgen_unmix/end_hybgen_unmix(CS)
This subroutine deallocates memory in the control structure for the hybgen unmixing module.
- Parameters:
cs :: Coordinate control structure
- Called from:
- subroutine mom_hybgen_unmix/set_hybgen_unmix_params(CS, min_thickness)
This subroutine can be used to set the parameters for the hybgen module.
- Parameters:
cs :: Coordinate unmixing control structure
min_thickness ::
min_thickness[in] Minimum allowed thickness [H ~> m or kg m-2]
- Call to:
- subroutine mom_hybgen_unmix/hybgen_unmix(G, GV, US, CS, tv, Reg, ntr, h)
Unmix the properties in the lowest layer with mass if it is too light, and make any other changes to the water column to prepare for regridding.
- Parameters:
g :: [in] Ocean grid structure
gv :: [in] Ocean vertical grid structure
us :: [in] A dimensional unit scaling type
cs :: [in] hybgen control structure
tv ::
tv[inout] Thermodynamics structurereg :: Tracer registry structure
ntr ::
ntr[in] The number of tracers in the registry, or 0 if the registry is not in use.h ::
h[inout] Layer thicknesses [H ~> m or kg m-2]
- Call to:
mom_interface_heights::calc_derived_thermohybgen_column_unmixmom_error_handler::mom_error- Called from:
- subroutine mom_hybgen_unmix/hybgen_column_unmix(CS, nk, Rcv_tgt, temp, saln, Rcv, eqn_of_state, ntr, tracer, trcflg, fixlay, qhrlx, h_col, terrain_following, h_thin)
Unmix the properties in the lowest layer if it is too light.
- Parameters:
cs :: [in] hybgen unmixing control structure
nk ::
nk[in] The number of layersfixlay ::
fixlay[in] deepest fixed coordinate layerqhrlx ::
qhrlx[in] Relaxation fraction per timestep [nondim], < 1.rcv_tgt :: [in] Target potential density [R ~> kg m-3]
temp ::
temp[inout] A column of potential temperature [C ~> degC]saln ::
saln[inout] A column of salinity [S ~> ppt]rcv :: [inout] Coordinate potential density [R ~> kg m-3]
eqn_of_state ::
eqn_of_state[in] Equation of state structurentr ::
ntr[in] The number of registered passive tracerstracer ::
tracer[inout] Columns of the passive tracers [Conc]trcflg ::
trcflg[in] Hycom tracer type flag for each tracerh_col ::
h_col[inout] Layer thicknesses [H ~> m or kg m-2]terrain_following ::
terrain_following[in] True if this column is terrain followingh_thin ::
h_thin[in] A negligibly small thickness to identify essentially vanished layers [H ~> m or kg m-2]
- Called from: