MOM6
regrid_consts.F90
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
7 use mom_string_functions, only : uppercase
8 
9 implicit none ; public
10 
11 ! List of regridding types. These should be consecutive and starting at 1.
12 ! This allows them to be used as array indices.
13 integer, parameter :: regridding_layer = 1 !< Layer mode identifier
14 integer, parameter :: regridding_zstar = 2 !< z* coordinates identifier
15 integer, parameter :: regridding_rho = 3 !< Density coordinates identifier
16 integer, parameter :: regridding_sigma = 4 !< Sigma coordinates identifier
17 integer, parameter :: regridding_arbitrary = 5 !< Arbitrary coordinates identifier
18 integer, parameter :: regridding_hycom1 = 6 !< Simple HyCOM coordinates without BBL
19 integer, parameter :: regridding_slight = 7 !< Identifier for stretched coordinates in the
20  !! lightest water, isopycnal below
21 integer, parameter :: regridding_sigma_shelf_zstar = 8 !< Identifiered for z* coordinates at the bottom,
22  !! sigma-near the top
23 integer, parameter :: regridding_adaptive = 9 !< Adaptive coordinate mode identifier
24 
25 character(len=*), parameter :: regridding_layer_string = "LAYER" !< Layer string
26 character(len=*), parameter :: regridding_zstar_string_old = "Z*" !< z* string (legacy name)
27 character(len=*), parameter :: regridding_zstar_string = "ZSTAR" !< z* string
28 character(len=*), parameter :: regridding_rho_string = "RHO" !< Rho string
29 character(len=*), parameter :: regridding_sigma_string = "SIGMA" !< Sigma string
30 character(len=*), parameter :: regridding_arbitrary_string = "ARB" !< Arbitrary coordinates
31 character(len=*), parameter :: regridding_hycom1_string = "HYCOM1" !< Hycom string
32 character(len=*), parameter :: regridding_slight_string = "SLIGHT" !< Hybrid S-rho string
33 character(len=*), parameter :: regridding_sigma_shelf_zstar_string = "SIGMA_SHELF_ZSTAR" !< Hybrid z*/sigma
34 character(len=*), parameter :: regridding_adaptive_string = "ADAPTIVE" !< Adaptive coordinate string
35 character(len=*), parameter :: default_coordinate_mode = regridding_layer_string !< Default coordinate mode
36 
37 !> Returns a string with the coordinate units associated with the coordinate mode.
38 interface coordinateunits
39  module procedure coordinateunitsi
40  module procedure coordinateunitss
41 end interface
42 
43 !> Returns true if the coordinate is dependent on the state density, returns false otherwise.
44 interface state_dependent
45  module procedure state_dependent_char
46  module procedure state_dependent_int
47 end interface
48 
49 contains
50 
51 !> Parse a string parameter specifying the coordinate mode and
52 !! return the appropriate enumerated integer
53 function coordinatemode(string)
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
70 end function coordinatemode
71 
72 !> Returns a string with the coordinate units associated with the
73 !! enumerated integer,
74 function coordinateunitsi(coordMode)
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
89 end function coordinateunitsi
90 
91 !> Returns a string with the coordinate units associated with the
92 !! string defining the coordinate mode.
93 function coordinateunitss(string)
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)
99 end function coordinateunitss
100 
101 !> Returns true if the coordinate is dependent on the state density, returns false otherwise.
102 logical function state_dependent_char(string)
103  character(len=*), intent(in) :: string !< String to indicate coordinate mode
104 
105  state_dependent_char = state_dependent_int( coordinatemode(string) )
106 
107 end function state_dependent_char
108 
109 !> Returns true if the coordinate is dependent on the state density, returns false otherwise.
110 logical function state_dependent_int(mode)
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
124 end function state_dependent_int
125 
126 end module regrid_consts
regrid_consts::state_dependent
Returns true if the coordinate is dependent on the state density, returns false otherwise.
Definition: regrid_consts.F90:44
mom_string_functions
Handy functions for manipulating strings.
Definition: MOM_string_functions.F90:2
regrid_consts::coordinateunits
Returns a string with the coordinate units associated with the coordinate mode.
Definition: regrid_consts.F90:38
regrid_consts
Contains constants for interpreting input parameters that control regridding.
Definition: regrid_consts.F90:2
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2