MOM6
|
Contains constants for interpreting input parameters that control regridding.
Data Types | |
interface | coordinateunits |
Returns a string with the coordinate units associated with the coordinate mode. More... | |
interface | state_dependent |
Returns true if the coordinate is dependent on the state density, returns false otherwise. More... | |
Functions/Subroutines | |
integer function | coordinatemode (string) |
Parse a string parameter specifying the coordinate mode and return the appropriate enumerated integer. More... | |
character(len=16) function | coordinateunitsi (coordMode) |
Returns a string with the coordinate units associated with the enumerated integer,. More... | |
character(len=16) function | coordinateunitss (string) |
Returns a string with the coordinate units associated with the string defining the coordinate mode. More... | |
logical function | state_dependent_char (string) |
Returns true if the coordinate is dependent on the state density, returns false otherwise. More... | |
logical function | state_dependent_int (mode) |
Returns true if the coordinate is dependent on the state density, returns false otherwise. More... | |
Variables | |
integer, parameter | regridding_layer = 1 |
Layer mode identifier. More... | |
integer, parameter | regridding_zstar = 2 |
z* coordinates identifier More... | |
integer, parameter | regridding_rho = 3 |
Density coordinates identifier. More... | |
integer, parameter | regridding_sigma = 4 |
Sigma coordinates identifier. More... | |
integer, parameter | regridding_arbitrary = 5 |
Arbitrary coordinates identifier. More... | |
integer, parameter | regridding_hycom1 = 6 |
Simple HyCOM coordinates without BBL. More... | |
integer, parameter | regridding_slight = 7 |
Identifier for stretched coordinates in the lightest water, isopycnal below. More... | |
integer, parameter | regridding_sigma_shelf_zstar = 8 |
Identifiered for z* coordinates at the bottom, sigma-near the top. More... | |
integer, parameter | regridding_adaptive = 9 |
Adaptive coordinate mode identifier. More... | |
character(len= *), parameter | regridding_layer_string = "LAYER" |
Layer string. More... | |
character(len= *), parameter | regridding_zstar_string_old = "Z*" |
z* string (legacy name) More... | |
character(len= *), parameter | regridding_zstar_string = "ZSTAR" |
z* string More... | |
character(len= *), parameter | regridding_rho_string = "RHO" |
Rho string. More... | |
character(len= *), parameter | regridding_sigma_string = "SIGMA" |
Sigma string. More... | |
character(len= *), parameter | regridding_arbitrary_string = "ARB" |
Arbitrary coordinates. More... | |
character(len= *), parameter | regridding_hycom1_string = "HYCOM1" |
Hycom string. More... | |
character(len= *), parameter | regridding_slight_string = "SLIGHT" |
Hybrid S-rho string. More... | |
character(len= *), parameter | regridding_sigma_shelf_zstar_string = "SIGMA_SHELF_ZSTAR" |
Hybrid z*/sigma. More... | |
character(len= *), parameter | regridding_adaptive_string = "ADAPTIVE" |
Adaptive coordinate string. More... | |
character(len= *), parameter | default_coordinate_mode = REGRIDDING_LAYER_STRING |
Default coordinate mode. More... | |
integer function regrid_consts::coordinatemode | ( | character(len=*), intent(in) | string | ) |
Parse a string parameter specifying the coordinate mode and return the appropriate enumerated integer.
[in] | string | String to indicate coordinate mode |
Definition at line 54 of file regrid_consts.F90.
References mom_error_handler::mom_error(), regridding_adaptive, regridding_adaptive_string, regridding_arbitrary, regridding_arbitrary_string, regridding_hycom1, regridding_hycom1_string, regridding_layer, regridding_layer_string, regridding_rho, regridding_rho_string, regridding_sigma, regridding_sigma_shelf_zstar, regridding_sigma_shelf_zstar_string, regridding_sigma_string, regridding_slight, regridding_slight_string, regridding_zstar, regridding_zstar_string, regridding_zstar_string_old, and mom_string_functions::uppercase().
Referenced by coordinateunitss(), mom_diag_remap::diag_remap_calc_hmask(), mom_diag_remap::diag_remap_configure_axes(), mom_diag_remap::diag_remap_init(), mom_diag_remap::diag_remap_update(), and state_dependent_char().
character(len=16) function regrid_consts::coordinateunitsi | ( | integer, intent(in) | coordMode | ) |
Returns a string with the coordinate units associated with the enumerated integer,.
[in] | coordmode | Coordinate mode |
Definition at line 75 of file regrid_consts.F90.
References mom_error_handler::mom_error(), regridding_adaptive, regridding_hycom1, regridding_layer, regridding_rho, regridding_sigma, regridding_sigma_shelf_zstar, regridding_slight, and regridding_zstar.
Referenced by coordinateunitss().
Returns a string with the coordinate units associated with the string defining the coordinate mode.
[in] | string | Coordinate mode |
Definition at line 94 of file regrid_consts.F90.
References coordinatemode(), and coordinateunitsi().
logical function regrid_consts::state_dependent_char | ( | character(len=*), intent(in) | string | ) |
Returns true if the coordinate is dependent on the state density, returns false otherwise.
[in] | string | String to indicate coordinate mode |
Definition at line 103 of file regrid_consts.F90.
References coordinatemode(), and state_dependent_int().
logical function regrid_consts::state_dependent_int | ( | integer, intent(in) | mode | ) |
Returns true if the coordinate is dependent on the state density, returns false otherwise.
[in] | mode | Coordinate mode |
Definition at line 111 of file regrid_consts.F90.
References mom_error_handler::mom_error(), regridding_adaptive, regridding_hycom1, regridding_layer, regridding_rho, regridding_sigma, regridding_sigma_shelf_zstar, regridding_slight, and regridding_zstar.
Referenced by state_dependent_char().
character(len=*), parameter regrid_consts::default_coordinate_mode = REGRIDDING_LAYER_STRING |
Default coordinate mode.
Definition at line 35 of file regrid_consts.F90.
integer, parameter regrid_consts::regridding_adaptive = 9 |
Adaptive coordinate mode identifier.
Definition at line 23 of file regrid_consts.F90.
Referenced by coordinatemode(), coordinateunitsi(), mom_regridding::getcoordinateshortname(), mom_regridding::getcoordinateunits(), mom_regridding::getstaticthickness(), mom_regridding::initcoord(), mom_regridding::initialize_regridding(), mom_regridding::regridding_main(), mom_regridding::set_regrid_params(), state_dependent_int(), and mom_regridding::uniformresolution().
character(len=*), parameter regrid_consts::regridding_adaptive_string = "ADAPTIVE" |
Adaptive coordinate string.
Definition at line 34 of file regrid_consts.F90.
Referenced by coordinatemode().
integer, parameter regrid_consts::regridding_arbitrary = 5 |
Arbitrary coordinates identifier.
Definition at line 17 of file regrid_consts.F90.
Referenced by coordinatemode().
character(len=*), parameter regrid_consts::regridding_arbitrary_string = "ARB" |
Arbitrary coordinates.
Definition at line 30 of file regrid_consts.F90.
Referenced by coordinatemode().
integer, parameter regrid_consts::regridding_hycom1 = 6 |
Simple HyCOM coordinates without BBL.
Definition at line 18 of file regrid_consts.F90.
Referenced by coordinatemode(), coordinateunitsi(), mom_regridding::getcoordinateshortname(), mom_regridding::getcoordinateunits(), mom_regridding::getstaticthickness(), mom_regridding::initcoord(), mom_regridding::initialize_regridding(), mom_regridding::regridding_main(), mom_regridding::set_regrid_max_depths(), mom_regridding::set_regrid_max_thickness(), mom_regridding::set_regrid_params(), state_dependent_int(), and mom_regridding::uniformresolution().
character(len=*), parameter regrid_consts::regridding_hycom1_string = "HYCOM1" |
integer, parameter regrid_consts::regridding_layer = 1 |
Layer mode identifier.
Definition at line 13 of file regrid_consts.F90.
Referenced by coordinatemode(), coordinateunitsi(), and state_dependent_int().
character(len=*), parameter regrid_consts::regridding_layer_string = "LAYER" |
integer, parameter regrid_consts::regridding_rho = 3 |
Density coordinates identifier.
Definition at line 15 of file regrid_consts.F90.
Referenced by adjustment_initialization::adjustment_initialize_temperature_salinity(), adjustment_initialization::adjustment_initialize_thickness(), coordinatemode(), coordinateunitsi(), dome2d_initialization::dome2d_initialize_temperature_salinity(), dome2d_initialization::dome2d_initialize_thickness(), dumbbell_initialization::dumbbell_initialize_thickness(), rossby_front_2d_initialization::rossby_front_initialize_thickness(), seamount_initialization::seamount_initialize_temperature_salinity(), seamount_initialization::seamount_initialize_thickness(), and state_dependent_int().
character(len=*), parameter regrid_consts::regridding_rho_string = "RHO" |
integer, parameter regrid_consts::regridding_sigma = 4 |
Sigma coordinates identifier.
Definition at line 16 of file regrid_consts.F90.
Referenced by adjustment_initialization::adjustment_initialize_temperature_salinity(), adjustment_initialization::adjustment_initialize_thickness(), coordinatemode(), coordinateunitsi(), dome2d_initialization::dome2d_initialize_temperature_salinity(), dome2d_initialization::dome2d_initialize_thickness(), dumbbell_initialization::dumbbell_initialize_thickness(), rossby_front_2d_initialization::rossby_front_initialize_thickness(), seamount_initialization::seamount_initialize_temperature_salinity(), seamount_initialization::seamount_initialize_thickness(), and state_dependent_int().
integer, parameter regrid_consts::regridding_sigma_shelf_zstar = 8 |
Identifiered for z* coordinates at the bottom, sigma-near the top.
Definition at line 21 of file regrid_consts.F90.
Referenced by coordinatemode(), coordinateunitsi(), isomip_initialization::isomip_initialize_sponges(), isomip_initialization::isomip_initialize_temperature_salinity(), isomip_initialization::isomip_initialize_thickness(), and state_dependent_int().
character(len=*), parameter regrid_consts::regridding_sigma_shelf_zstar_string = "SIGMA_SHELF_ZSTAR" |
character(len=*), parameter regrid_consts::regridding_sigma_string = "SIGMA" |
integer, parameter regrid_consts::regridding_slight = 7 |
Identifier for stretched coordinates in the lightest water, isopycnal below.
Definition at line 19 of file regrid_consts.F90.
Referenced by coordinatemode(), coordinateunitsi(), mom_regridding::getcoordinateshortname(), mom_regridding::getcoordinateunits(), mom_regridding::getstaticthickness(), mom_regridding::initcoord(), mom_regridding::initialize_regridding(), mom_regridding::regridding_main(), mom_regridding::set_regrid_max_depths(), mom_regridding::set_regrid_max_thickness(), mom_regridding::set_regrid_params(), state_dependent_int(), and mom_regridding::uniformresolution().
character(len=*), parameter regrid_consts::regridding_slight_string = "SLIGHT" |
Hybrid S-rho string.
Definition at line 32 of file regrid_consts.F90.
Referenced by coordinatemode().
integer, parameter regrid_consts::regridding_zstar = 2 |
z* coordinates identifier
Definition at line 14 of file regrid_consts.F90.
Referenced by coordinatemode(), coordinateunitsi(), and state_dependent_int().
character(len=*), parameter regrid_consts::regridding_zstar_string = "ZSTAR" |
character(len=*), parameter regrid_consts::regridding_zstar_string_old = "Z*" |
z* string (legacy name)
Definition at line 26 of file regrid_consts.F90.
Referenced by coordinatemode().