24 implicit none ;
private
26 public mom_initialize_coord
33 character(len=40) :: mdl =
"MOM_coord_initialization"
39 subroutine mom_initialize_coord(GV, US, PF, write_geom, output_dir, tv, max_depth)
44 logical,
intent(in) :: write_geom
45 character(len=*),
intent(in) :: output_dir
47 real,
intent(in) :: max_depth
49 character(len=200) :: config
52 #include "version_variable.h"
54 type(
eos_type),
pointer :: eos => null()
56 if (
associated(tv%eqn_of_state)) eos => tv%eqn_of_state
60 call calltree_enter(
"MOM_initialize_coord(), MOM_coord_initialization.F90")
62 call get_param(pf, mdl,
"DEBUG", debug, default=.false.)
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.", &
84 select case ( trim(config) )
86 call set_coord_from_gprime(gv%Rlay, gv%g_prime, gv, us, pf)
88 call set_coord_from_layer_density(gv%Rlay, gv%g_prime, gv, us, pf)
90 call set_coord_linear(gv%Rlay, gv%g_prime, gv, us, pf)
92 call set_coord_from_ts_ref(gv%Rlay, gv%g_prime, gv, us, pf, eos, tv%P_Ref)
94 call set_coord_from_ts_profile(gv%Rlay, gv%g_prime, gv, us, pf, eos, tv%P_Ref)
96 call set_coord_from_ts_range(gv%Rlay, gv%g_prime, gv, us, pf, eos, tv%P_Ref)
98 call set_coord_from_file(gv%Rlay, gv%g_prime, gv, us, pf)
100 call user_set_coord(gv%Rlay, gv%g_prime, gv, us, pf, eos)
102 call bfb_set_coord(gv%Rlay, gv%g_prime, gv, us, pf, eos)
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))
110 gv%g_prime(nz+1) = 10.0*gv%g_Earth
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 )
117 gv%max_depth = max_depth
120 if (write_geom)
call write_vertgrid_file(gv, us, pf, output_dir)
122 call calltree_leave(
'MOM_initialize_coord()')
124 end subroutine mom_initialize_coord
129 subroutine set_coord_from_gprime(Rlay, g_prime, GV, US, param_file)
131 real,
dimension(GV%ke),
intent(out) :: Rlay
133 real,
dimension(GV%ke+1),
intent(out) :: g_prime
140 character(len=40) :: mdl =
"set_coord_from_gprime"
144 call calltree_enter(trim(mdl)//
"(), MOM_coord_initialization.F90")
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)
154 do k=2,nz ; g_prime(k) = g_int ;
enddo
156 do k=2,nz ; rlay(k) = rlay(k-1) + g_prime(k)*(gv%Rho0/gv%g_Earth) ;
enddo
158 call calltree_leave(trim(mdl)//
'()')
160 end subroutine set_coord_from_gprime
163 subroutine set_coord_from_layer_density(Rlay, g_prime, GV, US, param_file)
165 real,
dimension(GV%ke),
intent(out) :: Rlay
167 real,
dimension(GV%ke+1),
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", &
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)
194 rlay(k) = rlay(k-1) + rlay_range/(real(nz-1))
199 g_prime(k) = (gv%g_Earth/(gv%Rho0)) * (rlay(k) - rlay(k-1))
202 call calltree_leave(trim(mdl)//
'()')
203 end subroutine set_coord_from_layer_density
206 subroutine set_coord_from_ts_ref(Rlay, g_prime, GV, US, param_file, eqn_of_state, P_Ref)
208 real,
dimension(GV%ke),
intent(out) :: Rlay
210 real,
dimension(GV%ke+1),
intent(out) :: g_prime
214 type(
eos_type),
pointer :: eqn_of_state
215 real,
intent(in) :: P_Ref
223 character(len=40) :: mdl =
"set_coord_from_TS_ref"
227 call calltree_enter(trim(mdl)//
"(), MOM_coord_initialization.F90")
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)
243 do k=2,nz ; g_prime(k) = g_int ;
enddo
251 do k=2,nz ; rlay(k) = rlay(k-1) + g_prime(k)*(gv%Rho0/gv%g_Earth) ;
enddo
253 call calltree_leave(trim(mdl)//
'()')
254 end subroutine set_coord_from_ts_ref
257 subroutine set_coord_from_ts_profile(Rlay, g_prime, GV, US, param_file, eqn_of_state, P_Ref)
259 real,
dimension(GV%ke),
intent(out) :: Rlay
260 real,
dimension(GV%ke+1),
intent(out) :: g_prime
264 type(
eos_type),
pointer :: eqn_of_state
265 real,
intent(in) :: P_Ref
269 real,
dimension(GV%ke) :: T0, S0, Pref
272 character(len=40) :: mdl =
"set_coord_from_TS_profile"
273 character(len=200) :: filename, coord_file, inputdir
276 call calltree_enter(trim(mdl)//
"(), MOM_coord_initialization.F90")
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.)
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)
292 if (.not.
file_exists(filename))
call mom_error(fatal, &
293 " set_coord_from_TS_profile: Unable to open " //trim(filename))
296 do k=1,nz ; pref(k) = p_ref ;
enddo
298 do k=2,nz; g_prime(k) = (gv%g_Earth/(gv%Rho0)) * (rlay(k) - rlay(k-1)) ;
enddo
300 call calltree_leave(trim(mdl)//
'()')
301 end subroutine set_coord_from_ts_profile
304 subroutine set_coord_from_ts_range(Rlay, g_prime, GV, US, param_file, eqn_of_state, P_Ref)
306 real,
dimension(GV%ke),
intent(out) :: Rlay
307 real,
dimension(GV%ke+1),
intent(out) :: g_prime
311 type(
eos_type),
pointer :: eqn_of_state
312 real,
intent(in) :: P_Ref
316 real,
dimension(GV%ke) :: T0, S0, Pref
317 real :: S_Ref, S_Light, S_Dense
318 real :: T_Ref, T_Light, T_Dense
324 real :: a1, frac_dense, k_frac
325 integer :: k, nz, k_light
326 character(len=40) :: mdl =
"set_coord_from_TS_range"
327 character(len=200) :: filename, coord_file, inputdir
330 call calltree_enter(trim(mdl)//
"(), MOM_coord_initialization.F90")
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)
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")
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")
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)
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.")
364 k_light = gv%nk_rho_varies + 1
367 t0(k_light) = t_light ; s0(k_light) = s_light
368 a1 = 2.0 * res_rat / (1.0 + res_rat)
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
377 do k=1,nz ; pref(k) = p_ref ;
enddo
381 rlay(k) = 2.0*rlay(k+1) - rlay(k+2)
383 do k=2,nz ; g_prime(k) = (gv%g_Earth/gv%Rho0) * (rlay(k) - rlay(k-1)) ;
enddo
385 call calltree_leave(trim(mdl)//
'()')
386 end subroutine set_coord_from_ts_range
389 subroutine set_coord_from_file(Rlay, g_prime, GV, US, param_file)
391 real,
dimension(GV%ke),
intent(out) :: Rlay
392 real,
dimension(GV%ke+1),
intent(out) :: g_prime
400 character(len=40) :: mdl =
"set_coord_from_file"
401 character(len=40) :: coord_var
402 character(len=200) :: filename,coord_file,inputdir
405 call calltree_enter(trim(mdl)//
"(), MOM_coord_initialization.F90")
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))
423 call read_axis_data(filename, coord_var, rlay)
424 do k=1,nz ; rlay(k) = us%kg_m3_to_R*rlay(k) ;
enddo
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 "//&
433 call calltree_leave(trim(mdl)//
'()')
434 end subroutine set_coord_from_file
441 subroutine set_coord_linear(Rlay, g_prime, GV, US, param_file)
443 real,
dimension(GV%ke),
intent(out) :: Rlay
444 real,
dimension(GV%ke+1),
intent(out) :: g_prime
450 character(len=40) :: mdl =
"set_coord_linear"
451 real :: Rlay_ref, Rlay_range, g_fs
455 call calltree_enter(trim(mdl)//
"(), MOM_coord_initialization.F90")
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)
471 rlay(k) = rlay_ref + rlay_range*((real(k)-0.5)/real(nz))
476 g_prime(k) = (gv%g_Earth/(gv%Rho0)) * (rlay(k) - rlay(k-1))
479 call calltree_leave(trim(mdl)//
'()')
480 end subroutine set_coord_linear
485 subroutine set_coord_to_none(Rlay, g_prime, GV, US, param_file)
487 real,
dimension(GV%ke),
intent(out) :: Rlay
488 real,
dimension(GV%ke+1),
intent(out) :: g_prime
494 character(len=40) :: mdl =
"set_coord_to_none"
498 call calltree_enter(trim(mdl)//
"(), MOM_coord_initialization.F90")
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)
505 do k=2,nz ; g_prime(k) = 0. ;
enddo
507 do k=2,nz ; rlay(k) = rlay(k-1) + g_prime(k)*(gv%Rho0/gv%g_Earth) ;
enddo
509 call calltree_leave(trim(mdl)//
'()')
511 end subroutine set_coord_to_none
515 subroutine write_vertgrid_file(GV, US, param_file, directory)
519 character(len=*),
intent(in) :: directory
521 character(len=240) :: filepath
523 type(fieldtype) :: fields(2)
526 filepath = trim(directory) // trim(
"Vertical_coordinate")
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')
531 call create_file(unit, trim(filepath), vars, 2, fields, single_file, gv=gv)
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(:))
536 call close_file(unit)
538 end subroutine write_vertgrid_file