MOM6
regrid_consts Module Reference

Detailed Description

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...
 

Function/Subroutine Documentation

◆ coordinatemode()

integer function regrid_consts::coordinatemode ( character(len=*), intent(in)  string)

Parse a string parameter specifying the coordinate mode and return the appropriate enumerated integer.

Returns
Enumerated integer indicating coordinate mode
Parameters
[in]stringString to indicate coordinate mode

Definition at line 54 of file regrid_consts.F90.

54  integer :: coordinateMode !< Enumerated integer indicating coordinate mode
55  character(len=*), intent(in) :: string !< String to indicate coordinate mode
56  select case ( uppercase(trim(string)) )
57  case (trim(regridding_layer_string)); coordinatemode = regridding_layer
58  case (trim(regridding_zstar_string)); coordinatemode = regridding_zstar
59  case (trim(regridding_zstar_string_old)); coordinatemode = regridding_zstar
60  case (trim(regridding_rho_string)); coordinatemode = regridding_rho
61  case (trim(regridding_sigma_string)); coordinatemode = regridding_sigma
62  case (trim(regridding_hycom1_string)); coordinatemode = regridding_hycom1
63  case (trim(regridding_slight_string)); coordinatemode = regridding_slight
64  case (trim(regridding_arbitrary_string)); coordinatemode = regridding_arbitrary
65  case (trim(regridding_sigma_shelf_zstar_string)); coordinatemode = regridding_sigma_shelf_zstar
66  case (trim(regridding_adaptive_string)); coordinatemode = regridding_adaptive
67  case default ; call mom_error(fatal, "coordinateMode: "//&
68  "Unrecognized choice of coordinate ("//trim(string)//").")
69  end select

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().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ coordinateunitsi()

character(len=16) function regrid_consts::coordinateunitsi ( integer, intent(in)  coordMode)

Returns a string with the coordinate units associated with the enumerated integer,.

Returns
Units of coordinate
Parameters
[in]coordmodeCoordinate mode

Definition at line 75 of file regrid_consts.F90.

75  character(len=16) :: coordinateUnitsI !< Units of coordinate
76  integer, intent(in) :: coordMode !< Coordinate mode
77  select case ( coordmode )
78  case (regridding_layer); coordinateunitsi = "kg m^-3"
79  case (regridding_zstar); coordinateunitsi = "m"
80  case (regridding_sigma_shelf_zstar); coordinateunitsi = "m"
81  case (regridding_rho); coordinateunitsi = "kg m^-3"
82  case (regridding_sigma); coordinateunitsi = "Non-dimensional"
83  case (regridding_hycom1); coordinateunitsi = "m"
84  case (regridding_slight); coordinateunitsi = "m"
85  case (regridding_adaptive); coordinateunitsi = "m"
86  case default ; call mom_error(fatal, "coordinateUnts: "//&
87  "Unrecognized coordinate mode.")
88  end select

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().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ coordinateunitss()

character(len=16) function regrid_consts::coordinateunitss ( character(len=*), intent(in)  string)

Returns a string with the coordinate units associated with the string defining the coordinate mode.

Returns
Units of coordinate
Parameters
[in]stringCoordinate mode

Definition at line 94 of file regrid_consts.F90.

94  character(len=16) :: coordinateUnitsS !< Units of coordinate
95  character(len=*), intent(in) :: string !< Coordinate mode
96  integer :: coordMode
97  coordmode = coordinatemode(string)
98  coordinateunitss = coordinateunitsi(coordmode)

References coordinatemode(), and coordinateunitsi().

Here is the call graph for this function:

◆ state_dependent_char()

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.

Parameters
[in]stringString to indicate coordinate mode

Definition at line 103 of file regrid_consts.F90.

103  character(len=*), intent(in) :: string !< String to indicate coordinate mode
104 
105  state_dependent_char = state_dependent_int( coordinatemode(string) )
106 

References coordinatemode(), and state_dependent_int().

Here is the call graph for this function:

◆ 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.

Parameters
[in]modeCoordinate mode

Definition at line 111 of file regrid_consts.F90.

111  integer, intent(in) :: mode !< Coordinate mode
112  select case ( mode )
113  case (regridding_layer); state_dependent_int = .true.
114  case (regridding_zstar); state_dependent_int = .false.
115  case (regridding_sigma_shelf_zstar); state_dependent_int = .false.
116  case (regridding_rho); state_dependent_int = .true.
117  case (regridding_sigma); state_dependent_int = .false.
118  case (regridding_hycom1); state_dependent_int = .true.
119  case (regridding_slight); state_dependent_int = .true.
120  case (regridding_adaptive); state_dependent_int = .true.
121  case default ; call mom_error(fatal, "state_dependent: "//&
122  "Unrecognized choice of coordinate.")
123  end select

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().

Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ default_coordinate_mode

character(len=*), parameter regrid_consts::default_coordinate_mode = REGRIDDING_LAYER_STRING

Default coordinate mode.

Definition at line 35 of file regrid_consts.F90.

35 character(len=*), parameter :: DEFAULT_COORDINATE_MODE = REGRIDDING_LAYER_STRING !< Default coordinate mode

◆ regridding_adaptive

integer, parameter regrid_consts::regridding_adaptive = 9

◆ regridding_adaptive_string

character(len=*), parameter regrid_consts::regridding_adaptive_string = "ADAPTIVE"

Adaptive coordinate string.

Definition at line 34 of file regrid_consts.F90.

34 character(len=*), parameter :: REGRIDDING_ADAPTIVE_STRING = "ADAPTIVE" !< Adaptive coordinate string

Referenced by coordinatemode().

◆ regridding_arbitrary

integer, parameter regrid_consts::regridding_arbitrary = 5

Arbitrary coordinates identifier.

Definition at line 17 of file regrid_consts.F90.

17 integer, parameter :: REGRIDDING_ARBITRARY = 5 !< Arbitrary coordinates identifier

Referenced by coordinatemode().

◆ regridding_arbitrary_string

character(len=*), parameter regrid_consts::regridding_arbitrary_string = "ARB"

Arbitrary coordinates.

Definition at line 30 of file regrid_consts.F90.

30 character(len=*), parameter :: REGRIDDING_ARBITRARY_STRING = "ARB" !< Arbitrary coordinates

Referenced by coordinatemode().

◆ regridding_hycom1

◆ regridding_hycom1_string

character(len=*), parameter regrid_consts::regridding_hycom1_string = "HYCOM1"

Hycom string.

Definition at line 31 of file regrid_consts.F90.

31 character(len=*), parameter :: REGRIDDING_HYCOM1_STRING = "HYCOM1" !< Hycom string

Referenced by coordinatemode().

◆ regridding_layer

integer, parameter regrid_consts::regridding_layer = 1

Layer mode identifier.

Definition at line 13 of file regrid_consts.F90.

13 integer, parameter :: REGRIDDING_LAYER = 1 !< Layer mode identifier

Referenced by coordinatemode(), coordinateunitsi(), and state_dependent_int().

◆ regridding_layer_string

character(len=*), parameter regrid_consts::regridding_layer_string = "LAYER"

Layer string.

Definition at line 25 of file regrid_consts.F90.

25 character(len=*), parameter :: REGRIDDING_LAYER_STRING = "LAYER" !< Layer string

Referenced by coordinatemode().

◆ regridding_rho

◆ regridding_rho_string

character(len=*), parameter regrid_consts::regridding_rho_string = "RHO"

Rho string.

Definition at line 28 of file regrid_consts.F90.

28 character(len=*), parameter :: REGRIDDING_RHO_STRING = "RHO" !< Rho string

Referenced by coordinatemode().

◆ regridding_sigma

◆ regridding_sigma_shelf_zstar

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.

21 integer, parameter :: REGRIDDING_SIGMA_SHELF_ZSTAR = 8 !< Identifiered for z* coordinates at the bottom,

Referenced by coordinatemode(), coordinateunitsi(), isomip_initialization::isomip_initialize_sponges(), isomip_initialization::isomip_initialize_temperature_salinity(), isomip_initialization::isomip_initialize_thickness(), and state_dependent_int().

◆ regridding_sigma_shelf_zstar_string

character(len=*), parameter regrid_consts::regridding_sigma_shelf_zstar_string = "SIGMA_SHELF_ZSTAR"

Hybrid z*/sigma.

Definition at line 33 of file regrid_consts.F90.

33 character(len=*), parameter :: REGRIDDING_SIGMA_SHELF_ZSTAR_STRING = "SIGMA_SHELF_ZSTAR" !< Hybrid z*/sigma

Referenced by coordinatemode().

◆ regridding_sigma_string

character(len=*), parameter regrid_consts::regridding_sigma_string = "SIGMA"

Sigma string.

Definition at line 29 of file regrid_consts.F90.

29 character(len=*), parameter :: REGRIDDING_SIGMA_STRING = "SIGMA" !< Sigma string

Referenced by coordinatemode().

◆ regridding_slight

integer, parameter regrid_consts::regridding_slight = 7

◆ regridding_slight_string

character(len=*), parameter regrid_consts::regridding_slight_string = "SLIGHT"

Hybrid S-rho string.

Definition at line 32 of file regrid_consts.F90.

32 character(len=*), parameter :: REGRIDDING_SLIGHT_STRING = "SLIGHT" !< Hybrid S-rho string

Referenced by coordinatemode().

◆ regridding_zstar

integer, parameter regrid_consts::regridding_zstar = 2

z* coordinates identifier

Definition at line 14 of file regrid_consts.F90.

14 integer, parameter :: REGRIDDING_ZSTAR = 2 !< z* coordinates identifier

Referenced by coordinatemode(), coordinateunitsi(), and state_dependent_int().

◆ regridding_zstar_string

character(len=*), parameter regrid_consts::regridding_zstar_string = "ZSTAR"

z* string

Definition at line 27 of file regrid_consts.F90.

27 character(len=*), parameter :: REGRIDDING_ZSTAR_STRING = "ZSTAR" !< z* string

Referenced by coordinatemode().

◆ regridding_zstar_string_old

character(len=*), parameter regrid_consts::regridding_zstar_string_old = "Z*"

z* string (legacy name)

Definition at line 26 of file regrid_consts.F90.

26 character(len=*), parameter :: REGRIDDING_ZSTAR_STRING_OLD = "Z*" !< z* string (legacy name)

Referenced by coordinatemode().