MOM6
mom_coord_initialization Module Reference

Detailed Description

Initializes fixed aspects of the related to its vertical coordinate.

Functions/Subroutines

subroutine, public mom_initialize_coord (GV, PF, write_geom, output_dir, tv, max_depth)
 MOM_initialize_coord sets up time-invariant quantities related to MOM6's vertical coordinate. More...
 
subroutine set_coord_from_gprime (Rlay, g_prime, GV, param_file)
 
subroutine set_coord_from_layer_density (Rlay, g_prime, GV, param_file)
 
subroutine set_coord_from_ts_ref (Rlay, g_prime, GV, param_file, eqn_of_state, P_Ref)
 
subroutine set_coord_from_ts_profile (Rlay, g_prime, GV, param_file, eqn_of_state, P_Ref)
 
subroutine set_coord_from_ts_range (Rlay, g_prime, GV, param_file, eqn_of_state, P_Ref)
 
subroutine set_coord_from_file (Rlay, g_prime, GV, param_file)
 
subroutine set_coord_linear (Rlay, g_prime, GV, param_file)
 
subroutine write_vertgrid_file (GV, param_file, directory)
 This subroutine writes out a file containing any available data related to the vertical grid used by the MOM ocean model. More...
 

Variables

character(len=40) mdl = "MOM_coord_initialization"
 

Function/Subroutine Documentation

◆ mom_initialize_coord()

subroutine, public mom_coord_initialization::mom_initialize_coord ( type(verticalgrid_type), intent(inout)  GV,
type(param_file_type), intent(in)  PF,
logical, intent(in)  write_geom,
character(len=*), intent(in)  output_dir,
type(thermo_var_ptrs), intent(inout)  tv,
real, intent(in)  max_depth 
)

MOM_initialize_coord sets up time-invariant quantities related to MOM6's vertical coordinate.

Parameters
[in,out]gvOcean vertical grid structure.
[in]pfA structure indicating the open file to parse for model parameter values.
[in]write_geomIf true, write grid geometry files.
[in]output_dirThe directory into which to write files.
[in,out]tvThe thermodynamic variable structure.
[in]max_depthThe ocean's maximum depth, in m.

Definition at line 35 of file MOM_coord_initialization.F90.

References bfb_initialization::bfb_set_coord(), mom_error_handler::calltree_enter(), mom_error_handler::calltree_leave(), mdl, set_coord_from_file(), set_coord_from_gprime(), set_coord_from_layer_density(), set_coord_from_ts_profile(), set_coord_from_ts_range(), set_coord_from_ts_ref(), set_coord_linear(), mom_verticalgrid::setverticalgridaxes(), user_initialization::user_set_coord(), and write_vertgrid_file().

Referenced by mom::initialize_mom().

35  type(verticalgrid_type), intent(inout) :: gv !< Ocean vertical grid structure.
36  type(param_file_type), intent(in) :: pf !< A structure indicating the open file
37  !! to parse for model parameter values.
38  logical, intent(in) :: write_geom !< If true, write grid geometry files.
39  character(len=*), intent(in) :: output_dir !< The directory into which to write files.
40  type(thermo_var_ptrs), intent(inout) :: tv !< The thermodynamic variable structure.
41  real, intent(in) :: max_depth !< The ocean's maximum depth, in m.
42  ! Local
43  character(len=200) :: config
44  logical :: debug
45 ! This include declares and sets the variable "version".
46 #include "version_variable.h"
47  integer :: nz
48  type(eos_type), pointer :: eos => null()
49 
50  if (associated(tv%eqn_of_state)) eos => tv%eqn_of_state
51 
52  nz = gv%ke
53 
54  call calltree_enter("MOM_initialize_coord(), MOM_coord_initialization.F90")
55  call log_version(pf, mdl, version, "")
56  call get_param(pf, mdl, "DEBUG", debug, default=.false.)
57 
58 ! Set-up the layer densities, GV%Rlay, and reduced gravities, GV%g_prime.
59  call get_param(pf, mdl, "COORD_CONFIG", config, &
60  "This specifies how layers are to be defined: \n"//&
61  " \t file - read coordinate information from the file \n"//&
62  " \t\t specified by (COORD_FILE).\n"//&
63  " \t BFB - Custom coords for buoyancy-forced basin case \n"//&
64  " \t\t based on SST_S, T_BOT and DRHO_DT.\n"//&
65  " \t linear - linear based on interfaces not layers \n"//&
66  " \t layer_ref - linear based on layer densities \n"//&
67  " \t ts_ref - use reference temperature and salinity \n"//&
68  " \t ts_range - use range of temperature and salinity \n"//&
69  " \t\t (T_REF and S_REF) to determine surface density \n"//&
70  " \t\t and GINT calculate internal densities. \n"//&
71  " \t gprime - use reference density (RHO_0) for surface \n"//&
72  " \t\t density and GINT calculate internal densities. \n"//&
73  " \t ts_profile - use temperature and salinity profiles \n"//&
74  " \t\t (read from COORD_FILE) to set layer densities. \n"//&
75  " \t USER - call a user modified routine.", &
76  fail_if_missing=.true.)
77  select case ( trim(config) )
78  case ("gprime")
79  call set_coord_from_gprime(gv%Rlay, gv%g_prime, gv, pf)
80  case ("layer_ref")
81  call set_coord_from_layer_density(gv%Rlay, gv%g_prime, gv, pf)
82  case ("linear")
83  call set_coord_linear(gv%Rlay, gv%g_prime, gv, pf)
84  case ("ts_ref")
85  call set_coord_from_ts_ref(gv%Rlay, gv%g_prime, gv, pf, eos, tv%P_Ref)
86  case ("ts_profile")
87  call set_coord_from_ts_profile(gv%Rlay, gv%g_prime, gv, pf, eos, tv%P_Ref)
88  case ("ts_range")
89  call set_coord_from_ts_range(gv%Rlay, gv%g_prime, gv, pf, eos, tv%P_Ref)
90  case ("file")
91  call set_coord_from_file(gv%Rlay, gv%g_prime, gv, pf)
92  case ("USER")
93  call user_set_coord(gv%Rlay, gv%g_prime, gv, pf, eos)
94  case ("BFB")
95  call bfb_set_coord(gv%Rlay, gv%g_prime, gv, pf, eos)
96  case ("none")
97  case default ; call mom_error(fatal,"MOM_initialize_coord: "// &
98  "Unrecognized coordinate setup"//trim(config))
99  end select
100  if (debug) call chksum(gv%Rlay, "MOM_initialize_coord: Rlay ", 1, nz)
101  if (debug) call chksum(gv%g_prime, "MOM_initialize_coord: g_prime ", 1, nz)
102  call setverticalgridaxes( gv%Rlay, gv )
103 
104 ! Copy the maximum depth across from the input argument
105  gv%max_depth = max_depth
106 
107 ! Write out all of the grid data used by this run.
108  if (write_geom) call write_vertgrid_file(gv, pf, output_dir)
109 
110  call calltree_leave('MOM_initialize_coord()')
111 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_coord_from_file()

subroutine mom_coord_initialization::set_coord_from_file ( real, dimension(:), intent(out)  Rlay,
real, dimension(:), intent(out)  g_prime,
type(verticalgrid_type), intent(in)  GV,
type(param_file_type), intent(in)  param_file 
)
private
Parameters
[out]rlayThe layers' target coordinate values (potential density).
[out]g_primeThe reduced gravity across the interfaces, in m s-2.
[in]gvThe ocean's vertical grid structure.
[in]param_fileA structure to parse for run-time parameters

Definition at line 420 of file MOM_coord_initialization.F90.

References mom_error_handler::calltree_enter(), and mom_error_handler::calltree_leave().

Referenced by mom_initialize_coord().

420  real, dimension(:), intent(out) :: rlay !< The layers' target coordinate values
421  !! (potential density).
422  real, dimension(:), intent(out) :: g_prime !< The reduced gravity across the interfaces,
423  !! in m s-2.
424  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
425  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
426 ! Arguments: Rlay - the layers' target coordinate values (potential density).
427 ! (out) g_prime - the reduced gravity across the interfaces, in m s-2.
428 ! (in) GV - The ocean's vertical grid structure.
429 ! (in) param_file - A structure indicating the open file to parse for
430 ! model parameter values.
431 
432 ! This subroutine sets the layer densities (Rlay) and the interface !
433 ! reduced gravities (g). !
434  real :: g_fs ! Reduced gravity across the free surface, in m s-2.
435  integer :: k, nz
436  character(len=40) :: mdl = "set_coord_from_file" ! This subroutine's name.
437  character(len=40) :: coord_var
438  character(len=200) :: filename,coord_file,inputdir ! Strings for file/path
439  nz = gv%ke
440 
441  call calltree_enter(trim(mdl)//"(), MOM_coord_initialization.F90")
442 
443  call get_param(param_file, mdl, "GFS", g_fs, &
444  "The reduced gravity at the free surface.", units="m s-2", &
445  default=gv%g_Earth)
446  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
447  inputdir = slasher(inputdir)
448  call get_param(param_file, mdl, "COORD_FILE", coord_file, &
449  "The file from which the coordinate densities are read.", &
450  fail_if_missing=.true.)
451  call get_param(param_file, mdl, "COORD_VAR", coord_var, &
452  "The variable in COORD_FILE that is to be used for the \n"//&
453  "coordinate densities.", default="Layer")
454  filename = trim(inputdir)//trim(coord_file)
455  call log_param(param_file, mdl, "INPUTDIR/COORD_FILE", filename)
456  if (.not.file_exists(filename)) call mom_error(fatal, &
457  " set_coord_from_file: Unable to open "//trim(filename))
458 
459  call read_axis_data(filename, coord_var, rlay)
460  g_prime(1) = g_fs
461  do k=2,nz ; g_prime(k) = (gv%g_Earth/gv%Rho0) * (rlay(k) - rlay(k-1)) ; enddo
462  do k=1,nz ; if (g_prime(k) <= 0.0) then
463  call mom_error(fatal, "MOM_initialization set_coord_from_file: "//&
464  "Zero or negative g_primes read from variable "//"Layer"//" in file "//&
465  trim(filename))
466  endif ; enddo
467 
468  call calltree_leave(trim(mdl)//'()')
Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_coord_from_gprime()

subroutine mom_coord_initialization::set_coord_from_gprime ( real, dimension(:), intent(out)  Rlay,
real, dimension(:), intent(out)  g_prime,
type(verticalgrid_type), intent(in)  GV,
type(param_file_type), intent(in)  param_file 
)
private
Parameters
[out]rlayThe layers' target coordinate values (potential density).
[out]g_primeA structure indicating the open file to parse for model parameter values.
[in]gvThe ocean's vertical grid structure.
[in]param_fileA structure to parse for run-time parameters

Definition at line 118 of file MOM_coord_initialization.F90.

References mom_error_handler::calltree_enter(), and mom_error_handler::calltree_leave().

Referenced by mom_initialize_coord().

118  real, dimension(:), intent(out) :: rlay !< The layers' target coordinate values
119  !! (potential density).
120  real, dimension(:), intent(out) :: g_prime !< A structure indicating the open file to
121  !! parse for model parameter values.
122  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
123  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
124 ! Arguments: Rlay - the layers' target coordinate values (potential density).
125 ! (out) g_prime - the reduced gravity across the interfaces, in m s-2.
126 ! (in) GV - The ocean's vertical grid structure.
127 ! (in) param_file - A structure indicating the open file to parse for
128 ! model parameter values.
129 
130 ! This subroutine sets the layer densities (Rlay) and the interface !
131 ! reduced gravities (g). !
132  real :: g_int ! Reduced gravities across the internal interfaces, in m s-2.
133  real :: g_fs ! Reduced gravity across the free surface, in m s-2.
134  character(len=40) :: mdl = "set_coord_from_gprime" ! This subroutine's name.
135  integer :: k, nz
136  nz = gv%ke
137 
138  call calltree_enter(trim(mdl)//"(), MOM_coord_initialization.F90")
139 
140  call get_param(param_file, mdl, "GFS" , g_fs, &
141  "The reduced gravity at the free surface.", units="m s-2", &
142  default=gv%g_Earth)
143  call get_param(param_file, mdl, "GINT", g_int, &
144  "The reduced gravity across internal interfaces.", &
145  units="m s-2", fail_if_missing=.true.)
146 
147  g_prime(1) = g_fs
148  do k=2,nz ; g_prime(k) = g_int ; enddo
149  rlay(1) = gv%Rho0
150  do k=2,nz ; rlay(k) = rlay(k-1) + g_prime(k)*(gv%Rho0/gv%g_Earth) ; enddo
151 
152  call calltree_leave(trim(mdl)//'()')
153 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_coord_from_layer_density()

subroutine mom_coord_initialization::set_coord_from_layer_density ( real, dimension(:), intent(out)  Rlay,
real, dimension(:), intent(out)  g_prime,
type(verticalgrid_type), intent(in)  GV,
type(param_file_type), intent(in)  param_file 
)
private
Parameters
[out]rlayThe layers' target coordinate values (potential density).
[out]g_primeThe reduced gravity across the interfaces, in m s-2.
[in]gvThe ocean's vertical grid structure
[in]param_fileA structure to parse for run-time parameters

Definition at line 159 of file MOM_coord_initialization.F90.

References mom_error_handler::calltree_enter(), and mom_error_handler::calltree_leave().

Referenced by mom_initialize_coord().

159  real, dimension(:), intent(out) :: rlay !< The layers' target coordinate values
160  !! (potential density).
161  real, dimension(:), intent(out) :: g_prime !< The reduced gravity across the interfaces,
162  !! in m s-2.
163  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
164  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
165 ! Arguments: Rlay - the layers' target coordinate values (potential density).
166 ! (out) g_prime - the reduced gravity across the interfaces, in m s-2.
167 ! (in) GV - The ocean's vertical grid structure.
168 ! (in) param_file - A structure indicating the open file to parse for
169 ! model parameter values.
170 
171 ! This subroutine sets the layer densities (Rlay) and the interface !
172 ! reduced gravities (g). !
173  real :: g_fs ! Reduced gravity across the free surface, in m s-2.
174  real :: rlay_ref! The surface layer's target density, in kg m-3.
175  real :: rlay_range ! The range of densities, in kg m-3.
176  character(len=40) :: mdl = "set_coord_from_layer_density" ! This subroutine's name.
177  integer :: k, nz
178  nz = gv%ke
179 
180  call calltree_enter(trim(mdl)//"(), MOM_coord_initialization.F90")
181 
182  call get_param(param_file, mdl, "GFS", g_fs, &
183  "The reduced gravity at the free surface.", units="m s-2", &
184  default=gv%g_Earth)
185  call get_param(param_file, mdl, "LIGHTEST_DENSITY", rlay_ref, &
186  "The reference potential density used for layer 1.", &
187  units="kg m-3", default=gv%Rho0)
188  call get_param(param_file, mdl, "DENSITY_RANGE", rlay_range, &
189  "The range of reference potential densities in the layers.", &
190  units="kg m-3", default=2.0)
191 
192  g_prime(1) = g_fs
193  rlay(1) = rlay_ref
194  do k=2,nz
195  rlay(k) = rlay(k-1) + rlay_range/(real(nz-1))
196  enddo
197 ! These statements set the interface reduced gravities. !
198  do k=2,nz
199  g_prime(k) = (gv%g_Earth/gv%Rho0) * (rlay(k) - rlay(k-1))
200  enddo
201 
202  call calltree_leave(trim(mdl)//'()')
Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_coord_from_ts_profile()

subroutine mom_coord_initialization::set_coord_from_ts_profile ( real, dimension(:), intent(out)  Rlay,
real, dimension(:), intent(out)  g_prime,
type(verticalgrid_type), intent(in)  GV,
type(param_file_type), intent(in)  param_file,
type(eos_type), pointer  eqn_of_state,
real, intent(in)  P_Ref 
)
private
Parameters
[out]rlayThe layers' target coordinate values (potential density).
[out]g_primeThe reduced gravity across the interfaces, in m s-2.
[in]gvThe ocean's vertical grid structure.
[in]param_fileA structure to parse for run-time parameters
eqn_of_stateinteger that selects equation of state.
[in]p_refThe coordinate-density reference pressure in Pa.

Definition at line 270 of file MOM_coord_initialization.F90.

References mom_error_handler::calltree_enter(), and mom_error_handler::calltree_leave().

Referenced by mom_initialize_coord().

270  real, dimension(:), intent(out) :: rlay !< The layers' target coordinate values
271  !! (potential density).
272  real, dimension(:), intent(out) :: g_prime !< The reduced gravity across the interfaces,
273  !! in m s-2.
274  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
275  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time
276  !! parameters
277  type(eos_type), pointer :: eqn_of_state !< integer that selects equation of state.
278  real, intent(in) :: p_ref !< The coordinate-density reference pressure
279  !! in Pa.
280 ! Arguments: Rlay - the layers' target coordinate values (potential density).
281 ! (out) g_prime - the reduced gravity across the interfaces, in m s-2.
282 ! (in) GV - The ocean's vertical grid structure.
283 ! (in) param_file - A structure indicating the open file to parse for
284 ! model parameter values.
285 ! (in) eqn_of_state - integer that selects equation of state
286 ! (in) P_Ref - The coordinate-density reference pressure in Pa.
287 
288 ! This subroutine sets the layer densities (Rlay) and the interface !
289 ! reduced gravities (g). !
290  real, dimension(GV%ke) :: t0, s0, pref
291  real :: g_fs ! Reduced gravity across the free surface, in m s-2.
292  integer :: k, nz
293  character(len=40) :: mdl = "set_coord_from_TS_profile" ! This subroutine's name.
294  character(len=200) :: filename, coord_file, inputdir ! Strings for file/path
295  nz = gv%ke
296 
297  call calltree_enter(trim(mdl)//"(), MOM_coord_initialization.F90")
298 
299  call get_param(param_file, mdl, "GFS", g_fs, &
300  "The reduced gravity at the free surface.", units="m s-2", &
301  default=gv%g_Earth)
302  call get_param(param_file, mdl, "COORD_FILE", coord_file, &
303  "The file from which the coordinate temperatures and \n"//&
304  "salnities are read.", fail_if_missing=.true.)
305 
306  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
307  filename = trim(slasher(inputdir))//trim(coord_file)
308  call log_param(param_file, mdl, "INPUTDIR/COORD_FILE", filename)
309 
310  call read_data(filename,"PTEMP",t0(:))
311  call read_data(filename,"SALT",s0(:))
312 
313  if (.not.file_exists(filename)) call mom_error(fatal, &
314  " set_coord_from_TS_profile: Unable to open " //trim(filename))
315 ! These statements set the interface reduced gravities. !
316  g_prime(1) = g_fs
317  do k=1,nz ; pref(k) = p_ref ; enddo
318  call calculate_density(t0, s0, pref, rlay, 1,nz,eqn_of_state)
319  do k=2,nz; g_prime(k) = (gv%g_Earth/gv%Rho0) * (rlay(k) - rlay(k-1)); enddo
320 
321  call calltree_leave(trim(mdl)//'()')
Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_coord_from_ts_range()

subroutine mom_coord_initialization::set_coord_from_ts_range ( real, dimension(:), intent(out)  Rlay,
real, dimension(:), intent(out)  g_prime,
type(verticalgrid_type), intent(in)  GV,
type(param_file_type), intent(in)  param_file,
type(eos_type), pointer  eqn_of_state,
real, intent(in)  P_Ref 
)
private
Parameters
[out]rlayThe layers' target coordinate values (potential density).
[out]g_primeThe reduced gravity across the interfaces, in m s-2.
[in]gvThe ocean's vertical grid structure.
[in]param_fileA structure to parse for run-time parameters
eqn_of_stateinteger that selects equation of state
[in]p_refThe coordinate-density reference pressure in Pa.

Definition at line 328 of file MOM_coord_initialization.F90.

References mom_error_handler::calltree_enter(), and mom_error_handler::calltree_leave().

Referenced by mom_initialize_coord().

328  real, dimension(:), intent(out) :: rlay !< The layers' target coordinate values
329  !! (potential density).
330  real, dimension(:), intent(out) :: g_prime !< The reduced gravity across the interfaces,
331  !! in m s-2.
332  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
333  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time
334  !! parameters
335  type(eos_type), pointer :: eqn_of_state !< integer that selects equation of state
336  real, intent(in) :: p_ref !< The coordinate-density reference pressure
337  !! in Pa.
338 ! Arguments: Rlay - the layers' target coordinate values (potential density).
339 ! (out) g_prime - the reduced gravity across the interfaces, in m s-2.
340 ! (in) GV - The ocean's vertical grid structure.
341 ! (in) param_file - A structure indicating the open file to parse for
342 ! model parameter values.
343 ! (in) eqn_of_state - integer that selects equation of state
344 ! (in) P_Ref - The coordinate-density reference pressure in Pa.
345 
346 ! This subroutine sets the layer densities (Rlay) and the interface !
347 ! reduced gravities (g). !
348  real, dimension(GV%ke) :: t0, s0, pref
349  real :: s_ref, s_light, s_dense ! Salnity range parameters in PSU.
350  real :: t_ref, t_light, t_dense ! Temperature range parameters in dec C.
351  real :: res_rat ! The ratio of density space resolution in the denser part
352  ! of the range to that in the lighter part of the range.
353  ! Setting this greater than 1 increases the resolution for
354  ! the denser water.
355  real :: g_fs ! Reduced gravity across the free surface, in m s-2.
356  real :: a1, frac_dense, k_frac
357  integer :: k, nz, k_light
358  character(len=40) :: mdl = "set_coord_from_TS_range" ! This subroutine's name.
359  character(len=200) :: filename, coord_file, inputdir ! Strings for file/path
360  nz = gv%ke
361 
362  call calltree_enter(trim(mdl)//"(), MOM_coord_initialization.F90")
363 
364  call get_param(param_file, mdl, "T_REF", t_ref, &
365  "The default initial temperatures.", units="degC", default=10.0)
366  call get_param(param_file, mdl, "TS_RANGE_T_LIGHT", t_light, &
367  "The initial temperature of the lightest layer when \n"//&
368  "COORD_CONFIG is set to ts_range.", units="degC", default=t_ref)
369  call get_param(param_file, mdl, "TS_RANGE_T_DENSE", t_dense, &
370  "The initial temperature of the densest layer when \n"//&
371  "COORD_CONFIG is set to ts_range.", units="degC", default=t_ref)
372 
373  call get_param(param_file, mdl, "S_REF", s_ref, &
374  "The default initial salinities.", units="PSU", default=35.0)
375  call get_param(param_file, mdl, "TS_RANGE_S_LIGHT", s_light, &
376  "The initial lightest salinities when COORD_CONFIG \n"//&
377  "is set to ts_range.", default = s_ref, units="PSU")
378  call get_param(param_file, mdl, "TS_RANGE_S_DENSE", s_dense, &
379  "The initial densest salinities when COORD_CONFIG \n"//&
380  "is set to ts_range.", default = s_ref, units="PSU")
381 
382  call get_param(param_file, mdl, "TS_RANGE_RESOLN_RATIO", res_rat, &
383  "The ratio of density space resolution in the densest \n"//&
384  "part of the range to that in the lightest part of the \n"//&
385  "range when COORD_CONFIG is set to ts_range. Values \n"//&
386  "greater than 1 increase the resolution of the denser water.",&
387  default=1.0, units="nondim")
388 
389  call get_param(param_file, mdl, "GFS", g_fs, &
390  "The reduced gravity at the free surface.", units="m s-2", &
391  default=gv%g_Earth)
392 
393  k_light = gv%nk_rho_varies + 1
394 
395  ! Set T0(k) to range from T_LIGHT to T_DENSE, and simliarly for S0(k).
396  t0(k_light) = t_light ; s0(k_light) = s_light
397  a1 = 2.0 * res_rat / (1.0 + res_rat)
398  do k=k_light+1,nz
399  k_frac = real(k-k_light)/real(nz-k_light)
400  frac_dense = a1 * k_frac + (1.0 - a1) * k_frac**2
401  t0(k) = frac_dense * (t_dense - t_light) + t_light
402  s0(k) = frac_dense * (s_dense - s_light) + s_light
403  enddo
404 
405  g_prime(1) = g_fs
406  do k=1,nz ; pref(k) = p_ref ; enddo
407  call calculate_density(t0, s0, pref, rlay, k_light,nz-k_light+1,eqn_of_state)
408  ! Extrapolate target densities for the variable density mixed and buffer layers.
409  do k=k_light-1,1,-1
410  rlay(k) = 2.0*rlay(k+1) - rlay(k+2)
411  enddo
412  do k=2,nz; g_prime(k) = (gv%g_Earth/gv%Rho0) * (rlay(k) - rlay(k-1)); enddo
413 
414  call calltree_leave(trim(mdl)//'()')
Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_coord_from_ts_ref()

subroutine mom_coord_initialization::set_coord_from_ts_ref ( real, dimension(:), intent(out)  Rlay,
real, dimension(:), intent(out)  g_prime,
type(verticalgrid_type), intent(in)  GV,
type(param_file_type), intent(in)  param_file,
type(eos_type), pointer  eqn_of_state,
real, intent(in)  P_Ref 
)
private
Parameters
[out]rlayThe layers' target coordinate values (potential density).
[out]g_primethe reduced gravity across the interfaces, in m s-2.
[in]gvThe ocean's vertical grid structure.
[in]param_fileA structure to parse for run-time parameters
eqn_of_stateinteger selecting the equation of state.
[in]p_refThe coordinate-density reference pressure in Pa.

Definition at line 209 of file MOM_coord_initialization.F90.

References mom_error_handler::calltree_enter(), and mom_error_handler::calltree_leave().

Referenced by mom_initialize_coord().

209  real, dimension(:), intent(out) :: rlay !< The layers' target coordinate values
210  !! (potential density).
211  real, dimension(:), intent(out) :: g_prime !< the reduced gravity across the interfaces,
212  !! in m s-2.
213  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
214  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time
215  !! parameters
216  type(eos_type), pointer :: eqn_of_state !< integer selecting the equation of state.
217  real, intent(in) :: p_ref !< The coordinate-density reference pressure
218  !! in Pa.
219 ! Arguments: Rlay - the layers' target coordinate values (potential density).
220 ! (out) g_prime - the reduced gravity across the interfaces, in m s-2.
221 ! (in) GV - The ocean's vertical grid structure.
222 ! (in) param_file - A structure indicating the open file to parse for
223 ! model parameter values.
224 ! (in) eqn_of_state - integer selecting the equation of state
225 ! (in) P_Ref - The coordinate-density reference pressure in Pa.
226 
227 ! This subroutine sets the layer densities (Rlay) and the interface !
228 ! reduced gravities (g). !
229  real :: t_ref ! Reference temperature
230  real :: s_ref ! Reference salinity
231  real :: g_int ! Reduced gravities across the internal interfaces, in m s-2.
232  real :: g_fs ! Reduced gravity across the free surface, in m s-2.
233  character(len=40) :: mdl = "set_coord_from_TS_ref" ! This subroutine's name.
234  integer :: k, nz
235  nz = gv%ke
236 
237  call calltree_enter(trim(mdl)//"(), MOM_coord_initialization.F90")
238 
239  call get_param(param_file, mdl, "T_REF", t_ref, &
240  "The initial temperature of the lightest layer.", units="degC", &
241  fail_if_missing=.true.)
242  call get_param(param_file, mdl, "S_REF", s_ref, &
243  "The initial salinities.", units="PSU", default=35.0)
244  call get_param(param_file, mdl, "GFS", g_fs, &
245  "The reduced gravity at the free surface.", units="m s-2", &
246  default=gv%g_Earth)
247  call get_param(param_file, mdl, "GINT", g_int, &
248  "The reduced gravity across internal interfaces.", &
249  units="m s-2", fail_if_missing=.true.)
250  !
251 ! These statements set the interface reduced gravities. !
252  g_prime(1) = g_fs
253  do k=2,nz ; g_prime(k) = g_int ; enddo
254 
255 ! The uppermost layer's density is set here. Subsequent layers' !
256 ! densities are determined from this value and the g values. !
257 ! T0 = 28.228 ; S0 = 34.5848 ; Pref = P_Ref
258  call calculate_density(t_ref, s_ref, p_ref, rlay(1), eqn_of_state)
259 
260 ! These statements set the layer densities. !
261  do k=2,nz ; rlay(k) = rlay(k-1) + g_prime(k)*(gv%Rho0/gv%g_Earth) ; enddo
262 
263  call calltree_leave(trim(mdl)//'()')
Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_coord_linear()

subroutine mom_coord_initialization::set_coord_linear ( real, dimension(:), intent(out)  Rlay,
real, dimension(:), intent(out)  g_prime,
type(verticalgrid_type), intent(in)  GV,
type(param_file_type), intent(in)  param_file 
)
private
Parameters
[out]rlayThe layers' target coordinate values (potential density).
[out]g_primeThe reduced gravity across the interfaces, in m s-2.
[in]gvThe ocean's vertical grid structure.
[in]param_fileA structure to parse for run-time parameters

Definition at line 474 of file MOM_coord_initialization.F90.

References mom_error_handler::calltree_enter(), and mom_error_handler::calltree_leave().

Referenced by mom_initialize_coord().

474  real, dimension(:), intent(out) :: rlay !< The layers' target coordinate values
475  !! (potential density).
476  real, dimension(:), intent(out) :: g_prime !< The reduced gravity across the interfaces,
477  !! in m s-2.
478  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
479  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
480 ! Arguments: Rlay - the layers' target coordinate values (potential density).
481 ! (out) g_prime - the reduced gravity across the interfaces, in m s-2.
482 ! (in) GV - The ocean's vertical grid structure.
483 ! (in) param_file - A structure indicating the open file to parse for
484 ! model parameter values.
485 
486 ! This subroutine sets the layer densities (Rlay) and the interface
487 ! reduced gravities (g) according to a linear profile starting at a
488 ! reference surface layer density and spanning a range of densities
489 ! to the bottom defined by the parameter RLAY_RANGE
490 ! (defaulting to 2.0 if not defined)
491  character(len=40) :: mdl = "set_coord_linear" ! This subroutine
492  real :: rlay_ref, rlay_range, g_fs
493  integer :: k, nz
494  nz = gv%ke
495 
496  call calltree_enter(trim(mdl)//"(), MOM_coord_initialization.F90")
497 
498  call get_param(param_file, mdl, "LIGHTEST_DENSITY", rlay_ref, &
499  "The reference potential density used for the surface \n"// &
500  "interface.", units="kg m-3", default=gv%Rho0)
501  call get_param(param_file, mdl, "DENSITY_RANGE", rlay_range, &
502  "The range of reference potential densities across \n"// &
503  "all interfaces.", units="kg m-3", default=2.0)
504  call get_param(param_file, mdl, "GFS", g_fs, &
505  "The reduced gravity at the free surface.", units="m s-2", &
506  default=gv%g_Earth)
507 
508  ! This following sets the target layer densities such that a the
509  ! surface interface has density Rlay_ref and the bottom
510  ! is Rlay_range larger
511  do k=1,nz
512  rlay(k) = rlay_ref + rlay_range*((real(k)-0.5)/real(nz))
513  enddo
514  ! These statements set the interface reduced gravities.
515  g_prime(1) = g_fs
516  do k=2,nz
517  g_prime(k) = (gv%g_Earth/gv%Rho0) * (rlay(k) - rlay(k-1))
518  enddo
519 
520  call calltree_leave(trim(mdl)//'()')
Here is the call graph for this function:
Here is the caller graph for this function:

◆ write_vertgrid_file()

subroutine mom_coord_initialization::write_vertgrid_file ( type(verticalgrid_type), intent(in)  GV,
type(param_file_type), intent(in)  param_file,
character(len=*), intent(in)  directory 
)
private

This subroutine writes out a file containing any available data related to the vertical grid used by the MOM ocean model.

Parameters
[in]gvThe ocean's vertical grid structure
[in]param_fileA structure to parse for run-time parameters
[in]directoryThe directory into which to place the file.

Definition at line 528 of file MOM_coord_initialization.F90.

References mom_io::var_desc().

Referenced by mom_initialize_coord().

528  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
529  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
530  character(len=*), intent(in) :: directory !< The directory into which to place the file.
531 ! This subroutine writes out a file containing any available data related
532 ! to the vertical grid used by the MOM ocean model.
533 ! Arguments: GV - The container for the vertical grid data.
534 ! (in) param_file - A structure indicating the open file to parse for
535 ! model parameter values.
536 ! (in) directory - The directory into which to place the file.
537  character(len=240) :: filepath
538  type(vardesc) :: vars(2)
539  type(fieldtype) :: fields(2)
540  integer :: unit
541 
542  filepath = trim(directory) // trim("Vertical_coordinate")
543 
544  vars(1) = var_desc("R","kilogram meter-3","Target Potential Density",'1','L','1')
545  vars(2) = var_desc("g","meter second-2","Reduced gravity",'1','L','1')
546 
547  call create_file(unit, trim(filepath), vars, 2, fields, single_file, gv=gv)
548 
549  call write_field(unit, fields(1), gv%Rlay)
550  call write_field(unit, fields(2), gv%g_prime)
551 
552  call close_file(unit)
553 
Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ mdl

character(len=40) mom_coord_initialization::mdl = "MOM_coord_initialization"
private

Definition at line 27 of file MOM_coord_initialization.F90.

Referenced by mom_initialize_coord().

27 character(len=40) :: mdl = "MOM_coord_initialization" ! This module's name.