MOM6
regrid_consts.F90
Go to the documentation of this file.
1 !> Contains constants for interpreting input parameters that control regridding.
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 use mom_error_handler, only : mom_error, fatal
8 
9 implicit none ; public
10 
11 integer, parameter :: regridding_num_types = 2
12 
13 ! List of regridding types. These should be consecutive and starting at 1.
14 ! This allows them to be used as array indices.
15 integer, parameter :: regridding_layer = 1 !< Layer mode
16 integer, parameter :: regridding_zstar = 2 !< z* coordinates
17 integer, parameter :: regridding_rho = 3 !< Target interface densities
18 integer, parameter :: regridding_sigma = 4 !< Sigma coordinates
19 integer, parameter :: regridding_arbitrary = 5 !< Arbitrary coordinates
20 integer, parameter :: regridding_hycom1 = 6 !< Simple HyCOM coordinates without BBL
21 integer, parameter :: regridding_slight = 7 !< Stretched coordinates in the
22 integer, parameter :: regridding_sigma_shelf_zstar = 8 !< z* coordinates at the bottom, sigma-near the top
23  !! lightest water, isopycnal below
24 integer, parameter :: regridding_adaptive = 9
25 
26 character(len=*), parameter :: regridding_layer_string = "LAYER" !< Layer string
27 character(len=*), parameter :: regridding_zstar_string_old = "Z*" !< z* string (legacy name)
28 character(len=*), parameter :: regridding_zstar_string = "ZSTAR" !< z* string
29 character(len=*), parameter :: regridding_rho_string = "RHO" !< Rho string
30 character(len=*), parameter :: regridding_sigma_string = "SIGMA" !< Sigma string
31 character(len=*), parameter :: regridding_arbitrary_string = "ARB" !< Arbitrary coordinates
32 character(len=*), parameter :: regridding_hycom1_string = "HYCOM1" !< Hycom string
33 character(len=*), parameter :: regridding_slight_string = "SLIGHT" !< Hybrid S-rho string
34 character(len=*), parameter :: regridding_sigma_shelf_zstar_string = "SIGMA_SHELF_ZSTAR" !< Hybrid z*/sigma
35 character(len=*), parameter :: regridding_adaptive_string = "ADAPTIVE"
36 character(len=*), parameter :: default_coordinate_mode = regridding_layer_string !< Default coordinate mode
37 
38 integer, dimension(REGRIDDING_NUM_TYPES), parameter :: vertical_coords = &
40  !(/ REGRIDDING_LAYER, REGRIDDING_ZSTAR, REGRIDDING_RHO, &
41  ! REGRIDDING_SIGMA, REGRIDDING_ARBITRARY, &
42  ! REGRIDDING_HYCOM1, REGRIDDING_SLIGHT /)
43 
44 character(len=*), dimension(REGRIDDING_NUM_TYPES), parameter :: vertical_coord_strings = &
46  !(/ REGRIDDING_LAYER_STRING, REGRIDDING_ZSTAR_STRING, REGRIDDING_RHO_STRING, &
47  ! REGRIDDING_SIGMA_STRING, REGRIDDING_ARBITRARY_STRING, &
48  ! REGRIDDING_HYCOM1_STRING, REGRIDDING_SLIGHT_STRING /)
49 
50 interface coordinateunits
51  module procedure coordinateunitsi
52  module procedure coordinateunitss
53 end interface
54 
55 interface state_dependent
56  module procedure state_dependent_char
57  module procedure state_dependent_int
58 end interface
59 
60 contains
61 
62 !> Parse a string parameter specifying the coordinate mode and
63 !! return the appropriate enumerated integer
64 function coordinatemode(string)
65  integer :: coordinateMode !< Enumerated integer indicating coordinate mode
66  character(len=*), intent(in) :: string !< String to indicate coordinate mode
67  select case ( uppercase(trim(string)) )
68  case (trim(regridding_layer_string)); coordinatemode = regridding_layer
69  case (trim(regridding_zstar_string)); coordinatemode = regridding_zstar
70  case (trim(regridding_zstar_string_old)); coordinatemode = regridding_zstar
71  case (trim(regridding_rho_string)); coordinatemode = regridding_rho
72  case (trim(regridding_sigma_string)); coordinatemode = regridding_sigma
73  case (trim(regridding_hycom1_string)); coordinatemode = regridding_hycom1
74  case (trim(regridding_slight_string)); coordinatemode = regridding_slight
75  case (trim(regridding_arbitrary_string)); coordinatemode = regridding_arbitrary
77  case (trim(regridding_adaptive_string)); coordinatemode = regridding_adaptive
78  case default ; call mom_error(fatal, "coordinateMode: "//&
79  "Unrecognized choice of coordinate ("//trim(string)//").")
80  end select
81 end function coordinatemode
82 
83 !> Returns a string with the coordinate units associated with the
84 !! enumerated integer,
85 function coordinateunitsi(coordMode)
86  character(len=16) :: coordinateUnitsI !< Units of coordinate
87  integer, intent(in) :: coordMode !< Coordinate mode
88  select case ( coordmode )
89  case (regridding_layer); coordinateunitsi = "kg m^-3"
90  case (regridding_zstar); coordinateunitsi = "m"
91  case (regridding_sigma_shelf_zstar); coordinateunitsi = "m"
92  case (regridding_rho); coordinateunitsi = "kg m^-3"
93  case (regridding_sigma); coordinateunitsi = "Non-dimensional"
94  case (regridding_hycom1); coordinateunitsi = "m"
95  case (regridding_slight); coordinateunitsi = "m"
96  case (regridding_adaptive); coordinateunitsi = "m"
97  case default ; call mom_error(fatal, "coordinateUnts: "//&
98  "Unrecognized coordinate mode.")
99  end select
100 end function coordinateunitsi
101 
102 !> Returns a string with the coordinate units associated with the
103 !! string defining the coordinate mode.
104 function coordinateunitss(string)
105  character(len=16) :: coordinateUnitsS !< Units of coordinate
106  character(len=*), intent(in) :: string !< Coordinate mode
107  integer :: coordMode
108  coordmode = coordinatemode(string)
109  coordinateunitss = coordinateunitsi(coordmode)
110 end function coordinateunitss
111 
112 !> Returns true if the coordinate is dependent on the state density, returns false otherwise.
113 logical function state_dependent_char(string)
114  character(len=*), intent(in) :: string !< String to indicate coordinate mode
115 
117 
118 end function state_dependent_char
119 
120 !> Returns true if the coordinate is dependent on the state density, returns false otherwise.
121 logical function state_dependent_int(mode)
122  integer, intent(in) :: mode !< Coordinate mode
123  select case ( mode )
124  case (regridding_layer); state_dependent_int = .true.
125  case (regridding_zstar); state_dependent_int = .false.
127  case (regridding_rho); state_dependent_int = .true.
128  case (regridding_sigma); state_dependent_int = .false.
129  case (regridding_hycom1); state_dependent_int = .true.
130  case (regridding_slight); state_dependent_int = .true.
132  case default ; call mom_error(fatal, "state_dependent: "//&
133  "Unrecognized choice of coordinate.")
134  end select
135 end function state_dependent_int
136 
137 end module regrid_consts
integer, parameter regridding_layer
Layer mode.
integer function coordinatemode(string)
Parse a string parameter specifying the coordinate mode and return the appropriate enumerated integer...
integer, parameter regridding_sigma
Sigma coordinates.
character(len= *), parameter regridding_sigma_string
Sigma string.
character(len= *), parameter regridding_rho_string
Rho string.
character(len= *), parameter regridding_zstar_string
z* string
logical function state_dependent_char(string)
Returns true if the coordinate is dependent on the state density, returns false otherwise.
integer, dimension(regridding_num_types), parameter vertical_coords
character(len= *), parameter default_coordinate_mode
Default coordinate mode.
character(len= *), parameter regridding_layer_string
Layer string.
integer, parameter regridding_hycom1
Simple HyCOM coordinates without BBL.
character(len=16) function coordinateunitsi(coordMode)
Returns a string with the coordinate units associated with the enumerated integer,.
character(len=len(input_string)) function, public uppercase(input_string)
logical function state_dependent_int(mode)
Returns true if the coordinate is dependent on the state density, returns false otherwise.
integer, parameter regridding_adaptive
character(len=16) function coordinateunitss(string)
Returns a string with the coordinate units associated with the string defining the coordinate mode...
character(len= *), parameter regridding_adaptive_string
integer, parameter regridding_zstar
z* coordinates
character(len= *), parameter regridding_sigma_shelf_zstar_string
Hybrid z*/sigma.
integer, parameter regridding_num_types
integer, parameter regridding_sigma_shelf_zstar
z* coordinates at the bottom, sigma-near the top lightest water, isopycnal below
character(len= *), parameter regridding_arbitrary_string
Arbitrary coordinates.
Contains constants for interpreting input parameters that control regridding.
character(len= *), parameter regridding_hycom1_string
Hycom string.
integer, parameter regridding_arbitrary
Arbitrary coordinates.
character(len= *), parameter regridding_slight_string
Hybrid S-rho string.
subroutine, public mom_error(level, message, all_print)
integer, parameter regridding_rho
Target interface densities.
character(len= *), parameter regridding_zstar_string_old
z* string (legacy name)
character(len= *), dimension(regridding_num_types), parameter vertical_coord_strings
integer, parameter regridding_slight
Stretched coordinates in the.