MOM6
mom_verticalgrid Module Reference

Data Types

type  verticalgrid_type
 

Functions/Subroutines

subroutine, public verticalgridinit (param_file, GV)
 Allocates and initializes the model's vertical grid structure. 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

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

Parameters
[in]gvThe ocean's vertical grid structure

Definition at line 186 of file MOM_verticalGrid.F90.

Referenced by mom::register_diags().

186  character(len=48) :: get_flux_units
187  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
188 ! This subroutine returns the appropriate units for thickness fluxes,
189 ! depending on whether the model is Boussinesq or not and the scaling for
190 ! the vertical thickness.
191 
192 ! Arguments: G - The ocean's grid structure.
193 ! (ret) get_flux_units - The model's thickness flux units.
194 
195  if (gv%Boussinesq) then
196  get_flux_units = "meter3 second-1"
197  else
198  get_flux_units = "kilogram second-1"
199  endif
Here is the caller graph for this function:

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

Parameters
[in]gvThe ocean's vertical grid structure

Definition at line 168 of file MOM_verticalGrid.F90.

Referenced by mom::register_diags().

168  character(len=48) :: get_thickness_units
169  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
170 ! This subroutine returns the appropriate units for thicknesses,
171 ! depending on whether the model is Boussinesq or not and the scaling for
172 ! the vertical thickness.
173 
174 ! Arguments: G - The ocean's grid structure.
175 ! (ret) get_thickness_units - The model's vertical thickness units.
176 
177  if (gv%Boussinesq) then
178  get_thickness_units = "meter"
179  else
180  get_thickness_units = "kilogram meter-2"
181  endif
Here is the caller graph for this function:

◆ 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, forexample 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 204 of file MOM_verticalGrid.F90.

References mom_error_handler::mom_error().

Referenced by mom::register_diags().

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

◆ 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 271 of file MOM_verticalGrid.F90.

Referenced by mom_coord_initialization::mom_initialize_coord().

271  ! Arguments
272  type(verticalgrid_type), intent(inout) :: gv !< The container for vertical grid data
273  real, dimension(GV%ke), intent(in) :: rlay !< The layer target density
274  ! Local variables
275  integer :: k, nk
276 
277  nk = gv%ke
278 
279  gv%zAxisLongName = 'Target Potential Density'
280  gv%zAxisUnits = 'kg m-3'
281  do k=1,nk ; gv%sLayer(k) = rlay(k) ; enddo
282  if (nk > 1) then
283  gv%sInterface(1) = 1.5*rlay(1) - 0.5*rlay(2)
284  do k=2,nk ; gv%sInterface(k) = 0.5*( rlay(k-1) + rlay(k) ) ; enddo
285  gv%sInterface(nk+1) = 1.5*rlay(nk) - 0.5*rlay(nk-1)
286  else
287  gv%sInterface(1) = 0.0 ; gv%sInterface(nk+1) = 2.0*rlay(nk)
288  endif
289 
Here is the caller graph for this function:

◆ 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 294 of file MOM_verticalGrid.F90.

294 ! Arguments: G - The ocean's grid structure.
295  type(verticalgrid_type), pointer :: gv !< The ocean's vertical grid structure
296 
297  dealloc_(gv%g_prime) ; dealloc_(gv%Rlay)
298  dealloc_( gv%sInterface )
299  dealloc_( gv%sLayer )
300  deallocate( gv )
301 

◆ verticalgridinit()

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

Allocates and initializes the model's vertical grid structure.

Definition at line 87 of file MOM_verticalGrid.F90.

References mom_error_handler::mom_error().

87 ! This routine initializes the verticalGrid_type structure (GV).
88 ! All memory is allocated but not necessarily set to meaningful values until later.
89  type(param_file_type), intent(in) :: param_file ! Parameter file handle/type
90  type(verticalgrid_type), pointer :: gv ! The container for vertical grid data
91 ! This include declares and sets the variable "version".
92 #include "version_variable.h"
93  integer :: nk
94  character(len=16) :: mdl = 'MOM_verticalGrid'
95 
96  if (associated(gv)) call mom_error(fatal, &
97  'verticalGridInit: called with an associated GV pointer.')
98  allocate(gv)
99 
100  ! Read all relevant parameters and write them to the model log.
101  call log_version(param_file, mdl, version, &
102  "Parameters providing information about the vertical grid.")
103  call get_param(param_file, mdl, "G_EARTH", gv%g_Earth, &
104  "The gravitational acceleration of the Earth.", &
105  units="m s-2", default = 9.80)
106  call get_param(param_file, mdl, "RHO_0", gv%Rho0, &
107  "The mean ocean density used with BOUSSINESQ true to \n"//&
108  "calculate accelerations and the mass for conservation \n"//&
109  "properties, or with BOUSSINSEQ false to convert some \n"//&
110  "parameters from vertical units of m to kg m-2.", &
111  units="kg m-3", default=1035.0)
112  call get_param(param_file, mdl, "BOUSSINESQ", gv%Boussinesq, &
113  "If true, make the Boussinesq approximation.", default=.true.)
114  call get_param(param_file, mdl, "ANGSTROM", gv%Angstrom_z, &
115  "The minumum layer thickness, usually one-Angstrom.", &
116  units="m", default=1.0e-10)
117  if (.not.gv%Boussinesq) then
118  call get_param(param_file, mdl, "H_TO_KG_M2", gv%H_to_kg_m2,&
119  "A constant that translates thicknesses from the model's \n"//&
120  "internal units of thickness to kg m-2.", units="kg m-2 H-1", &
121  default=1.0)
122  else
123  call get_param(param_file, mdl, "H_TO_M", gv%H_to_m, &
124  "A constant that translates the model's internal \n"//&
125  "units of thickness into m.", units="m H-1", default=1.0)
126  endif
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 = gv%m_to_H * gv%Angstrom_z
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 = gv%Angstrom_z*1000.0*gv%kg_m2_to_H
151  endif
152  gv%H_subroundoff = 1e-20 * max(gv%Angstrom,gv%m_to_H*1e-17)
153  gv%H_to_Pa = gv%g_Earth * gv%H_to_kg_m2
154 
155 ! Log derivative values.
156  call log_param(param_file, mdl, "M to THICKNESS", gv%m_to_H)
157 
158  alloc_( gv%sInterface(nk+1) )
159  alloc_( gv%sLayer(nk) )
160  alloc_( gv%g_prime(nk+1) ) ; gv%g_prime(:) = 0.0
161  ! The extent of Rlay should be changed to nk?
162  alloc_( gv%Rlay(nk+1) ) ; gv%Rlay(:) = 0.0
163 
Here is the call graph for this function: