MOM6
MOM_coord_initialization.F90
Go to the documentation of this file.
1 !> Initializes fixed aspects of the related to its vertical coordinate.
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 use mom_debugging, only : chksum
8 use mom_error_handler, only : mom_mesg, mom_error, fatal, warning, is_root_pe
11 use mom_file_parser, only : log_version
12 use mom_io, only : close_file, create_file, fieldtype, file_exists
13 use mom_io, only : open_file, mom_read_data, read_axis_data, single_file, multiple
14 use mom_io, only : slasher, vardesc, write_field, var_desc
21 
22 use netcdf
23 
24 implicit none ; private
25 
27 
28 ! A note on unit descriptions in comments: MOM6 uses units that can be rescaled for dimensional
29 ! consistency testing. These are noted in comments with units like Z, H, L, and T, along with
30 ! their mks counterparts with notation like "a velocity [Z T-1 ~> m s-1]". If the units
31 ! vary with the Boussinesq approximation, the Boussinesq variant is given first.
32 
33 character(len=40) :: mdl = "MOM_coord_initialization" !< This module's name.
34 
35 contains
36 
37 !> MOM_initialize_coord sets up time-invariant quantities related to MOM6's
38 !! vertical coordinate.
39 subroutine mom_initialize_coord(GV, US, PF, write_geom, output_dir, tv, max_depth)
40  type(verticalgrid_type), intent(inout) :: gv !< Ocean vertical grid structure.
41  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
42  type(param_file_type), intent(in) :: pf !< A structure indicating the open file
43  !! to parse for model parameter values.
44  logical, intent(in) :: write_geom !< If true, write grid geometry files.
45  character(len=*), intent(in) :: output_dir !< The directory into which to write files.
46  type(thermo_var_ptrs), intent(inout) :: tv !< The thermodynamic variable structure.
47  real, intent(in) :: max_depth !< The ocean's maximum depth [Z ~> m].
48  ! Local
49  character(len=200) :: config
50  logical :: debug
51 ! This include declares and sets the variable "version".
52 #include "version_variable.h"
53  integer :: nz
54  type(eos_type), pointer :: eos => null()
55 
56  if (associated(tv%eqn_of_state)) eos => tv%eqn_of_state
57 
58  nz = gv%ke
59 
60  call calltree_enter("MOM_initialize_coord(), MOM_coord_initialization.F90")
61  call log_version(pf, mdl, version, "")
62  call get_param(pf, mdl, "DEBUG", debug, default=.false.)
63 
64 ! Set-up the layer densities, GV%Rlay, and reduced gravities, GV%g_prime.
65  call get_param(pf, mdl, "COORD_CONFIG", config, &
66  "This specifies how layers are to be defined: \n"//&
67  " \t ALE or none - used to avoid defining layers in ALE mode \n"//&
68  " \t file - read coordinate information from the file \n"//&
69  " \t\t specified by (COORD_FILE).\n"//&
70  " \t BFB - Custom coords for buoyancy-forced basin case \n"//&
71  " \t\t based on SST_S, T_BOT and DRHO_DT.\n"//&
72  " \t linear - linear based on interfaces not layers \n"//&
73  " \t layer_ref - linear based on layer densities \n"//&
74  " \t ts_ref - use reference temperature and salinity \n"//&
75  " \t ts_range - use range of temperature and salinity \n"//&
76  " \t\t (T_REF and S_REF) to determine surface density \n"//&
77  " \t\t and GINT calculate internal densities. \n"//&
78  " \t gprime - use reference density (RHO_0) for surface \n"//&
79  " \t\t density and GINT calculate internal densities. \n"//&
80  " \t ts_profile - use temperature and salinity profiles \n"//&
81  " \t\t (read from COORD_FILE) to set layer densities. \n"//&
82  " \t USER - call a user modified routine.", &
83  fail_if_missing=.true.)
84  select case ( trim(config) )
85  case ("gprime")
86  call set_coord_from_gprime(gv%Rlay, gv%g_prime, gv, us, pf)
87  case ("layer_ref")
88  call set_coord_from_layer_density(gv%Rlay, gv%g_prime, gv, us, pf)
89  case ("linear")
90  call set_coord_linear(gv%Rlay, gv%g_prime, gv, us, pf)
91  case ("ts_ref")
92  call set_coord_from_ts_ref(gv%Rlay, gv%g_prime, gv, us, pf, eos, tv%P_Ref)
93  case ("ts_profile")
94  call set_coord_from_ts_profile(gv%Rlay, gv%g_prime, gv, us, pf, eos, tv%P_Ref)
95  case ("ts_range")
96  call set_coord_from_ts_range(gv%Rlay, gv%g_prime, gv, us, pf, eos, tv%P_Ref)
97  case ("file")
98  call set_coord_from_file(gv%Rlay, gv%g_prime, gv, us, pf)
99  case ("USER")
100  call user_set_coord(gv%Rlay, gv%g_prime, gv, us, pf, eos)
101  case ("BFB")
102  call bfb_set_coord(gv%Rlay, gv%g_prime, gv, us, pf, eos)
103  case ("none", "ALE")
104  call set_coord_to_none(gv%Rlay, gv%g_prime, gv, us, pf)
105  case default ; call mom_error(fatal,"MOM_initialize_coord: "// &
106  "Unrecognized coordinate setup"//trim(config))
107  end select
108  if (debug) call chksum(us%R_to_kg_m3*gv%Rlay(:), "MOM_initialize_coord: Rlay ", 1, nz)
109  if (debug) call chksum(us%m_to_Z*us%L_to_m**2*us%s_to_T**2*gv%g_prime(:), "MOM_initialize_coord: g_prime ", 1, nz)
110  call setverticalgridaxes( gv%Rlay, gv, scale=us%R_to_kg_m3 )
111 
112 ! Copy the maximum depth across from the input argument
113  gv%max_depth = max_depth
114 
115 ! Write out all of the grid data used by this run.
116  if (write_geom) call write_vertgrid_file(gv, us, pf, output_dir)
117 
118  call calltree_leave('MOM_initialize_coord()')
119 
120 end subroutine mom_initialize_coord
121 
122 ! The set_coord routines deal with initializing aspects of the vertical grid.
123 
124 !> Sets the layer densities (Rlay) and the interface reduced gravities (g).
125 subroutine set_coord_from_gprime(Rlay, g_prime, GV, US, param_file)
126  real, dimension(:), intent(out) :: Rlay !< The layers' target coordinate values
127  !! (potential density) [R ~> kg m-3].
128  real, dimension(:), intent(out) :: g_prime !< The reduced gravity across the interfaces
129  !! [L2 Z-1 T-2 ~> m s-2].
130  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
131  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
132  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
133  ! Local variables
134  real :: g_int ! Reduced gravities across the internal interfaces [m s-2].
135  real :: g_fs ! Reduced gravity across the free surface [m s-2].
136  character(len=40) :: mdl = "set_coord_from_gprime" ! This subroutine's name.
137  integer :: k, nz
138  nz = gv%ke
139 
140  call calltree_enter(trim(mdl)//"(), MOM_coord_initialization.F90")
141 
142  call get_param(param_file, mdl, "GFS" , g_fs, &
143  "The reduced gravity at the free surface.", units="m s-2", &
144  default=gv%mks_g_Earth, scale=us%m_s_to_L_T**2*us%Z_to_m)
145  call get_param(param_file, mdl, "GINT", g_int, &
146  "The reduced gravity across internal interfaces.", &
147  units="m s-2", fail_if_missing=.true., scale=us%m_s_to_L_T**2*us%Z_to_m)
148 
149  g_prime(1) = g_fs
150  do k=2,nz ; g_prime(k) = g_int ; enddo
151  rlay(1) = gv%Rho0
152  do k=2,nz ; rlay(k) = rlay(k-1) + g_prime(k)*(gv%Rho0/gv%g_Earth) ; enddo
153 
154  call calltree_leave(trim(mdl)//'()')
155 
156 end subroutine set_coord_from_gprime
157 
158 !> Sets the layer densities (Rlay) and the interface reduced gravities (g).
159 subroutine set_coord_from_layer_density(Rlay, g_prime, GV, US, param_file)
160  real, dimension(:), intent(out) :: Rlay !< The layers' target coordinate values
161  !! (potential density) [R ~> kg m-3].
162  real, dimension(:), intent(out) :: g_prime !< The reduced gravity across the interfaces
163  !! [L2 Z-1 T-2 ~> m s-2].
164  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
165  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
166  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
167  ! Local variables
168  real :: g_fs ! Reduced gravity across the free surface [m s-2].
169  real :: Rlay_Ref! The surface layer's target density [kg m-3].
170  real :: RLay_range ! The range of densities [kg m-3].
171  character(len=40) :: mdl = "set_coord_from_layer_density" ! This subroutine's name.
172  integer :: k, nz
173  nz = gv%ke
174 
175  call calltree_enter(trim(mdl)//"(), MOM_coord_initialization.F90")
176 
177  call get_param(param_file, mdl, "GFS", g_fs, &
178  "The reduced gravity at the free surface.", units="m s-2", &
179  default=gv%mks_g_Earth, scale=us%m_s_to_L_T**2*us%Z_to_m)
180  call get_param(param_file, mdl, "LIGHTEST_DENSITY", rlay_ref, &
181  "The reference potential density used for layer 1.", &
182  units="kg m-3", default=us%R_to_kg_m3*gv%Rho0, scale=us%kg_m3_to_R)
183  call get_param(param_file, mdl, "DENSITY_RANGE", rlay_range, &
184  "The range of reference potential densities in the layers.", &
185  units="kg m-3", default=2.0, scale=us%kg_m3_to_R)
186 
187  g_prime(1) = g_fs
188  rlay(1) = rlay_ref
189  do k=2,nz
190  rlay(k) = rlay(k-1) + rlay_range/(real(nz-1))
191  enddo
192 ! These statements set the interface reduced gravities. !
193  do k=2,nz
194  g_prime(k) = (gv%g_Earth/(gv%Rho0)) * (rlay(k) - rlay(k-1))
195  enddo
196 
197  call calltree_leave(trim(mdl)//'()')
198 end subroutine set_coord_from_layer_density
199 
200 !> Sets the layer densities (Rlay) and the interface reduced gravities (g) from a profile of g'.
201 subroutine set_coord_from_ts_ref(Rlay, g_prime, GV, US, param_file, eqn_of_state, &
202  P_Ref)
203  real, dimension(:), intent(out) :: Rlay !< The layers' target coordinate values
204  !! (potential density) [R ~> kg m-3].
205  real, dimension(:), intent(out) :: g_prime !< The reduced gravity across the interfaces
206  !! [L2 Z-1 T-2 ~> m s-2].
207  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
208  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
209  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time
210  !! parameters
211  type(eos_type), pointer :: eqn_of_state !< integer selecting the equation of state.
212  real, intent(in) :: P_Ref !< The coordinate-density reference pressure [Pa].
213  ! Local variables
214  real :: T_ref ! Reference temperature
215  real :: S_ref ! Reference salinity
216  real :: g_int ! Reduced gravities across the internal interfaces [m s-2].
217  real :: g_fs ! Reduced gravity across the free surface [m s-2].
218  character(len=40) :: mdl = "set_coord_from_TS_ref" ! This subroutine's name.
219  integer :: k, nz
220  nz = gv%ke
221 
222  call calltree_enter(trim(mdl)//"(), MOM_coord_initialization.F90")
223 
224  call get_param(param_file, mdl, "T_REF", t_ref, &
225  "The initial temperature of the lightest layer.", units="degC", &
226  fail_if_missing=.true.)
227  call get_param(param_file, mdl, "S_REF", s_ref, &
228  "The initial salinities.", units="PSU", default=35.0)
229  call get_param(param_file, mdl, "GFS", g_fs, &
230  "The reduced gravity at the free surface.", units="m s-2", &
231  default=gv%mks_g_Earth, scale=us%m_s_to_L_T**2*us%Z_to_m)
232  call get_param(param_file, mdl, "GINT", g_int, &
233  "The reduced gravity across internal interfaces.", &
234  units="m s-2", fail_if_missing=.true., scale=us%m_s_to_L_T**2*us%Z_to_m)
235  !
236 ! These statements set the interface reduced gravities. !
237  g_prime(1) = g_fs
238  do k=2,nz ; g_prime(k) = g_int ; enddo
239 
240 ! The uppermost layer's density is set here. Subsequent layers' !
241 ! densities are determined from this value and the g values. !
242 ! T0 = 28.228 ; S0 = 34.5848 ; Pref = P_Ref
243  call calculate_density(t_ref, s_ref, p_ref, rlay(1), eqn_of_state, scale=us%kg_m3_to_R)
244 
245 ! These statements set the layer densities. !
246  do k=2,nz ; rlay(k) = rlay(k-1) + g_prime(k)*(gv%Rho0/gv%g_Earth) ; enddo
247 
248  call calltree_leave(trim(mdl)//'()')
249 end subroutine set_coord_from_ts_ref
250 
251 !> Sets the layer densities (Rlay) and the interface reduced gravities (g) from a T-S profile.
252 subroutine set_coord_from_ts_profile(Rlay, g_prime, GV, US, param_file, &
253  eqn_of_state, P_Ref)
254  real, dimension(:), intent(out) :: Rlay !< The layers' target coordinate values
255  !! (potential density) [R ~> kg m-3].
256  real, dimension(:), intent(out) :: g_prime !< The reduced gravity across the interfaces
257  !! [L2 Z-1 T-2 ~> m s-2].
258  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
259  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
260  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time
261  !! parameters
262  type(eos_type), pointer :: eqn_of_state !< integer that selects equation of state.
263  real, intent(in) :: P_Ref !< The coordinate-density reference pressure [Pa].
264  ! Local variables
265  real, dimension(GV%ke) :: T0, S0, Pref
266  real :: g_fs ! Reduced gravity across the free surface [m s-2].
267  integer :: k, nz
268  character(len=40) :: mdl = "set_coord_from_TS_profile" ! This subroutine's name.
269  character(len=200) :: filename, coord_file, inputdir ! Strings for file/path
270  nz = gv%ke
271 
272  call calltree_enter(trim(mdl)//"(), MOM_coord_initialization.F90")
273 
274  call get_param(param_file, mdl, "GFS", g_fs, &
275  "The reduced gravity at the free surface.", units="m s-2", &
276  default=gv%mks_g_Earth, scale=us%m_s_to_L_T**2*us%Z_to_m)
277  call get_param(param_file, mdl, "COORD_FILE", coord_file, &
278  "The file from which the coordinate temperatures and "//&
279  "salinities are read.", fail_if_missing=.true.)
280 
281  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
282  filename = trim(slasher(inputdir))//trim(coord_file)
283  call log_param(param_file, mdl, "INPUTDIR/COORD_FILE", filename)
284 
285  call mom_read_data(filename,"PTEMP",t0(:))
286  call mom_read_data(filename,"SALT",s0(:))
287 
288  if (.not.file_exists(filename)) call mom_error(fatal, &
289  " set_coord_from_TS_profile: Unable to open " //trim(filename))
290 ! These statements set the interface reduced gravities. !
291  g_prime(1) = g_fs
292  do k=1,nz ; pref(k) = p_ref ; enddo
293  call calculate_density(t0, s0, pref, rlay, 1, nz, eqn_of_state, scale=us%kg_m3_to_R)
294  do k=2,nz; g_prime(k) = (gv%g_Earth/(gv%Rho0)) * (rlay(k) - rlay(k-1)) ; enddo
295 
296  call calltree_leave(trim(mdl)//'()')
297 end subroutine set_coord_from_ts_profile
298 
299 !> Sets the layer densities (Rlay) and the interface reduced gravities (g) from a linear T-S profile.
300 subroutine set_coord_from_ts_range(Rlay, g_prime, GV, US, param_file, &
301  eqn_of_state, P_Ref)
302  real, dimension(:), intent(out) :: Rlay !< The layers' target coordinate values
303  !! (potential density) [R ~> kg m-3].
304  real, dimension(:), intent(out) :: g_prime !< The reduced gravity across the interfaces
305  !! [L2 Z-1 T-2 ~> m s-2].
306  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
307  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
308  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time
309  !! parameters
310  type(eos_type), pointer :: eqn_of_state !< integer that selects equation of state
311  real, intent(in) :: P_Ref !< The coordinate-density reference pressure [Pa]
312 
313  ! Local variables
314  real, dimension(GV%ke) :: T0, S0, Pref
315  real :: S_Ref, S_Light, S_Dense ! Salinity range parameters [ppt].
316  real :: T_Ref, T_Light, T_Dense ! Temperature range parameters [decC].
317  real :: res_rat ! The ratio of density space resolution in the denser part
318  ! of the range to that in the lighter part of the range.
319  ! Setting this greater than 1 increases the resolution for
320  ! the denser water.
321  real :: g_fs ! Reduced gravity across the free surface [m s-2].
322  real :: a1, frac_dense, k_frac
323  integer :: k, nz, k_light
324  character(len=40) :: mdl = "set_coord_from_TS_range" ! This subroutine's name.
325  character(len=200) :: filename, coord_file, inputdir ! Strings for file/path
326  nz = gv%ke
327 
328  call calltree_enter(trim(mdl)//"(), MOM_coord_initialization.F90")
329 
330  call get_param(param_file, mdl, "T_REF", t_ref, &
331  "The default initial temperatures.", units="degC", default=10.0)
332  call get_param(param_file, mdl, "TS_RANGE_T_LIGHT", t_light, &
333  "The initial temperature of the lightest layer when "//&
334  "COORD_CONFIG is set to ts_range.", units="degC", default=t_ref)
335  call get_param(param_file, mdl, "TS_RANGE_T_DENSE", t_dense, &
336  "The initial temperature of the densest layer when "//&
337  "COORD_CONFIG is set to ts_range.", units="degC", default=t_ref)
338 
339  call get_param(param_file, mdl, "S_REF", s_ref, &
340  "The default initial salinities.", units="PSU", default=35.0)
341  call get_param(param_file, mdl, "TS_RANGE_S_LIGHT", s_light, &
342  "The initial lightest salinities when COORD_CONFIG "//&
343  "is set to ts_range.", default = s_ref, units="PSU")
344  call get_param(param_file, mdl, "TS_RANGE_S_DENSE", s_dense, &
345  "The initial densest salinities when COORD_CONFIG "//&
346  "is set to ts_range.", default = s_ref, units="PSU")
347 
348  call get_param(param_file, mdl, "TS_RANGE_RESOLN_RATIO", res_rat, &
349  "The ratio of density space resolution in the densest "//&
350  "part of the range to that in the lightest part of the "//&
351  "range when COORD_CONFIG is set to ts_range. Values "//&
352  "greater than 1 increase the resolution of the denser water.",&
353  default=1.0, units="nondim")
354 
355  call get_param(param_file, mdl, "GFS", g_fs, &
356  "The reduced gravity at the free surface.", units="m s-2", &
357  default=gv%mks_g_Earth, scale=us%m_s_to_L_T**2*us%Z_to_m)
358 
359  k_light = gv%nk_rho_varies + 1
360 
361  ! Set T0(k) to range from T_LIGHT to T_DENSE, and simliarly for S0(k).
362  t0(k_light) = t_light ; s0(k_light) = s_light
363  a1 = 2.0 * res_rat / (1.0 + res_rat)
364  do k=k_light+1,nz
365  k_frac = real(k-k_light)/real(nz-k_light)
366  frac_dense = a1 * k_frac + (1.0 - a1) * k_frac**2
367  t0(k) = frac_dense * (t_dense - t_light) + t_light
368  s0(k) = frac_dense * (s_dense - s_light) + s_light
369  enddo
370 
371  g_prime(1) = g_fs
372  do k=1,nz ; pref(k) = p_ref ; enddo
373  call calculate_density(t0, s0, pref, rlay, k_light, nz-k_light+1, eqn_of_state, scale=us%kg_m3_to_R)
374  ! Extrapolate target densities for the variable density mixed and buffer layers.
375  do k=k_light-1,1,-1
376  rlay(k) = 2.0*rlay(k+1) - rlay(k+2)
377  enddo
378  do k=2,nz ; g_prime(k) = (gv%g_Earth/(gv%Rho0)) * (rlay(k) - rlay(k-1)) ; enddo
379 
380  call calltree_leave(trim(mdl)//'()')
381 end subroutine set_coord_from_ts_range
382 
383 ! Sets the layer densities (Rlay) and the interface reduced gravities (g) from data in file.
384 subroutine set_coord_from_file(Rlay, g_prime, GV, US, param_file)
385  real, dimension(:), intent(out) :: Rlay !< The layers' target coordinate values
386  !! (potential density) [R ~> kg m-3].
387  real, dimension(:), intent(out) :: g_prime !< The reduced gravity across the interfaces
388  !! [L2 Z-1 T-2 ~> m s-2].
389  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
390  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
391  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
392  ! Local variables
393  real :: g_fs ! Reduced gravity across the free surface [m s-2].
394  integer :: k, nz
395  character(len=40) :: mdl = "set_coord_from_file" ! This subroutine's name.
396  character(len=40) :: coord_var
397  character(len=200) :: filename,coord_file,inputdir ! Strings for file/path
398  nz = gv%ke
399 
400  call calltree_enter(trim(mdl)//"(), MOM_coord_initialization.F90")
401 
402  call get_param(param_file, mdl, "GFS", g_fs, &
403  "The reduced gravity at the free surface.", units="m s-2", &
404  default=gv%mks_g_Earth, scale=us%m_s_to_L_T**2*us%Z_to_m)
405  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
406  inputdir = slasher(inputdir)
407  call get_param(param_file, mdl, "COORD_FILE", coord_file, &
408  "The file from which the coordinate densities are read.", &
409  fail_if_missing=.true.)
410  call get_param(param_file, mdl, "COORD_VAR", coord_var, &
411  "The variable in COORD_FILE that is to be used for the "//&
412  "coordinate densities.", default="Layer")
413  filename = trim(inputdir)//trim(coord_file)
414  call log_param(param_file, mdl, "INPUTDIR/COORD_FILE", filename)
415  if (.not.file_exists(filename)) call mom_error(fatal, &
416  " set_coord_from_file: Unable to open "//trim(filename))
417 
418  call read_axis_data(filename, coord_var, rlay)
419  do k=1,nz ; rlay(k) = us%kg_m3_to_R*rlay(k) ; enddo
420  g_prime(1) = g_fs
421  do k=2,nz ; g_prime(k) = (gv%g_Earth/(gv%Rho0)) * (rlay(k) - rlay(k-1)) ; enddo
422  do k=1,nz ; if (g_prime(k) <= 0.0) then
423  call mom_error(fatal, "MOM_initialization set_coord_from_file: "//&
424  "Zero or negative g_primes read from variable "//"Layer"//" in file "//&
425  trim(filename))
426  endif ; enddo
427 
428  call calltree_leave(trim(mdl)//'()')
429 end subroutine set_coord_from_file
430 
431 !> Sets the layer densities (Rlay) and the interface
432 !! reduced gravities (g) according to a linear profile starting at a
433 !! reference surface layer density and spanning a range of densities
434 !! to the bottom defined by the parameter RLAY_RANGE
435 !! (defaulting to 2.0 if not defined)
436 subroutine set_coord_linear(Rlay, g_prime, GV, US, param_file)
437  real, dimension(:), intent(out) :: Rlay !< The layers' target coordinate values
438  !! (potential density) [R ~> kg m-3].
439  real, dimension(:), intent(out) :: g_prime !< The reduced gravity across the interfaces
440  !! [L2 Z-1 T-2 ~> m s-2].
441  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
442  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
443  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
444  ! Local variables
445  character(len=40) :: mdl = "set_coord_linear" ! This subroutine
446  real :: Rlay_ref, Rlay_range, g_fs
447  integer :: k, nz
448  nz = gv%ke
449 
450  call calltree_enter(trim(mdl)//"(), MOM_coord_initialization.F90")
451 
452  call get_param(param_file, mdl, "LIGHTEST_DENSITY", rlay_ref, &
453  "The reference potential density used for the surface interface.", &
454  units="kg m-3", default=us%R_to_kg_m3*gv%Rho0, scale=us%kg_m3_to_R)
455  call get_param(param_file, mdl, "DENSITY_RANGE", rlay_range, &
456  "The range of reference potential densities across all interfaces.", &
457  units="kg m-3", default=2.0, scale=us%kg_m3_to_R)
458  call get_param(param_file, mdl, "GFS", g_fs, &
459  "The reduced gravity at the free surface.", units="m s-2", &
460  default=gv%mks_g_Earth, scale=us%m_s_to_L_T**2*us%Z_to_m)
461 
462  ! This following sets the target layer densities such that a the
463  ! surface interface has density Rlay_ref and the bottom
464  ! is Rlay_range larger
465  do k=1,nz
466  rlay(k) = rlay_ref + rlay_range*((real(k)-0.5)/real(nz))
467  enddo
468  ! These statements set the interface reduced gravities.
469  g_prime(1) = g_fs
470  do k=2,nz
471  g_prime(k) = (gv%g_Earth/(gv%Rho0)) * (rlay(k) - rlay(k-1))
472  enddo
473 
474  call calltree_leave(trim(mdl)//'()')
475 end subroutine set_coord_linear
476 
477 !> Sets Rlay to Rho0 and g_prime to zero except for the free surface.
478 !! This is for use only in ALE mode where Rlay should not be used and g_prime(1) alone
479 !! might be used.
480 subroutine set_coord_to_none(Rlay, g_prime, GV, US, param_file)
481  real, dimension(:), intent(out) :: Rlay !< The layers' target coordinate values
482  !! (potential density) [R ~> kg m-3].
483  real, dimension(:), intent(out) :: g_prime !< The reduced gravity across the interfaces,
484  !! [L2 Z-1 T-2 ~> m s-2].
485  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
486  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
487  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
488  ! Local variables
489  real :: g_fs ! Reduced gravity across the free surface [m s-2].
490  character(len=40) :: mdl = "set_coord_to_none" ! This subroutine's name.
491  integer :: k, nz
492  nz = gv%ke
493 
494  call calltree_enter(trim(mdl)//"(), MOM_coord_initialization.F90")
495 
496  call get_param(param_file, mdl, "GFS" , g_fs, &
497  "The reduced gravity at the free surface.", units="m s-2", &
498  default=gv%mks_g_Earth, scale=us%m_s_to_L_T**2*us%Z_to_m)
499 
500  g_prime(1) = g_fs
501  do k=2,nz ; g_prime(k) = 0. ; enddo
502  rlay(1) = gv%Rho0
503  do k=2,nz ; rlay(k) = rlay(k-1) + g_prime(k)*(gv%Rho0/gv%g_Earth) ; enddo
504 
505  call calltree_leave(trim(mdl)//'()')
506 
507 end subroutine set_coord_to_none
508 
509 !> Writes out a file containing any available data related
510 !! to the vertical grid used by the MOM ocean model.
511 subroutine write_vertgrid_file(GV, US, param_file, directory)
512  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
513  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
514  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
515  character(len=*), intent(in) :: directory !< The directory into which to place the file.
516  ! Local variables
517  character(len=240) :: filepath
518  type(vardesc) :: vars(2)
519  type(fieldtype) :: fields(2)
520  integer :: unit
521 
522  filepath = trim(directory) // trim("Vertical_coordinate")
523 
524  vars(1) = var_desc("R","kilogram meter-3","Target Potential Density",'1','L','1')
525  vars(2) = var_desc("g","meter second-2","Reduced gravity",'1','L','1')
526 
527  call create_file(unit, trim(filepath), vars, 2, fields, single_file, gv=gv)
528 
529  call write_field(unit, fields(1), us%R_to_kg_m3*gv%Rlay(:))
530  call write_field(unit, fields(2), us%L_T_to_m_s**2*us%m_to_Z*gv%g_prime(:))
531 
532  call close_file(unit)
533 
534 end subroutine write_vertgrid_file
535 
536 end module mom_coord_initialization
mom_coord_initialization::set_coord_linear
subroutine set_coord_linear(Rlay, g_prime, GV, US, param_file)
Sets the layer densities (Rlay) and the interface reduced gravities (g) according to a linear profile...
Definition: MOM_coord_initialization.F90:437
mom_verticalgrid::setverticalgridaxes
subroutine, public setverticalgridaxes(Rlay, GV, scale)
This sets the coordinate data for the "layer mode" of the isopycnal model.
Definition: MOM_verticalGrid.F90:276
mom_io::var_desc
type(vardesc) function, public var_desc(name, units, longname, hor_grid, z_grid, t_grid, cmor_field_name, cmor_units, cmor_longname, conversion, caller)
Returns a vardesc type whose elements have been filled with the provided fields. The argument name is...
Definition: MOM_io.F90:600
mom_verticalgrid
Provides a transparent vertical ocean grid type and supporting routines.
Definition: MOM_verticalGrid.F90:2
mom_coord_initialization::write_vertgrid_file
subroutine write_vertgrid_file(GV, US, param_file, directory)
Writes out a file containing any available data related to the vertical grid used by the MOM ocean mo...
Definition: MOM_coord_initialization.F90:512
mom_file_parser::log_version
An overloaded interface to log version information about modules.
Definition: MOM_file_parser.F90:109
mom_variables::thermo_var_ptrs
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
Definition: MOM_variables.F90:78
mom_coord_initialization::set_coord_from_layer_density
subroutine set_coord_from_layer_density(Rlay, g_prime, GV, US, param_file)
Sets the layer densities (Rlay) and the interface reduced gravities (g).
Definition: MOM_coord_initialization.F90:160
mom_string_functions
Handy functions for manipulating strings.
Definition: MOM_string_functions.F90:2
mom_coord_initialization::mom_initialize_coord
subroutine, public mom_initialize_coord(GV, US, PF, write_geom, output_dir, tv, max_depth)
MOM_initialize_coord sets up time-invariant quantities related to MOM6's vertical coordinate.
Definition: MOM_coord_initialization.F90:40
mom_coord_initialization::set_coord_from_ts_ref
subroutine set_coord_from_ts_ref(Rlay, g_prime, GV, US, param_file, eqn_of_state, P_Ref)
Sets the layer densities (Rlay) and the interface reduced gravities (g) from a profile of g'.
Definition: MOM_coord_initialization.F90:203
mom_error_handler::mom_mesg
subroutine, public mom_mesg(message, verb, all_print)
This provides a convenient interface for writing an informative comment.
Definition: MOM_error_handler.F90:53
mom_eos
Provides subroutines for quantities specific to the equation of state.
Definition: MOM_EOS.F90:2
mom_file_parser::param_file_type
A structure that can be parsed to read and document run-time parameters.
Definition: MOM_file_parser.F90:54
mom_file_parser::get_param
An overloaded interface to read and log the values of various types of parameters.
Definition: MOM_file_parser.F90:102
mom_io
This module contains I/O framework code.
Definition: MOM_io.F90:2
mom_string_functions::uppercase
character(len=len(input_string)) function, public uppercase(input_string)
Return a string in which all uppercase letters have been replaced by their lowercase counterparts.
Definition: MOM_string_functions.F90:42
mom_io::read_axis_data
subroutine, public read_axis_data(filename, axis_name, var)
Read the data associated with a named axis in a file.
Definition: MOM_io.F90:438
mom_unit_scaling::unit_scale_type
Describes various unit conversion factors.
Definition: MOM_unit_scaling.F90:14
mom_coord_initialization::set_coord_from_file
subroutine set_coord_from_file(Rlay, g_prime, GV, US, param_file)
Definition: MOM_coord_initialization.F90:385
mom_eos::eos_type
A control structure for the equation of state.
Definition: MOM_EOS.F90:91
mom_coord_initialization::set_coord_to_none
subroutine set_coord_to_none(Rlay, g_prime, GV, US, param_file)
Sets Rlay to Rho0 and g_prime to zero except for the free surface. This is for use only in ALE mode w...
Definition: MOM_coord_initialization.F90:481
mom_coord_initialization
Initializes fixed aspects of the related to its vertical coordinate.
Definition: MOM_coord_initialization.F90:2
mom_verticalgrid::verticalgrid_type
Describes the vertical ocean grid, including unit conversion factors.
Definition: MOM_verticalGrid.F90:24
mom_coord_initialization::set_coord_from_ts_profile
subroutine set_coord_from_ts_profile(Rlay, g_prime, GV, US, param_file, eqn_of_state, P_Ref)
Sets the layer densities (Rlay) and the interface reduced gravities (g) from a T-S profile.
Definition: MOM_coord_initialization.F90:254
mom_variables
Provides transparent structures with groups of MOM6 variables and supporting routines.
Definition: MOM_variables.F90:2
user_initialization::user_set_coord
subroutine, public user_set_coord(Rlay, g_prime, GV, US, param_file, eqn_of_state)
Set vertical coordinates.
Definition: user_initialization.F90:41
bfb_initialization::bfb_set_coord
subroutine, public bfb_set_coord(Rlay, g_prime, GV, US, param_file, eqn_of_state)
This subroutine specifies the vertical coordinate in terms of temperature at the surface and at the b...
Definition: BFB_initialization.F90:39
mom_io::mom_read_data
Read a data field from a file.
Definition: MOM_io.F90:74
mom_coord_initialization::set_coord_from_ts_range
subroutine set_coord_from_ts_range(Rlay, g_prime, GV, US, param_file, eqn_of_state, P_Ref)
Sets the layer densities (Rlay) and the interface reduced gravities (g) from a linear T-S profile.
Definition: MOM_coord_initialization.F90:302
mom_file_parser
The MOM6 facility to parse input files for runtime parameters.
Definition: MOM_file_parser.F90:2
bfb_initialization
Initialization of the boundary-forced-basing configuration.
Definition: BFB_initialization.F90:2
mom_error_handler::calltree_leave
subroutine, public calltree_leave(mesg)
Writes a message about leaving a subroutine if call tree reporting is active.
Definition: MOM_error_handler.F90:151
mom_io::create_file
subroutine, public create_file(unit, filename, vars, novars, fields, threading, timeunit, G, dG, GV, checksums)
Routine creates a new NetCDF file. It also sets up structures that describe this file and variables t...
Definition: MOM_io.F90:93
mom_unit_scaling
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Definition: MOM_unit_scaling.F90:2
mom_error_handler::is_root_pe
logical function, public is_root_pe()
This returns .true. if the current PE is the root PE.
Definition: MOM_error_handler.F90:44
mom_debugging
Provides checksumming functions for debugging.
Definition: MOM_debugging.F90:7
mom_io::vardesc
Type for describing a variable, typically a tracer.
Definition: MOM_io.F90:53
mom_io::file_exists
Indicate whether a file exists, perhaps with domain decomposition.
Definition: MOM_io.F90:68
mom_coord_initialization::set_coord_from_gprime
subroutine set_coord_from_gprime(Rlay, g_prime, GV, US, param_file)
Sets the layer densities (Rlay) and the interface reduced gravities (g).
Definition: MOM_coord_initialization.F90:126
mom_error_handler::mom_error
subroutine, public mom_error(level, message, all_print)
This provides a convenient interface for writing an mpp_error message with run-time filter based on a...
Definition: MOM_error_handler.F90:72
mom_error_handler::calltree_waypoint
subroutine, public calltree_waypoint(mesg, n)
Writes a message about reaching a milestone if call tree reporting is active.
Definition: MOM_error_handler.F90:161
mom_file_parser::log_param
An overloaded interface to log the values of various types of parameters.
Definition: MOM_file_parser.F90:96
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2
mom_error_handler::calltree_enter
subroutine, public calltree_enter(mesg, n)
Writes a message about entering a subroutine if call tree reporting is active.
Definition: MOM_error_handler.F90:130
mom_coord_initialization::mdl
character(len=40) mdl
This module's name.
Definition: MOM_coord_initialization.F90:33
user_initialization
A template of a user to code up customized initial conditions.
Definition: user_initialization.F90:2
mom_eos::calculate_density
Calculates density of sea water from T, S and P.
Definition: MOM_EOS.F90:60
mom_file_parser::read_param
An overloaded interface to read various types of parameters.
Definition: MOM_file_parser.F90:90