23 implicit none ;
private 27 character(len=40) ::
mdl =
"MOM_coord_initialization" 38 logical,
intent(in) :: write_geom
39 character(len=*),
intent(in) :: output_dir
41 real,
intent(in) :: max_depth
43 character(len=200) :: config
46 #include "version_variable.h" 48 type(
eos_type),
pointer :: eos => null()
50 if (
associated(tv%eqn_of_state)) eos => tv%eqn_of_state
54 call calltree_enter(
"MOM_initialize_coord(), MOM_coord_initialization.F90")
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) )
97 case default ;
call mom_error(fatal,
"MOM_initialize_coord: "// &
98 "Unrecognized coordinate setup"//trim(config))
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)
105 gv%max_depth = max_depth
118 real,
dimension(:),
intent(out) :: Rlay
120 real,
dimension(:),
intent(out) :: g_prime
134 character(len=40) :: mdl =
"set_coord_from_gprime" 138 call calltree_enter(trim(mdl)//
"(), MOM_coord_initialization.F90")
140 call get_param(param_file, mdl,
"GFS" , g_fs, &
141 "The reduced gravity at the free surface.", units=
"m s-2", &
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.)
148 do k=2,nz ; g_prime(k) = g_int ;
enddo 150 do k=2,nz ; rlay(k) = rlay(k-1) + g_prime(k)*(gv%Rho0/gv%g_Earth) ;
enddo 159 real,
dimension(:),
intent(out) :: Rlay
161 real,
dimension(:),
intent(out) :: g_prime
176 character(len=40) :: mdl =
"set_coord_from_layer_density" 180 call calltree_enter(trim(mdl)//
"(), MOM_coord_initialization.F90")
182 call get_param(param_file, mdl,
"GFS", g_fs, &
183 "The reduced gravity at the free surface.", units=
"m s-2", &
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)
195 rlay(k) = rlay(k-1) + rlay_range/(
real(nz-1))
199 g_prime(k) = (gv%g_Earth/gv%Rho0) * (rlay(k) - rlay(k-1))
209 real,
dimension(:),
intent(out) :: Rlay
211 real,
dimension(:),
intent(out) :: g_prime
216 type(
eos_type),
pointer :: eqn_of_state
217 real,
intent(in) :: P_Ref
233 character(len=40) :: mdl =
"set_coord_from_TS_ref" 237 call calltree_enter(trim(mdl)//
"(), MOM_coord_initialization.F90")
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", &
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.)
253 do k=2,nz ; g_prime(k) = g_int ;
enddo 261 do k=2,nz ; rlay(k) = rlay(k-1) + g_prime(k)*(gv%Rho0/gv%g_Earth) ;
enddo 270 real,
dimension(:),
intent(out) :: Rlay
272 real,
dimension(:),
intent(out) :: g_prime
277 type(
eos_type),
pointer :: eqn_of_state
278 real,
intent(in) :: P_Ref
290 real,
dimension(GV%ke) :: T0, S0, Pref
293 character(len=40) :: mdl =
"set_coord_from_TS_profile" 294 character(len=200) :: filename, coord_file, inputdir
297 call calltree_enter(trim(mdl)//
"(), MOM_coord_initialization.F90")
299 call get_param(param_file, mdl,
"GFS", g_fs, &
300 "The reduced gravity at the free surface.", units=
"m s-2", &
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.)
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)
310 call read_data(filename,
"PTEMP",t0(:))
311 call read_data(filename,
"SALT",s0(:))
313 if (.not.
file_exists(filename))
call mom_error(fatal, &
314 " set_coord_from_TS_profile: Unable to open " //trim(filename))
317 do k=1,nz ; pref(k) = p_ref ;
enddo 319 do k=2,nz; g_prime(k) = (gv%g_Earth/gv%Rho0) * (rlay(k) - rlay(k-1));
enddo 328 real,
dimension(:),
intent(out) :: Rlay
330 real,
dimension(:),
intent(out) :: g_prime
335 type(
eos_type),
pointer :: eqn_of_state
336 real,
intent(in) :: P_Ref
348 real,
dimension(GV%ke) :: T0, S0, Pref
349 real :: S_Ref, S_Light, S_Dense
350 real :: T_Ref, T_Light, T_Dense
356 real :: a1, frac_dense, k_frac
357 integer :: k, nz, k_light
358 character(len=40) :: mdl =
"set_coord_from_TS_range" 359 character(len=200) :: filename, coord_file, inputdir
362 call calltree_enter(trim(mdl)//
"(), MOM_coord_initialization.F90")
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)
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")
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")
389 call get_param(param_file, mdl,
"GFS", g_fs, &
390 "The reduced gravity at the free surface.", units=
"m s-2", &
393 k_light = gv%nk_rho_varies + 1
396 t0(k_light) = t_light ; s0(k_light) = s_light
397 a1 = 2.0 * res_rat / (1.0 + res_rat)
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
406 do k=1,nz ; pref(k) = p_ref ;
enddo 410 rlay(k) = 2.0*rlay(k+1) - rlay(k+2)
412 do k=2,nz; g_prime(k) = (gv%g_Earth/gv%Rho0) * (rlay(k) - rlay(k-1));
enddo 420 real,
dimension(:),
intent(out) :: Rlay
422 real,
dimension(:),
intent(out) :: g_prime
436 character(len=40) :: mdl =
"set_coord_from_file" 437 character(len=40) :: coord_var
438 character(len=200) :: filename,coord_file,inputdir
441 call calltree_enter(trim(mdl)//
"(), MOM_coord_initialization.F90")
443 call get_param(param_file, mdl,
"GFS", g_fs, &
444 "The reduced gravity at the free surface.", units=
"m s-2", &
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))
459 call read_axis_data(filename, coord_var, rlay)
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 "//&
474 real,
dimension(:),
intent(out) :: Rlay
476 real,
dimension(:),
intent(out) :: g_prime
491 character(len=40) :: mdl =
"set_coord_linear" 492 real :: Rlay_ref, Rlay_range, g_fs
496 call calltree_enter(trim(mdl)//
"(), MOM_coord_initialization.F90")
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", &
512 rlay(k) = rlay_ref + rlay_range*((
real(k)-0.5)/
real(nz))
517 g_prime(k) = (gv%g_Earth/gv%Rho0) * (rlay(k) - rlay(k-1))
530 character(len=*),
intent(in) :: directory
537 character(len=240) :: filepath
539 type(fieldtype) :: fields(2)
542 filepath = trim(directory) // trim(
"Vertical_coordinate")
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')
547 call create_file(unit, trim(filepath), vars, 2, fields, single_file, gv=gv)
549 call write_field(unit, fields(1), gv%Rlay)
550 call write_field(unit, fields(2), gv%g_prime)
552 call close_file(unit)
subroutine set_coord_from_layer_density(Rlay, g_prime, GV, param_file)
subroutine, public read_axis_data(filename, axis_name, var)
type(vardesc) function, public var_desc(name, units, longname, hor_grid, z_grid, t_grid, cmor_field_name, cmor_units, conversion, caller)
Returns a vardesc type whose elements have been filled with the provided fields. The argument name is...
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 ...
subroutine set_coord_from_ts_profile(Rlay, g_prime, GV, param_file, eqn_of_state, P_Ref)
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...
subroutine set_coord_from_file(Rlay, g_prime, GV, param_file)
Calculates density of sea water from T, S and P.
subroutine set_coord_from_ts_range(Rlay, g_prime, GV, param_file, eqn_of_state, P_Ref)
subroutine, public calltree_leave(mesg)
Writes a message about leaving a subroutine if call tree reporting is active.
This module contains I/O framework code.
By Robert Hallberg, April 1994 - June 2002 *This subroutine initializes the fields for the simulation...
subroutine, public setverticalgridaxes(Rlay, GV)
This sets the coordinate data for the "layer mode" of the isopycnal model.
subroutine set_coord_from_gprime(Rlay, g_prime, GV, param_file)
subroutine set_coord_linear(Rlay, g_prime, GV, param_file)
subroutine, public calltree_waypoint(mesg, n)
Writes a message about reaching a milestone if call tree reporting is active.
character(len=len(input_string)) function, public uppercase(input_string)
logical function, public is_root_pe()
subroutine set_coord_from_ts_ref(Rlay, g_prime, GV, param_file, eqn_of_state, P_Ref)
subroutine, public mom_mesg(message, verb, all_print)
Provides subroutines for quantities specific to the equation of state.
Type for describing a variable, typically a tracer.
subroutine, public bfb_set_coord(Rlay, g_prime, GV, param_file, eqn_of_state)
Initializes fixed aspects of the related to its vertical coordinate.
The thermo_var_ptrs structure contains pointers to an assortment of thermodynamic fields that may be ...
subroutine, public user_set_coord(Rlay, g_prime, GV, param_file, eqn_of_state)
Set vertical coordinates.
subroutine, public create_file(unit, filename, vars, novars, fields, threading, timeunit, G, dG, GV)
Routine creates a new NetCDF file. It also sets up structures that describe this file and variables t...
subroutine, public mom_error(level, message, all_print)
A control structure for the equation of state.
subroutine, public calltree_enter(mesg, n)
Writes a message about entering a subroutine if call tree reporting is active.