MOM6
MOM_coord_initialization.F90
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
9 use mom_error_handler, only : calltree_enter, calltree_leave, calltree_waypoint
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
15 use mom_string_functions, only : uppercase
18 use mom_verticalgrid, only : verticalgrid_type, setverticalgridaxes
19 use user_initialization, only : user_set_coord
20 use bfb_initialization, only : bfb_set_coord
21 
22 use netcdf
23 
24 implicit none ; private
25 
26 public mom_initialize_coord
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  default="none")
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  ! There are nz+1 values of g_prime because it is an interface field, but the value at the bottom
109  ! should not matter. This is here just to avoid having an uninitialized value in some output.
110  gv%g_prime(nz+1) = 10.0*gv%g_Earth
111 
112  if (debug) call chksum(us%R_to_kg_m3*gv%Rlay(:), "MOM_initialize_coord: Rlay ", 1, nz)
113  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)
114  call setverticalgridaxes( gv%Rlay, gv, scale=us%R_to_kg_m3 )
115 
116 ! Copy the maximum depth across from the input argument
117  gv%max_depth = max_depth
118 
119 ! Write out all of the grid data used by this run.
120  if (write_geom) call write_vertgrid_file(gv, us, pf, output_dir)
121 
122  call calltree_leave('MOM_initialize_coord()')
123 
124 end subroutine mom_initialize_coord
125 
126 ! The set_coord routines deal with initializing aspects of the vertical grid.
127 
128 !> Sets the layer densities (Rlay) and the interface reduced gravities (g).
129 subroutine set_coord_from_gprime(Rlay, g_prime, GV, US, param_file)
130  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
131  real, dimension(GV%ke), intent(out) :: Rlay !< The layers' target coordinate values
132  !! (potential density) [R ~> kg m-3].
133  real, dimension(GV%ke+1), intent(out) :: g_prime !< The reduced gravity across the interfaces
134  !! [L2 Z-1 T-2 ~> m s-2].
135  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
136  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
137  ! Local variables
138  real :: g_int ! Reduced gravities across the internal interfaces [L2 Z-1 T-2 ~> m s-2].
139  real :: g_fs ! Reduced gravity across the free surface [L2 Z-1 T-2 ~> m s-2].
140  character(len=40) :: mdl = "set_coord_from_gprime" ! This subroutine's name.
141  integer :: k, nz
142  nz = gv%ke
143 
144  call calltree_enter(trim(mdl)//"(), MOM_coord_initialization.F90")
145 
146  call get_param(param_file, mdl, "GFS" , g_fs, &
147  "The reduced gravity at the free surface.", units="m s-2", &
148  default=gv%g_Earth*us%L_T_to_m_s**2*us%m_to_Z, scale=us%m_s_to_L_T**2*us%Z_to_m)
149  call get_param(param_file, mdl, "GINT", g_int, &
150  "The reduced gravity across internal interfaces.", &
151  units="m s-2", fail_if_missing=.true., scale=us%m_s_to_L_T**2*us%Z_to_m)
152 
153  g_prime(1) = g_fs
154  do k=2,nz ; g_prime(k) = g_int ; enddo
155  rlay(1) = gv%Rho0
156  do k=2,nz ; rlay(k) = rlay(k-1) + g_prime(k)*(gv%Rho0/gv%g_Earth) ; enddo
157 
158  call calltree_leave(trim(mdl)//'()')
159 
160 end subroutine set_coord_from_gprime
161 
162 !> Sets the layer densities (Rlay) and the interface reduced gravities (g).
163 subroutine set_coord_from_layer_density(Rlay, g_prime, GV, US, param_file)
164  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
165  real, dimension(GV%ke), intent(out) :: Rlay !< The layers' target coordinate values
166  !! (potential density) [R ~> kg m-3].
167  real, dimension(GV%ke+1), intent(out) :: g_prime !< The reduced gravity across the interfaces
168  !! [L2 Z-1 T-2 ~> m s-2].
169  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
170  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
171 
172  ! Local variables
173  real :: g_fs ! Reduced gravity across the free surface [L2 Z-1 T-2 ~> m s-2].
174  real :: Rlay_Ref! The surface layer's target density [R ~> kg m-3].
175  real :: RLay_range ! The range of densities [R ~> 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*us%L_T_to_m_s**2*us%m_to_Z, scale=us%m_s_to_L_T**2*us%Z_to_m)
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=us%R_to_kg_m3*gv%Rho0, scale=us%kg_m3_to_R)
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, scale=us%kg_m3_to_R)
191 
192  rlay(1) = rlay_ref
193  do k=2,nz
194  rlay(k) = rlay(k-1) + rlay_range/(real(nz-1))
195  enddo
196 ! These statements set the interface reduced gravities. !
197  g_prime(1) = g_fs
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)//'()')
203 end subroutine set_coord_from_layer_density
204 
205 !> Sets the layer densities (Rlay) and the interface reduced gravities (g) from a profile of g'.
206 subroutine set_coord_from_ts_ref(Rlay, g_prime, GV, US, param_file, eqn_of_state, P_Ref)
207  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
208  real, dimension(GV%ke), intent(out) :: Rlay !< The layers' target coordinate values
209  !! (potential density) [R ~> kg m-3].
210  real, dimension(GV%ke+1), intent(out) :: g_prime !< The reduced gravity across the interfaces
211  !! [L2 Z-1 T-2 ~> m s-2].
212  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
213  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
214  type(eos_type), pointer :: eqn_of_state !< Equation of state structure
215  real, intent(in) :: P_Ref !< The coordinate-density reference pressure
216  !! [R L2 T-2 ~> Pa].
217 
218  ! Local variables
219  real :: T_ref ! Reference temperature
220  real :: S_ref ! Reference salinity
221  real :: g_int ! Reduced gravities across the internal interfaces [L2 Z-1 T-2 ~> m s-2].
222  real :: g_fs ! Reduced gravity across the free surface [L2 Z-1 T-2 ~> m s-2].
223  character(len=40) :: mdl = "set_coord_from_TS_ref" ! This subroutine's name.
224  integer :: k, nz
225  nz = gv%ke
226 
227  call calltree_enter(trim(mdl)//"(), MOM_coord_initialization.F90")
228 
229  call get_param(param_file, mdl, "T_REF", t_ref, &
230  "The initial temperature of the lightest layer.", units="degC", &
231  fail_if_missing=.true.)
232  call get_param(param_file, mdl, "S_REF", s_ref, &
233  "The initial salinities.", units="PSU", default=35.0)
234  call get_param(param_file, mdl, "GFS", g_fs, &
235  "The reduced gravity at the free surface.", units="m s-2", &
236  default=gv%g_Earth*us%L_T_to_m_s**2*us%m_to_Z, scale=us%m_s_to_L_T**2*us%Z_to_m)
237  call get_param(param_file, mdl, "GINT", g_int, &
238  "The reduced gravity across internal interfaces.", &
239  units="m s-2", fail_if_missing=.true., scale=us%m_s_to_L_T**2*us%Z_to_m)
240 
241 ! These statements set the interface reduced gravities. !
242  g_prime(1) = g_fs
243  do k=2,nz ; g_prime(k) = g_int ; enddo
244 
245 ! The uppermost layer's density is set here. Subsequent layers' !
246 ! densities are determined from this value and the g values. !
247 ! T0 = 28.228 ; S0 = 34.5848 ; Pref = P_Ref
248  call calculate_density(t_ref, s_ref, p_ref, rlay(1), eqn_of_state)
249 
250 ! These statements set the layer densities. !
251  do k=2,nz ; rlay(k) = rlay(k-1) + g_prime(k)*(gv%Rho0/gv%g_Earth) ; enddo
252 
253  call calltree_leave(trim(mdl)//'()')
254 end subroutine set_coord_from_ts_ref
255 
256 !> Sets the layer densities (Rlay) and the interface reduced gravities (g) from a T-S profile.
257 subroutine set_coord_from_ts_profile(Rlay, g_prime, GV, US, param_file, eqn_of_state, P_Ref)
258  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
259  real, dimension(GV%ke), intent(out) :: Rlay !< Layer potential density [R ~> kg m-3].
260  real, dimension(GV%ke+1), intent(out) :: g_prime !< The reduced gravity at each
261  !! interface [L2 Z-1 T-2 ~> m s-2].
262  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
263  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
264  type(eos_type), pointer :: eqn_of_state !< Equation of state structure
265  real, intent(in) :: P_Ref !< The coordinate-density reference pressure
266  !! [R L2 T-2 ~> Pa].
267 
268  ! Local variables
269  real, dimension(GV%ke) :: T0, S0, Pref
270  real :: g_fs ! Reduced gravity across the free surface [L2 Z-1 T-2 ~> m s-2].
271  integer :: k, nz
272  character(len=40) :: mdl = "set_coord_from_TS_profile" ! This subroutine's name.
273  character(len=200) :: filename, coord_file, inputdir ! Strings for file/path
274  nz = gv%ke
275 
276  call calltree_enter(trim(mdl)//"(), MOM_coord_initialization.F90")
277 
278  call get_param(param_file, mdl, "GFS", g_fs, &
279  "The reduced gravity at the free surface.", units="m s-2", &
280  default=gv%g_Earth*us%L_T_to_m_s**2*us%m_to_Z, scale=us%m_s_to_L_T**2*us%Z_to_m)
281  call get_param(param_file, mdl, "COORD_FILE", coord_file, &
282  "The file from which the coordinate temperatures and "//&
283  "salinities are read.", fail_if_missing=.true.)
284 
285  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
286  filename = trim(slasher(inputdir))//trim(coord_file)
287  call log_param(param_file, mdl, "INPUTDIR/COORD_FILE", filename)
288 
289  call mom_read_data(filename,"PTEMP",t0(:))
290  call mom_read_data(filename,"SALT",s0(:))
291 
292  if (.not.file_exists(filename)) call mom_error(fatal, &
293  " set_coord_from_TS_profile: Unable to open " //trim(filename))
294 ! These statements set the interface reduced gravities. !
295  g_prime(1) = g_fs
296  do k=1,nz ; pref(k) = p_ref ; enddo
297  call calculate_density(t0, s0, pref, rlay, eqn_of_state, (/1,nz/) )
298  do k=2,nz; g_prime(k) = (gv%g_Earth/(gv%Rho0)) * (rlay(k) - rlay(k-1)) ; enddo
299 
300  call calltree_leave(trim(mdl)//'()')
301 end subroutine set_coord_from_ts_profile
302 
303 !> Sets the layer densities (Rlay) and the interface reduced gravities (g) from a linear T-S profile.
304 subroutine set_coord_from_ts_range(Rlay, g_prime, GV, US, param_file, eqn_of_state, P_Ref)
305  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
306  real, dimension(GV%ke), intent(out) :: Rlay !< Layer potential density [R ~> kg m-3].
307  real, dimension(GV%ke+1), intent(out) :: g_prime !< The reduced gravity at each
308  !! interface [L2 Z-1 T-2 ~> m s-2].
309  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
310  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
311  type(eos_type), pointer :: eqn_of_state !< Equation of state structure
312  real, intent(in) :: P_Ref !< The coordinate-density reference pressure
313  !! [R L2 T-2 ~> Pa].
314 
315  ! Local variables
316  real, dimension(GV%ke) :: T0, S0, Pref
317  real :: S_Ref, S_Light, S_Dense ! Salinity range parameters [ppt].
318  real :: T_Ref, T_Light, T_Dense ! Temperature range parameters [decC].
319  real :: res_rat ! The ratio of density space resolution in the denser part
320  ! of the range to that in the lighter part of the range.
321  ! Setting this greater than 1 increases the resolution for
322  ! the denser water.
323  real :: g_fs ! Reduced gravity across the free surface [L2 Z-1 T-2 ~> m s-2].
324  real :: a1, frac_dense, k_frac
325  integer :: k, nz, k_light
326  character(len=40) :: mdl = "set_coord_from_TS_range" ! This subroutine's name.
327  character(len=200) :: filename, coord_file, inputdir ! Strings for file/path
328  nz = gv%ke
329 
330  call calltree_enter(trim(mdl)//"(), MOM_coord_initialization.F90")
331 
332  call get_param(param_file, mdl, "T_REF", t_ref, &
333  "The default initial temperatures.", units="degC", default=10.0)
334  call get_param(param_file, mdl, "TS_RANGE_T_LIGHT", t_light, &
335  "The initial temperature of the lightest layer when "//&
336  "COORD_CONFIG is set to ts_range.", units="degC", default=t_ref)
337  call get_param(param_file, mdl, "TS_RANGE_T_DENSE", t_dense, &
338  "The initial temperature of the densest layer when "//&
339  "COORD_CONFIG is set to ts_range.", units="degC", default=t_ref)
340 
341  call get_param(param_file, mdl, "S_REF", s_ref, &
342  "The default initial salinities.", units="PSU", default=35.0)
343  call get_param(param_file, mdl, "TS_RANGE_S_LIGHT", s_light, &
344  "The initial lightest salinities when COORD_CONFIG "//&
345  "is set to ts_range.", default = s_ref, units="PSU")
346  call get_param(param_file, mdl, "TS_RANGE_S_DENSE", s_dense, &
347  "The initial densest salinities when COORD_CONFIG "//&
348  "is set to ts_range.", default = s_ref, units="PSU")
349 
350  call get_param(param_file, mdl, "TS_RANGE_RESOLN_RATIO", res_rat, &
351  "The ratio of density space resolution in the densest "//&
352  "part of the range to that in the lightest part of the "//&
353  "range when COORD_CONFIG is set to ts_range. Values "//&
354  "greater than 1 increase the resolution of the denser water.",&
355  default=1.0, units="nondim")
356 
357  call get_param(param_file, mdl, "GFS", g_fs, &
358  "The reduced gravity at the free surface.", units="m s-2", &
359  default=gv%g_Earth*us%L_T_to_m_s**2*us%m_to_Z, scale=us%m_s_to_L_T**2*us%Z_to_m)
360 
361  if ((gv%nk_rho_varies > 0) .and. (nz < gv%nk_rho_varies+2)) &
362  call mom_error(fatal, "set_coord_from_TS_range requires that NZ >= NKML+NKBL+2.")
363 
364  k_light = gv%nk_rho_varies + 1
365 
366  ! Set T0(k) to range from T_LIGHT to T_DENSE, and simliarly for S0(k).
367  t0(k_light) = t_light ; s0(k_light) = s_light
368  a1 = 2.0 * res_rat / (1.0 + res_rat)
369  do k=k_light+1,nz
370  k_frac = real(k-k_light)/real(nz-k_light)
371  frac_dense = a1 * k_frac + (1.0 - a1) * k_frac**2
372  t0(k) = frac_dense * (t_dense - t_light) + t_light
373  s0(k) = frac_dense * (s_dense - s_light) + s_light
374  enddo
375 
376  g_prime(1) = g_fs
377  do k=1,nz ; pref(k) = p_ref ; enddo
378  call calculate_density(t0, s0, pref, rlay, eqn_of_state, (/k_light,nz/) )
379  ! Extrapolate target densities for the variable density mixed and buffer layers.
380  do k=k_light-1,1,-1
381  rlay(k) = 2.0*rlay(k+1) - rlay(k+2)
382  enddo
383  do k=2,nz ; g_prime(k) = (gv%g_Earth/gv%Rho0) * (rlay(k) - rlay(k-1)) ; enddo
384 
385  call calltree_leave(trim(mdl)//'()')
386 end subroutine set_coord_from_ts_range
387 
388 ! Sets the layer densities (Rlay) and the interface reduced gravities (g) from data in file.
389 subroutine set_coord_from_file(Rlay, g_prime, GV, US, param_file)
390  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
391  real, dimension(GV%ke), intent(out) :: Rlay !< Layer potential density [R ~> kg m-3].
392  real, dimension(GV%ke+1), intent(out) :: g_prime !< The reduced gravity at each
393  !! interface [L2 Z-1 T-2 ~> m s-2].
394  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
395  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
396 
397  ! Local variables
398  real :: g_fs ! Reduced gravity across the free surface [L2 Z-1 T-2 ~> m s-2].
399  integer :: k, nz
400  character(len=40) :: mdl = "set_coord_from_file" ! This subroutine's name.
401  character(len=40) :: coord_var
402  character(len=200) :: filename,coord_file,inputdir ! Strings for file/path
403  nz = gv%ke
404 
405  call calltree_enter(trim(mdl)//"(), MOM_coord_initialization.F90")
406 
407  call get_param(param_file, mdl, "GFS", g_fs, &
408  "The reduced gravity at the free surface.", units="m s-2", &
409  default=gv%g_Earth*us%L_T_to_m_s**2*us%m_to_Z, scale=us%m_s_to_L_T**2*us%Z_to_m)
410  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
411  inputdir = slasher(inputdir)
412  call get_param(param_file, mdl, "COORD_FILE", coord_file, &
413  "The file from which the coordinate densities are read.", &
414  fail_if_missing=.true.)
415  call get_param(param_file, mdl, "COORD_VAR", coord_var, &
416  "The variable in COORD_FILE that is to be used for the "//&
417  "coordinate densities.", default="Layer")
418  filename = trim(inputdir)//trim(coord_file)
419  call log_param(param_file, mdl, "INPUTDIR/COORD_FILE", filename)
420  if (.not.file_exists(filename)) call mom_error(fatal, &
421  " set_coord_from_file: Unable to open "//trim(filename))
422 
423  call read_axis_data(filename, coord_var, rlay)
424  do k=1,nz ; rlay(k) = us%kg_m3_to_R*rlay(k) ; enddo
425  g_prime(1) = g_fs
426  do k=2,nz ; g_prime(k) = (gv%g_Earth/(gv%Rho0)) * (rlay(k) - rlay(k-1)) ; enddo
427  do k=1,nz ; if (g_prime(k) <= 0.0) then
428  call mom_error(fatal, "MOM_initialization set_coord_from_file: "//&
429  "Zero or negative g_primes read from variable "//"Layer"//" in file "//&
430  trim(filename))
431  endif ; enddo
432 
433  call calltree_leave(trim(mdl)//'()')
434 end subroutine set_coord_from_file
435 
436 !> Sets the layer densities (Rlay) and the interface
437 !! reduced gravities (g) according to a linear profile starting at a
438 !! reference surface layer density and spanning a range of densities
439 !! to the bottom defined by the parameter RLAY_RANGE
440 !! (defaulting to 2.0 if not defined)
441 subroutine set_coord_linear(Rlay, g_prime, GV, US, param_file)
442  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
443  real, dimension(GV%ke), intent(out) :: Rlay !< Layer potential density [R ~> kg m-3].
444  real, dimension(GV%ke+1), intent(out) :: g_prime !< The reduced gravity at each
445  !! interface [L2 Z-1 T-2 ~> m s-2].
446  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
447  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
448 
449  ! Local variables
450  character(len=40) :: mdl = "set_coord_linear" ! This subroutine
451  real :: Rlay_ref, Rlay_range, g_fs
452  integer :: k, nz
453  nz = gv%ke
454 
455  call calltree_enter(trim(mdl)//"(), MOM_coord_initialization.F90")
456 
457  call get_param(param_file, mdl, "LIGHTEST_DENSITY", rlay_ref, &
458  "The reference potential density used for the surface interface.", &
459  units="kg m-3", default=us%R_to_kg_m3*gv%Rho0, scale=us%kg_m3_to_R)
460  call get_param(param_file, mdl, "DENSITY_RANGE", rlay_range, &
461  "The range of reference potential densities across all interfaces.", &
462  units="kg m-3", default=2.0, scale=us%kg_m3_to_R)
463  call get_param(param_file, mdl, "GFS", g_fs, &
464  "The reduced gravity at the free surface.", units="m s-2", &
465  default=gv%g_Earth*us%L_T_to_m_s**2*us%m_to_Z, scale=us%m_s_to_L_T**2*us%Z_to_m)
466 
467  ! This following sets the target layer densities such that a the
468  ! surface interface has density Rlay_ref and the bottom
469  ! is Rlay_range larger
470  do k=1,nz
471  rlay(k) = rlay_ref + rlay_range*((real(k)-0.5)/real(nz))
472  enddo
473  ! These statements set the interface reduced gravities.
474  g_prime(1) = g_fs
475  do k=2,nz
476  g_prime(k) = (gv%g_Earth/(gv%Rho0)) * (rlay(k) - rlay(k-1))
477  enddo
478 
479  call calltree_leave(trim(mdl)//'()')
480 end subroutine set_coord_linear
481 
482 !> Sets Rlay to Rho0 and g_prime to zero except for the free surface.
483 !! This is for use only in ALE mode where Rlay should not be used and g_prime(1) alone
484 !! might be used.
485 subroutine set_coord_to_none(Rlay, g_prime, GV, US, param_file)
486  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
487  real, dimension(GV%ke), intent(out) :: Rlay !< Layer potential density [R ~> kg m-3].
488  real, dimension(GV%ke+1), intent(out) :: g_prime !< The reduced gravity at each
489  !! interface [L2 Z-1 T-2 ~> m s-2].
490  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
491  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
492  ! Local variables
493  real :: g_fs ! Reduced gravity across the free surface [L2 Z-1 T-2 ~> m s-2].
494  character(len=40) :: mdl = "set_coord_to_none" ! This subroutine's name.
495  integer :: k, nz
496  nz = gv%ke
497 
498  call calltree_enter(trim(mdl)//"(), MOM_coord_initialization.F90")
499 
500  call get_param(param_file, mdl, "GFS" , g_fs, &
501  "The reduced gravity at the free surface.", units="m s-2", &
502  default=gv%g_Earth*us%L_T_to_m_s**2*us%m_to_Z, scale=us%m_s_to_L_T**2*us%Z_to_m)
503 
504  g_prime(1) = g_fs
505  do k=2,nz ; g_prime(k) = 0. ; enddo
506  rlay(1) = gv%Rho0
507  do k=2,nz ; rlay(k) = rlay(k-1) + g_prime(k)*(gv%Rho0/gv%g_Earth) ; enddo
508 
509  call calltree_leave(trim(mdl)//'()')
510 
511 end subroutine set_coord_to_none
512 
513 !> Writes out a file containing any available data related
514 !! to the vertical grid used by the MOM ocean model.
515 subroutine write_vertgrid_file(GV, US, param_file, directory)
516  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
517  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
518  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
519  character(len=*), intent(in) :: directory !< The directory into which to place the file.
520  ! Local variables
521  character(len=240) :: filepath
522  type(vardesc) :: vars(2)
523  type(fieldtype) :: fields(2)
524  integer :: unit
525 
526  filepath = trim(directory) // trim("Vertical_coordinate")
527 
528  vars(1) = var_desc("R","kilogram meter-3","Target Potential Density",'1','L','1')
529  vars(2) = var_desc("g","meter second-2","Reduced gravity",'1','L','1')
530 
531  call create_file(unit, trim(filepath), vars, 2, fields, single_file, gv=gv)
532 
533  call write_field(unit, fields(1), us%R_to_kg_m3*gv%Rlay(:))
534  call write_field(unit, fields(2), us%L_T_to_m_s**2*us%m_to_Z*gv%g_prime(:))
535 
536  call close_file(unit)
537 
538 end subroutine write_vertgrid_file
539 
540 end module mom_coord_initialization
Calculates density of sea water from T, S and P.
Definition: MOM_EOS.F90:68
A structure that can be parsed to read and document run-time parameters.
This module contains I/O framework code.
Definition: MOM_io.F90:2
The MOM6 facility to parse input files for runtime parameters.
An overloaded interface to log the values of various types of parameters.
A template of a user to code up customized initial conditions.
Describes various unit conversion factors.
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Routines for error handling and I/O management.
Provides subroutines for quantities specific to the equation of state.
Definition: MOM_EOS.F90:2
Type for describing a variable, typically a tracer.
Definition: MOM_io.F90:53
An overloaded interface to log version information about modules.
Describes the vertical ocean grid, including unit conversion factors.
Indicate whether a file exists, perhaps with domain decomposition.
Definition: MOM_io.F90:68
Initializes fixed aspects of the related to its vertical coordinate.
Initialization of the boundary-forced-basing configuration.
An overloaded interface to read various types of parameters.
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
Handy functions for manipulating strings.
Provides a transparent vertical ocean grid type and supporting routines.
Provides transparent structures with groups of MOM6 variables and supporting routines.
Read a data field from a file.
Definition: MOM_io.F90:74
An overloaded interface to read and log the values of various types of parameters.
Provides checksumming functions for debugging.
A control structure for the equation of state.
Definition: MOM_EOS.F90:108