MOM6
mom_verticalgrid Module Reference

Detailed Description

Provides a transparent vertical ocean grid type and supporting routines.

Data Types

type  verticalgrid_type
 Describes the vertical ocean grid, including unit conversion factors. More...
 

Functions/Subroutines

subroutine, public verticalgridinit (param_file, GV, US)
 Allocates and initializes the ocean model vertical grid structure. More...
 
subroutine, public fix_restart_scaling (GV)
 Set the scaling factors for restart files to the scaling factors for this run. More...
 
character(len=48) function, public get_thickness_units (GV)
 Returns the model's thickness units, usually m or kg/m^2. More...
 
character(len=48) function, public get_flux_units (GV)
 Returns the model's thickness flux units, usually m^3/s or kg/s. More...
 
character(len=48) function, public get_tr_flux_units (GV, tr_units, tr_vol_conc_units, tr_mass_conc_units)
 Returns the model's tracer flux units. More...
 
subroutine, public setverticalgridaxes (Rlay, GV)
 This sets the coordinate data for the "layer mode" of the isopycnal model. More...
 
subroutine, public verticalgridend (GV)
 Deallocates the model's vertical grid structure. More...
 

Function/Subroutine Documentation

◆ fix_restart_scaling()

subroutine, public mom_verticalgrid::fix_restart_scaling ( type(verticalgrid_type), intent(inout)  GV)

Set the scaling factors for restart files to the scaling factors for this run.

Parameters
[in,out]gvThe ocean's vertical grid structure

Definition at line 174 of file MOM_verticalGrid.F90.

174  type(verticalGrid_type), intent(inout) :: GV !< The ocean's vertical grid structure
175 
176  gv%m_to_H_restart = gv%m_to_H

◆ get_flux_units()

character(len=48) function, public mom_verticalgrid::get_flux_units ( type(verticalgrid_type), intent(in)  GV)

Returns the model's thickness flux units, usually m^3/s or kg/s.

Returns
The thickness flux units
Parameters
[in]gvThe ocean's vertical grid structure

Definition at line 196 of file MOM_verticalGrid.F90.

196  character(len=48) :: get_flux_units !< The thickness flux units
197  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
198  ! This subroutine returns the appropriate units for thickness fluxes,
199  ! depending on whether the model is Boussinesq or not and the scaling for
200  ! the vertical thickness.
201 
202  if (gv%Boussinesq) then
203  get_flux_units = "m3 s-1"
204  else
205  get_flux_units = "kg s-1"
206  endif

◆ get_thickness_units()

character(len=48) function, public mom_verticalgrid::get_thickness_units ( type(verticalgrid_type), intent(in)  GV)

Returns the model's thickness units, usually m or kg/m^2.

Returns
The vertical thickness units
Parameters
[in]gvThe ocean's vertical grid structure

Definition at line 181 of file MOM_verticalGrid.F90.

181  character(len=48) :: get_thickness_units !< The vertical thickness units
182  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
183  ! This subroutine returns the appropriate units for thicknesses,
184  ! depending on whether the model is Boussinesq or not and the scaling for
185  ! the vertical thickness.
186 
187  if (gv%Boussinesq) then
188  get_thickness_units = "m"
189  else
190  get_thickness_units = "kg m-2"
191  endif

◆ get_tr_flux_units()

character(len=48) function, public mom_verticalgrid::get_tr_flux_units ( type(verticalgrid_type), intent(in)  GV,
character(len=*), intent(in), optional  tr_units,
character(len=*), intent(in), optional  tr_vol_conc_units,
character(len=*), intent(in), optional  tr_mass_conc_units 
)

Returns the model's tracer flux units.

Returns
The model's flux units for a tracer.
Parameters
[in]gvThe ocean's vertical grid structure.
[in]tr_unitsUnits for a tracer, for example Celsius or PSU.
[in]tr_vol_conc_unitsThe concentration units per unit volume, for example if the units are umol m-3, tr_vol_conc_units would be umol.
[in]tr_mass_conc_unitsThe concentration units per unit mass of sea water, for example if the units are mol kg-1, tr_vol_conc_units would be mol.

Definition at line 211 of file MOM_verticalGrid.F90.

211  character(len=48) :: get_tr_flux_units !< The model's flux units
212  !! for a tracer.
213  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical
214  !! grid structure.
215  character(len=*), optional, intent(in) :: tr_units !< Units for a tracer, for example
216  !! Celsius or PSU.
217  character(len=*), optional, intent(in) :: tr_vol_conc_units !< The concentration units per unit
218  !! volume, for example if the units are
219  !! umol m-3, tr_vol_conc_units would
220  !! be umol.
221  character(len=*), optional, intent(in) :: tr_mass_conc_units !< The concentration units per unit
222  !! mass of sea water, for example if
223  !! the units are mol kg-1,
224  !! tr_vol_conc_units would be mol.
225 
226  ! This subroutine returns the appropriate units for thicknesses and fluxes,
227  ! depending on whether the model is Boussinesq or not and the scaling for
228  ! the vertical thickness.
229  integer :: cnt
230 
231  cnt = 0
232  if (present(tr_units)) cnt = cnt+1
233  if (present(tr_vol_conc_units)) cnt = cnt+1
234  if (present(tr_mass_conc_units)) cnt = cnt+1
235 
236  if (cnt == 0) call mom_error(fatal, "get_tr_flux_units: One of the three "//&
237  "arguments tr_units, tr_vol_conc_units, or tr_mass_conc_units "//&
238  "must be present.")
239  if (cnt > 1) call mom_error(fatal, "get_tr_flux_units: Only one of "//&
240  "tr_units, tr_vol_conc_units, and tr_mass_conc_units may be present.")
241  if (present(tr_units)) then
242  if (gv%Boussinesq) then
243  get_tr_flux_units = trim(tr_units)//" m3 s-1"
244  else
245  get_tr_flux_units = trim(tr_units)//" kg s-1"
246  endif
247  endif
248  if (present(tr_vol_conc_units)) then
249  if (gv%Boussinesq) then
250  get_tr_flux_units = trim(tr_vol_conc_units)//" s-1"
251  else
252  get_tr_flux_units = trim(tr_vol_conc_units)//" m-3 kg s-1"
253  endif
254  endif
255  if (present(tr_mass_conc_units)) then
256  if (gv%Boussinesq) then
257  get_tr_flux_units = trim(tr_mass_conc_units)//" kg-1 m3 s-1"
258  else
259  get_tr_flux_units = trim(tr_mass_conc_units)//" s-1"
260  endif
261  endif
262 

◆ setverticalgridaxes()

subroutine, public mom_verticalgrid::setverticalgridaxes ( real, dimension(gv%ke), intent(in)  Rlay,
type(verticalgrid_type), intent(inout)  GV 
)

This sets the coordinate data for the "layer mode" of the isopycnal model.

Parameters
[in,out]gvThe container for vertical grid data
[in]rlayThe layer target density

Definition at line 267 of file MOM_verticalGrid.F90.

267  type(verticalGrid_type), intent(inout) :: GV !< The container for vertical grid data
268  real, dimension(GV%ke), intent(in) :: Rlay !< The layer target density
269  ! Local variables
270  integer :: k, nk
271 
272  nk = gv%ke
273 
274  gv%zAxisLongName = 'Target Potential Density'
275  gv%zAxisUnits = 'kg m-3'
276  do k=1,nk ; gv%sLayer(k) = rlay(k) ; enddo
277  if (nk > 1) then
278  gv%sInterface(1) = 1.5*rlay(1) - 0.5*rlay(2)
279  do k=2,nk ; gv%sInterface(k) = 0.5*( rlay(k-1) + rlay(k) ) ; enddo
280  gv%sInterface(nk+1) = 1.5*rlay(nk) - 0.5*rlay(nk-1)
281  else
282  gv%sInterface(1) = 0.0 ; gv%sInterface(nk+1) = 2.0*rlay(nk)
283  endif
284 

◆ verticalgridend()

subroutine, public mom_verticalgrid::verticalgridend ( type(verticalgrid_type), pointer  GV)

Deallocates the model's vertical grid structure.

Parameters
gvThe ocean's vertical grid structure

Definition at line 289 of file MOM_verticalGrid.F90.

289  type(verticalGrid_type), pointer :: GV !< The ocean's vertical grid structure
290 
291  deallocate( gv%g_prime, gv%Rlay )
292  deallocate( gv%sInterface , gv%sLayer )
293  deallocate( gv )
294 

◆ verticalgridinit()

subroutine, public mom_verticalgrid::verticalgridinit ( type(param_file_type), intent(in)  param_file,
type(verticalgrid_type), pointer  GV,
type(unit_scale_type), intent(in)  US 
)

Allocates and initializes the ocean model vertical grid structure.

Parameters
[in]param_fileParameter file handle/type
gvThe container for vertical grid data
[in]usA dimensional unit scaling type

Definition at line 72 of file MOM_verticalGrid.F90.

72  type(param_file_type), intent(in) :: param_file !< Parameter file handle/type
73  type(verticalGrid_type), pointer :: GV !< The container for vertical grid data
74  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
75  ! This routine initializes the verticalGrid_type structure (GV).
76  ! All memory is allocated but not necessarily set to meaningful values until later.
77 
78  ! Local variables
79  integer :: nk, H_power
80  real :: H_rescale_factor
81  ! This include declares and sets the variable "version".
82 # include "version_variable.h"
83  character(len=16) :: mdl = 'MOM_verticalGrid'
84 
85  if (associated(gv)) call mom_error(fatal, &
86  'verticalGridInit: called with an associated GV pointer.')
87  allocate(gv)
88 
89  ! Read all relevant parameters and write them to the model log.
90  call log_version(param_file, mdl, version, &
91  "Parameters providing information about the vertical grid.")
92  call get_param(param_file, mdl, "G_EARTH", gv%mks_g_Earth, &
93  "The gravitational acceleration of the Earth.", &
94  units="m s-2", default = 9.80)
95  call get_param(param_file, mdl, "RHO_0", gv%Rho0, &
96  "The mean ocean density used with BOUSSINESQ true to "//&
97  "calculate accelerations and the mass for conservation "//&
98  "properties, or with BOUSSINSEQ false to convert some "//&
99  "parameters from vertical units of m to kg m-2.", &
100  units="kg m-3", default=1035.0)
101  call get_param(param_file, mdl, "BOUSSINESQ", gv%Boussinesq, &
102  "If true, make the Boussinesq approximation.", default=.true.)
103  call get_param(param_file, mdl, "ANGSTROM", gv%Angstrom_m, &
104  "The minimum layer thickness, usually one-Angstrom.", &
105  units="m", default=1.0e-10)
106  call get_param(param_file, mdl, "H_RESCALE_POWER", h_power, &
107  "An integer power of 2 that is used to rescale the model's "//&
108  "intenal units of thickness. Valid values range from -300 to 300.", &
109  units="nondim", default=0, debuggingparam=.true.)
110  if (abs(h_power) > 300) call mom_error(fatal, "verticalGridInit: "//&
111  "H_RESCALE_POWER is outside of the valid range of -300 to 300.")
112  h_rescale_factor = 1.0
113  if (h_power /= 0) h_rescale_factor = 2.0**h_power
114  if (.not.gv%Boussinesq) then
115  call get_param(param_file, mdl, "H_TO_KG_M2", gv%H_to_kg_m2,&
116  "A constant that translates thicknesses from the model's "//&
117  "internal units of thickness to kg m-2.", units="kg m-2 H-1", &
118  default=1.0)
119  gv%H_to_kg_m2 = gv%H_to_kg_m2 * h_rescale_factor
120  else
121  call get_param(param_file, mdl, "H_TO_M", gv%H_to_m, &
122  "A constant that translates the model's internal "//&
123  "units of thickness into m.", units="m H-1", default=1.0)
124  gv%H_to_m = gv%H_to_m * h_rescale_factor
125  endif
126  gv%g_Earth = us%m_to_L**2*us%Z_to_m*us%T_to_s**2 * gv%mks_g_Earth
127 #ifdef STATIC_MEMORY_
128  ! Here NK_ is a macro, while nk is a variable.
129  call get_param(param_file, mdl, "NK", nk, &
130  "The number of model layers.", units="nondim", &
131  static_value=nk_)
132  if (nk /= nk_) call mom_error(fatal, "verticalGridInit: " // &
133  "Mismatched number of layers NK_ between MOM_memory.h and param_file")
134 
135 #else
136  call get_param(param_file, mdl, "NK", nk, &
137  "The number of model layers.", units="nondim", fail_if_missing=.true.)
138 #endif
139  gv%ke = nk
140 
141  if (gv%Boussinesq) then
142  gv%H_to_kg_m2 = gv%Rho0 * gv%H_to_m
143  gv%kg_m2_to_H = 1.0 / gv%H_to_kg_m2
144  gv%m_to_H = 1.0 / gv%H_to_m
145  gv%Angstrom_H = gv%m_to_H * gv%Angstrom_m
146  else
147  gv%kg_m2_to_H = 1.0 / gv%H_to_kg_m2
148  gv%m_to_H = gv%Rho0 * gv%kg_m2_to_H
149  gv%H_to_m = gv%H_to_kg_m2 / gv%Rho0
150  gv%Angstrom_H = gv%Angstrom_m*1000.0*gv%kg_m2_to_H
151  endif
152  gv%H_subroundoff = 1e-20 * max(gv%Angstrom_H,gv%m_to_H*1e-17)
153  gv%H_to_Pa = gv%mks_g_Earth * gv%H_to_kg_m2
154 
155  gv%H_to_Z = gv%H_to_m * us%m_to_Z
156  gv%Z_to_H = us%Z_to_m * gv%m_to_H
157  gv%Angstrom_Z = us%m_to_Z * gv%Angstrom_m
158 
159 ! Log derivative values.
160  call log_param(param_file, mdl, "M to THICKNESS", gv%m_to_H*h_rescale_factor)
161  call log_param(param_file, mdl, "M to THICKNESS rescaled by 2^-n", gv%m_to_H)
162  call log_param(param_file, mdl, "THICKNESS to M rescaled by 2^n", gv%H_to_m)
163 
164  allocate( gv%sInterface(nk+1) )
165  allocate( gv%sLayer(nk) )
166  allocate( gv%g_prime(nk+1) ) ; gv%g_prime(:) = 0.0
167  ! The extent of Rlay should be changed to nk?
168  allocate( gv%Rlay(nk+1) ) ; gv%Rlay(:) = 0.0
169