10 use mom_cpu_clock,
only : cpu_clock_id, cpu_clock_begin, cpu_clock_end
13 use mom_domains,
only : root_pe, to_all, scalar_pair, cgrid_ne, agrid
32 use mom_sponge,
only : set_up_sponge_field, set_up_sponge_ml_density
42 use mom_eos,
only : convert_temp_salt_for_teos10
52 use rgc_initialization,
only : rgc_initialize_sponges
88 use mom_ale,
only : ale_initregridding,
ale_cs, ale_initthicknesstocoord
89 use mom_ale,
only : ale_remap_scalar, ale_build_grid, ale_regrid_accelerated
90 use mom_ale,
only : ts_plm_edge_values
97 implicit none ;
private
99 #include <MOM_memory.h>
101 public mom_initialize_state
108 character(len=40) :: mdl =
"MOM_state_initialization"
114 subroutine mom_initialize_state(u, v, h, tv, Time, G, GV, US, PF, dirs, &
115 restart_CS, ALE_CSp, tracer_Reg, sponge_CSp, &
116 ALE_sponge_CSp, OBC, Time_in)
120 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
123 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
126 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
130 type(time_type),
intent(inout) :: time
137 type(
ale_cs),
pointer :: ale_csp
142 type(time_type),
optional,
intent(in) :: time_in
145 character(len=200) :: filename
146 character(len=200) :: filename2
147 character(len=200) :: inputdir
148 character(len=200) :: config
154 logical :: from_z_file, useale
156 integer :: write_geom
157 logical :: use_temperature, use_sponge, use_obc
159 logical :: depress_sfc
161 logical :: trim_ic_for_p_surf
163 logical :: regrid_accelerate
164 integer :: regrid_iterations
170 type(
eos_type),
pointer :: eos => null()
173 logical :: debug_layers = .false.
174 character(len=80) :: mesg
176 #include "version_variable.h"
177 integer :: i, j, k, is, ie, js, je, isq, ieq, jsq, jeq, nz
178 integer :: isd, ied, jsd, jed, isdb, iedb, jsdb, jedb
180 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
181 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
182 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
183 isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
185 call calltree_enter(
"MOM_initialize_state(), MOM_state_initialization.F90")
187 call get_param(pf, mdl,
"DEBUG", debug, default=.false.)
188 call get_param(pf, mdl,
"DEBUG_OBC", debug_obc, default=.false.)
190 new_sim = determine_is_new_run(dirs%input_filename, dirs%restart_input_dir, &
192 just_read = .not.new_sim
194 call get_param(pf, mdl,
"INPUTDIR", inputdir, &
195 "The directory in which input files are found.", default=
".")
196 inputdir = slasher(inputdir)
198 use_temperature =
associated(tv%T)
199 useale =
associated(ale_csp)
200 use_eos =
associated(tv%eqn_of_state)
201 use_obc =
associated(obc)
202 if (use_eos) eos => tv%eqn_of_state
210 call mom_mesg(
"Run initialized internally.", 3)
212 if (
present(time_in)) time = time_in
227 call get_param(pf, mdl,
"INIT_LAYERS_FROM_Z_FILE", from_z_file, &
228 "If true, initialize the layer thicknesses, temperatures, "//&
229 "and salinities from a Z-space file on a latitude-longitude "//&
230 "grid.", default=.false., do_not_log=just_read)
232 if (from_z_file)
then
234 if (.NOT.use_temperature)
call mom_error(fatal,
"MOM_initialize_state : "//&
235 "use_temperature must be true if INIT_LAYERS_FROM_Z_FILE is true")
237 call mom_temp_salt_initialize_from_z(h, tv, g, gv, us, pf, just_read_params=just_read)
241 call get_param(pf, mdl,
"THICKNESS_CONFIG", config, &
242 "A string that determines how the initial layer "//&
243 "thicknesses are specified for a new run: \n"//&
244 " \t file - read interface heights from the file specified \n"//&
245 " \t thickness_file - read thicknesses from the file specified \n"//&
246 " \t\t by (THICKNESS_FILE).\n"//&
247 " \t coord - determined by ALE coordinate.\n"//&
248 " \t uniform - uniform thickness layers evenly distributed \n"//&
249 " \t\t between the surface and MAXIMUM_DEPTH. \n"//&
250 " \t list - read a list of positive interface depths. \n"//&
251 " \t DOME - use a slope and channel configuration for the \n"//&
252 " \t\t DOME sill-overflow test case. \n"//&
253 " \t ISOMIP - use a configuration for the \n"//&
254 " \t\t ISOMIP test case. \n"//&
255 " \t benchmark - use the benchmark test case thicknesses. \n"//&
256 " \t Neverworld - use the Neverworld test case thicknesses. \n"//&
257 " \t search - search a density profile for the interface \n"//&
258 " \t\t densities. This is not yet implemented. \n"//&
259 " \t circle_obcs - the circle_obcs test case is used. \n"//&
260 " \t DOME2D - 2D version of DOME initialization. \n"//&
261 " \t adjustment2d - 2D lock exchange thickness ICs. \n"//&
262 " \t sloshing - sloshing gravity thickness ICs. \n"//&
263 " \t seamount - no motion test with seamount ICs. \n"//&
264 " \t dumbbell - sloshing channel ICs. \n"//&
265 " \t soliton - Equatorial Rossby soliton. \n"//&
266 " \t rossby_front - a mixed layer front in thermal wind balance.\n"//&
267 " \t USER - call a user modified routine.", &
268 default=
"uniform", do_not_log=just_read)
269 select case (trim(config))
271 call initialize_thickness_from_file(h, g, gv, us, pf, .false., just_read_params=just_read)
272 case (
"thickness_file")
273 call initialize_thickness_from_file(h, g, gv, us, pf, .true., just_read_params=just_read)
275 if (new_sim .and. useale)
then
276 call ale_initthicknesstocoord( ale_csp, g, gv, h )
277 elseif (new_sim)
then
278 call mom_error(fatal,
"MOM_initialize_state: USE_REGRIDDING must be True "//&
279 "for THICKNESS_CONFIG of 'coord'")
281 case (
"uniform");
call initialize_thickness_uniform(h, g, gv, pf, &
282 just_read_params=just_read)
283 case (
"list");
call initialize_thickness_list(h, g, gv, us, pf, &
284 just_read_params=just_read)
285 case (
"DOME");
call dome_initialize_thickness(h, g, gv, pf, &
286 just_read_params=just_read)
287 case (
"ISOMIP");
call isomip_initialize_thickness(h, g, gv, us, pf, tv, &
288 just_read_params=just_read)
289 case (
"benchmark");
call benchmark_initialize_thickness(h, g, gv, us, pf, &
290 tv%eqn_of_state, tv%P_Ref, just_read_params=just_read)
291 case (
"Neverwoorld",
"Neverland");
call neverworld_initialize_thickness(h, g, gv, us, pf, &
292 tv%eqn_of_state, tv%P_Ref)
293 case (
"search");
call initialize_thickness_search
294 case (
"circle_obcs");
call circle_obcs_initialize_thickness(h, g, gv, pf, &
295 just_read_params=just_read)
296 case (
"lock_exchange");
call lock_exchange_initialize_thickness(h, g, gv, us, &
297 pf, just_read_params=just_read)
298 case (
"external_gwave");
call external_gwave_initialize_thickness(h, g, gv, us, &
299 pf, just_read_params=just_read)
300 case (
"DOME2D");
call dome2d_initialize_thickness(h, g, gv, us, pf, &
301 just_read_params=just_read)
302 case (
"adjustment2d");
call adjustment_initialize_thickness(h, g, gv, us, &
303 pf, just_read_params=just_read)
304 case (
"sloshing");
call sloshing_initialize_thickness(h, g, gv, us, pf, &
305 just_read_params=just_read)
306 case (
"seamount");
call seamount_initialize_thickness(h, g, gv, us, pf, &
307 just_read_params=just_read)
308 case (
"dumbbell");
call dumbbell_initialize_thickness(h, g, gv, us, pf, &
309 just_read_params=just_read)
310 case (
"soliton");
call soliton_initialize_thickness(h, g, gv, us)
311 case (
"phillips");
call phillips_initialize_thickness(h, g, gv, us, pf, &
312 just_read_params=just_read)
313 case (
"rossby_front");
call rossby_front_initialize_thickness(h, g, gv, us, &
314 pf, just_read_params=just_read)
315 case (
"USER");
call user_initialize_thickness(h, g, gv, pf, &
316 just_read_params=just_read)
317 case default ;
call mom_error(fatal,
"MOM_initialize_state: "//&
318 "Unrecognized layer thickness configuration "//trim(config))
322 if ( use_temperature )
then
323 call get_param(pf, mdl,
"TS_CONFIG", config, &
324 "A string that determines how the initial tempertures "//&
325 "and salinities are specified for a new run: \n"//&
326 " \t file - read velocities from the file specified \n"//&
327 " \t\t by (TS_FILE). \n"//&
328 " \t fit - find the temperatures that are consistent with \n"//&
329 " \t\t the layer densities and salinity S_REF. \n"//&
330 " \t TS_profile - use temperature and salinity profiles \n"//&
331 " \t\t (read from TS_FILE) to set layer densities. \n"//&
332 " \t benchmark - use the benchmark test case T & S. \n"//&
333 " \t linear - linear in logical layer space. \n"//&
334 " \t DOME2D - 2D DOME initialization. \n"//&
335 " \t ISOMIP - ISOMIP initialization. \n"//&
336 " \t adjustment2d - 2d lock exchange T/S ICs. \n"//&
337 " \t sloshing - sloshing mode T/S ICs. \n"//&
338 " \t seamount - no motion test with seamount ICs. \n"//&
339 " \t dumbbell - sloshing channel ICs. \n"//&
340 " \t rossby_front - a mixed layer front in thermal wind balance.\n"//&
341 " \t SCM_CVMix_tests - used in the SCM CVMix tests.\n"//&
342 " \t USER - call a user modified routine.", &
343 fail_if_missing=new_sim, do_not_log=just_read)
345 select case (trim(config))
346 case (
"fit");
call initialize_temp_salt_fit(tv%T, tv%S, g, gv, us, pf, &
347 eos, tv%P_Ref, just_read_params=just_read)
348 case (
"file");
call initialize_temp_salt_from_file(tv%T, tv%S, g, &
349 pf, just_read_params=just_read)
350 case (
"benchmark");
call benchmark_init_temperature_salinity(tv%T, tv%S, &
351 g, gv, us, pf, eos, tv%P_Ref, just_read_params=just_read)
352 case (
"TS_profile") ;
call initialize_temp_salt_from_profile(tv%T, tv%S, &
353 g, pf, just_read_params=just_read)
354 case (
"linear");
call initialize_temp_salt_linear(tv%T, tv%S, g, pf, &
355 just_read_params=just_read)
356 case (
"DOME2D");
call dome2d_initialize_temperature_salinity ( tv%T, &
357 tv%S, h, g, gv, pf, eos, just_read_params=just_read)
358 case (
"ISOMIP");
call isomip_initialize_temperature_salinity ( tv%T, &
359 tv%S, h, g, gv, us, pf, eos, just_read_params=just_read)
360 case (
"adjustment2d");
call adjustment_initialize_temperature_salinity ( tv%T, &
361 tv%S, h, g, gv, pf, eos, just_read_params=just_read)
362 case (
"baroclinic_zone");
call baroclinic_zone_init_temperature_salinity( tv%T, &
363 tv%S, h, g, gv, us, pf, just_read_params=just_read)
364 case (
"sloshing");
call sloshing_initialize_temperature_salinity(tv%T, &
365 tv%S, h, g, gv, pf, eos, just_read_params=just_read)
366 case (
"seamount");
call seamount_initialize_temperature_salinity(tv%T, &
367 tv%S, h, g, gv, pf, eos, just_read_params=just_read)
368 case (
"dumbbell");
call dumbbell_initialize_temperature_salinity(tv%T, &
369 tv%S, h, g, gv, pf, eos, just_read_params=just_read)
370 case (
"rossby_front");
call rossby_front_initialize_temperature_salinity ( tv%T, &
371 tv%S, h, g, gv, pf, eos, just_read_params=just_read)
372 case (
"SCM_CVMix_tests");
call scm_cvmix_tests_ts_init(tv%T, tv%S, h, &
373 g, gv, us, pf, just_read_params=just_read)
374 case (
"dense");
call dense_water_initialize_ts(g, gv, pf, eos, tv%T, tv%S, &
375 h, just_read_params=just_read)
376 case (
"USER");
call user_init_temperature_salinity(tv%T, tv%S, g, pf, eos, &
377 just_read_params=just_read)
378 case default ;
call mom_error(fatal,
"MOM_initialize_state: "//&
379 "Unrecognized Temp & salt configuration "//trim(config))
383 if (use_temperature .and. use_obc) &
384 call fill_temp_salt_segments(g, obc, tv)
387 if (new_sim)
call pass_var(h, g%Domain)
390 call get_param(pf, mdl,
"VELOCITY_CONFIG", config, &
391 "A string that determines how the initial velocities "//&
392 "are specified for a new run: \n"//&
393 " \t file - read velocities from the file specified \n"//&
394 " \t\t by (VELOCITY_FILE). \n"//&
395 " \t zero - the fluid is initially at rest. \n"//&
396 " \t uniform - the flow is uniform (determined by\n"//&
397 " \t\t parameters INITIAL_U_CONST and INITIAL_V_CONST).\n"//&
398 " \t rossby_front - a mixed layer front in thermal wind balance.\n"//&
399 " \t soliton - Equatorial Rossby soliton.\n"//&
400 " \t USER - call a user modified routine.", default=
"zero", &
401 do_not_log=just_read)
402 select case (trim(config))
403 case (
"file");
call initialize_velocity_from_file(u, v, g, us, pf, &
404 just_read_params=just_read)
405 case (
"zero");
call initialize_velocity_zero(u, v, g, pf, &
406 just_read_params=just_read)
407 case (
"uniform");
call initialize_velocity_uniform(u, v, g, us, pf, &
408 just_read_params=just_read)
409 case (
"circular");
call initialize_velocity_circular(u, v, g, us, pf, &
410 just_read_params=just_read)
411 case (
"phillips");
call phillips_initialize_velocity(u, v, g, gv, us, pf, &
412 just_read_params=just_read)
413 case (
"rossby_front");
call rossby_front_initialize_velocity(u, v, h, &
414 g, gv, us, pf, just_read_params=just_read)
415 case (
"soliton");
call soliton_initialize_velocity(u, v, h, g, us)
416 case (
"USER");
call user_initialize_velocity(u, v, g, us, pf, &
417 just_read_params=just_read)
418 case default ;
call mom_error(fatal,
"MOM_initialize_state: "//&
419 "Unrecognized velocity configuration "//trim(config))
423 if (debug .and. new_sim)
then
424 call uvchksum(
"MOM_initialize_state [uv]", u, v, g%HI, haloshift=1, scale=us%m_s_to_L_T)
429 call get_param(pf, mdl,
"CONVERT_THICKNESS_UNITS", convert, &
430 "If true, convert the thickness initial conditions from "//&
431 "units of m to kg m-2 or vice versa, depending on whether "//&
432 "BOUSSINESQ is defined. This does not apply if a restart "//&
433 "file is read.", default=.not.gv%Boussinesq, do_not_log=just_read)
435 if (new_sim .and. convert .and. .not.gv%Boussinesq) &
437 call convert_thickness(h, g, gv, us, tv)
440 call get_param(pf, mdl,
"DEPRESS_INITIAL_SURFACE", depress_sfc, &
441 "If true, depress the initial surface to avoid huge "//&
442 "tsunamis when a large surface pressure is applied.", &
443 default=.false., do_not_log=just_read)
444 call get_param(pf, mdl,
"TRIM_IC_FOR_P_SURF", trim_ic_for_p_surf, &
445 "If true, cuts way the top of the column for initial conditions "//&
446 "at the depth where the hydrostatic pressure matches the imposed "//&
447 "surface pressure which is read from file.", default=.false., &
448 do_not_log=just_read)
449 if (depress_sfc .and. trim_ic_for_p_surf)
call mom_error(fatal,
"MOM_initialize_state: "//&
450 "DEPRESS_INITIAL_SURFACE and TRIM_IC_FOR_P_SURF are exclusive and cannot both be True")
451 if (new_sim .and. debug .and. (depress_sfc .or. trim_ic_for_p_surf)) &
452 call hchksum(h,
"Pre-depress: h ", g%HI, haloshift=1, scale=gv%H_to_m)
453 if (depress_sfc)
call depress_surface(h, g, gv, us, pf, tv, just_read_params=just_read)
454 if (trim_ic_for_p_surf)
call trim_for_ice(pf, g, gv, us, ale_csp, tv, h, just_read_params=just_read)
459 call get_param(pf, mdl,
"REGRID_ACCELERATE_INIT", regrid_accelerate, &
460 "If true, runs REGRID_ACCELERATE_ITERATIONS iterations of the regridding "//&
461 "algorithm to push the initial grid to be consistent with the initial "//&
462 "condition. Useful only for state-based and iterative coordinates.", &
463 default=.false., do_not_log=just_read)
464 if (regrid_accelerate)
then
465 call get_param(pf, mdl,
"REGRID_ACCELERATE_ITERATIONS", regrid_iterations, &
466 "The number of regridding iterations to perform to generate "//&
467 "an initial grid that is consistent with the initial conditions.", &
468 default=1, do_not_log=just_read)
470 call get_param(pf, mdl,
"DT", dt,
"Timestep", fail_if_missing=.true., scale=us%s_to_T)
472 if (new_sim .and. debug) &
473 call hchksum(h,
"Pre-ALE_regrid: h ", g%HI, haloshift=1, scale=gv%H_to_m)
474 call ale_regrid_accelerated(ale_csp, g, gv, h, tv, regrid_iterations, u, v, obc, tracer_reg, &
475 dt=dt, initial=.true.)
481 if (.not.new_sim)
then
484 call restore_state(dirs%input_filename, dirs%restart_input_dir, time, &
486 if (
present(time_in)) time = time_in
487 if ((gv%m_to_H_restart /= 0.0) .and. (gv%m_to_H_restart /= gv%m_to_H))
then
488 h_rescale = gv%m_to_H / gv%m_to_H_restart
489 do k=1,nz ;
do j=js,je ;
do i=is,ie ; h(i,j,k) = h_rescale * h(i,j,k) ;
enddo ;
enddo ;
enddo
491 if ( (us%s_to_T_restart * us%m_to_L_restart /= 0.0) .and. &
492 ((us%m_to_L * us%s_to_T_restart) /= (us%m_to_L_restart * us%s_to_T)) )
then
493 vel_rescale = (us%m_to_L * us%s_to_T_restart) / (us%m_to_L_restart * us%s_to_T)
494 do k=1,nz ;
do j=jsd,jed ;
do i=isdb,iedb ; u(i,j,k) = vel_rescale * u(i,j,k) ;
enddo ;
enddo ;
enddo
495 do k=1,nz ;
do j=jsdb,jedb ;
do i=isd,ied ; v(i,j,k) = vel_rescale * v(i,j,k) ;
enddo ;
enddo ;
enddo
499 if ( use_temperature )
then
500 call pass_var(tv%T, g%Domain, complete=.false.)
501 call pass_var(tv%S, g%Domain, complete=.false.)
506 call hchksum(h,
"MOM_initialize_state: h ", g%HI, haloshift=1, scale=gv%H_to_m)
507 if ( use_temperature )
call hchksum(tv%T,
"MOM_initialize_state: T ", g%HI, haloshift=1)
508 if ( use_temperature )
call hchksum(tv%S,
"MOM_initialize_state: S ", g%HI, haloshift=1)
509 if ( use_temperature .and. debug_layers)
then ;
do k=1,nz
510 write(mesg,
'("MOM_IS: T[",I2,"]")') k
511 call hchksum(tv%T(:,:,k), mesg, g%HI, haloshift=1)
512 write(mesg,
'("MOM_IS: S[",I2,"]")') k
513 call hchksum(tv%S(:,:,k), mesg, g%HI, haloshift=1)
517 call get_param(pf, mdl,
"SPONGE", use_sponge, &
518 "If true, sponges may be applied anywhere in the domain. "//&
519 "The exact location and properties of those sponges are "//&
520 "specified via SPONGE_CONFIG.", default=.false.)
521 if ( use_sponge )
then
522 call get_param(pf, mdl,
"SPONGE_CONFIG", config, &
523 "A string that sets how the sponges are configured: \n"//&
524 " \t file - read sponge properties from the file \n"//&
525 " \t\t specified by (SPONGE_FILE).\n"//&
526 " \t ISOMIP - apply ale sponge in the ISOMIP case \n"//&
527 " \t RGC - apply sponge in the rotating_gravity_current case \n"//&
528 " \t DOME - use a slope and channel configuration for the \n"//&
529 " \t\t DOME sill-overflow test case. \n"//&
530 " \t BFB - Sponge at the southern boundary of the domain\n"//&
531 " \t\t for buoyancy-forced basin case.\n"//&
532 " \t USER - call a user modified routine.", default=
"file")
533 select case (trim(config))
534 case (
"DOME");
call dome_initialize_sponges(g, gv, us, tv, pf, sponge_csp)
535 case (
"DOME2D");
call dome2d_initialize_sponges(g, gv, us, tv, pf, useale, &
536 sponge_csp, ale_sponge_csp)
537 case (
"ISOMIP");
call isomip_initialize_sponges(g, gv, us, tv, pf, useale, &
538 sponge_csp, ale_sponge_csp)
539 case(
"RGC");
call rgc_initialize_sponges(g, gv, us, tv, u, v, pf, useale, &
540 sponge_csp, ale_sponge_csp)
541 case (
"USER");
call user_initialize_sponges(g, gv, use_temperature, tv, pf, sponge_csp, h)
542 case (
"BFB");
call bfb_initialize_sponges_southonly(g, gv, us, use_temperature, tv, pf, &
544 case (
"DUMBBELL");
call dumbbell_initialize_sponges(g, gv, us, tv, pf, useale, &
545 sponge_csp, ale_sponge_csp)
546 case (
"phillips");
call phillips_initialize_sponges(g, gv, us, tv, pf, sponge_csp, h)
547 case (
"dense");
call dense_water_initialize_sponges(g, gv, us, tv, pf, useale, &
548 sponge_csp, ale_sponge_csp)
549 case (
"file");
call initialize_sponges_file(g, gv, us, use_temperature, tv, pf, &
550 sponge_csp, ale_sponge_csp, time)
551 case default ;
call mom_error(fatal,
"MOM_initialize_state: "//&
552 "Unrecognized sponge configuration "//trim(config))
557 call open_boundary_init(g, gv, us, pf, obc, restart_cs)
560 if (
associated(obc))
then
561 call initialize_segment_data(g, obc, pf)
564 if (.not. obc%needs_IO_for_data) &
565 call update_obc_segment_data(g, gv, us, obc, tv, h, time)
567 call get_param(pf, mdl,
"OBC_USER_CONFIG", config, &
568 "A string that sets how the user code is invoked to set open boundary data: \n"//&
569 " DOME - specified inflow on northern boundary\n"//&
570 " dyed_channel - supercritical with dye on the inflow boundary\n"//&
571 " dyed_obcs - circle_obcs with dyes on the open boundaries\n"//&
572 " Kelvin - barotropic Kelvin wave forcing on the western boundary\n"//&
573 " shelfwave - Flather with shelf wave forcing on western boundary\n"//&
574 " supercritical - now only needed here for the allocations\n"//&
575 " tidal_bay - Flather with tidal forcing on eastern boundary\n"//&
576 " USER - user specified", default=
"none")
577 if (trim(config) ==
"DOME")
then
578 call dome_set_obc_data(obc, tv, g, gv, us, pf, tracer_reg)
579 elseif (trim(config) ==
"dyed_channel")
then
580 call dyed_channel_set_obc_tracer_data(obc, g, gv, pf, tracer_reg)
581 obc%update_OBC = .true.
582 elseif (trim(config) ==
"dyed_obcs")
then
583 call dyed_obcs_set_obc_data(obc, g, gv, pf, tracer_reg)
584 elseif (trim(config) ==
"Kelvin")
then
585 obc%update_OBC = .true.
586 elseif (trim(config) ==
"shelfwave")
then
587 obc%update_OBC = .true.
588 elseif (lowercase(trim(config)) ==
"supercritical")
then
589 call supercritical_set_obc_data(obc, g, pf)
590 elseif (trim(config) ==
"tidal_bay")
then
591 obc%update_OBC = .true.
592 elseif (trim(config) ==
"USER")
then
593 call user_set_obc_data(obc, tv, g, pf, tracer_reg)
594 elseif (.not. trim(config) ==
"none")
then
595 call mom_error(fatal,
"The open boundary conditions specified by "//&
596 "OBC_USER_CONFIG = "//trim(config)//
" have not been fully implemented.")
598 if (open_boundary_query(obc, apply_open_obc=.true.))
then
599 call set_tracer_data(obc, tv, h, g, pf, tracer_reg)
606 if (debug.and.
associated(obc))
then
607 call hchksum(g%mask2dT,
'MOM_initialize_state: mask2dT ', g%HI)
608 call uvchksum(
'MOM_initialize_state: mask2dC[uv]', g%mask2dCu, &
610 call qchksum(g%mask2dBu,
'MOM_initialize_state: mask2dBu ', g%HI)
613 if (debug_obc)
call open_boundary_test_extern_h(g, gv, obc, h)
614 call calltree_leave(
'MOM_initialize_state()')
616 end subroutine mom_initialize_state
619 subroutine initialize_thickness_from_file(h, G, GV, US, param_file, file_has_thickness, &
624 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
628 logical,
intent(in) :: file_has_thickness
631 logical,
optional,
intent(in) :: just_read_params
635 real :: eta(SZI_(G),SZJ_(G),SZK_(G)+1)
636 integer :: inconsistent = 0
637 logical :: correct_thickness
639 character(len=40) :: mdl =
"initialize_thickness_from_file"
640 character(len=200) :: filename, thickness_file, inputdir, mesg
641 integer :: i, j, k, is, ie, js, je, nz
643 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
645 just_read = .false. ;
if (
present(just_read_params)) just_read = just_read_params
647 if (.not.just_read) &
648 call calltree_enter(trim(mdl)//
"(), MOM_state_initialization.F90")
650 call get_param(param_file, mdl,
"INPUTDIR", inputdir, default=
".", do_not_log=just_read)
651 inputdir = slasher(inputdir)
652 call get_param(param_file, mdl,
"THICKNESS_FILE", thickness_file, &
653 "The name of the thickness file.", &
654 fail_if_missing=.not.just_read, do_not_log=just_read)
656 filename = trim(inputdir)//trim(thickness_file)
657 if (.not.just_read)
call log_param(param_file, mdl,
"INPUTDIR/THICKNESS_FILE", filename)
659 if ((.not.just_read) .and. (.not.
file_exists(filename, g%Domain)))
call mom_error(fatal, &
660 " initialize_thickness_from_file: Unable to open "//trim(filename))
662 if (file_has_thickness)
then
664 if (just_read)
return
665 call mom_read_data(filename,
"h", h(:,:,:), g%Domain, scale=gv%m_to_H)
667 call get_param(param_file, mdl,
"ADJUST_THICKNESS", correct_thickness, &
668 "If true, all mass below the bottom removed if the "//&
669 "topography is shallower than the thickness input file "//&
670 "would indicate.", default=.false., do_not_log=just_read)
671 if (just_read)
return
673 call mom_read_data(filename,
"eta", eta(:,:,:), g%Domain, scale=us%m_to_Z)
675 if (correct_thickness)
then
676 call adjustetatofitbathymetry(g, gv, us, eta, h)
678 do k=nz,1,-1 ;
do j=js,je ;
do i=is,ie
679 if (eta(i,j,k) < (eta(i,j,k+1) + gv%Angstrom_Z))
then
680 eta(i,j,k) = eta(i,j,k+1) + gv%Angstrom_Z
681 h(i,j,k) = gv%Angstrom_H
683 h(i,j,k) = gv%Z_to_H * (eta(i,j,k) - eta(i,j,k+1))
685 enddo ;
enddo ;
enddo
687 do j=js,je ;
do i=is,ie
688 if (abs(eta(i,j,nz+1) + g%bathyT(i,j)) > 1.0*us%m_to_Z) &
689 inconsistent = inconsistent + 1
691 call sum_across_pes(inconsistent)
693 if ((inconsistent > 0) .and. (is_root_pe()))
then
694 write(mesg,
'("Thickness initial conditions are inconsistent ",'// &
695 '"with topography in ",I8," places.")') inconsistent
696 call mom_error(warning, mesg)
701 call calltree_leave(trim(mdl)//
'()')
702 end subroutine initialize_thickness_from_file
712 subroutine adjustetatofitbathymetry(G, GV, US, eta, h)
716 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1),
intent(inout) :: eta
717 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(inout) :: h
719 integer :: i, j, k, is, ie, js, je, nz, contractions, dilations
720 real :: hTolerance = 0.1
721 real :: hTmp, eTmp, dilate
722 character(len=100) :: mesg
724 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
725 htolerance = 0.1*us%m_to_Z
728 do j=js,je ;
do i=is,ie
729 if (-eta(i,j,nz+1) > g%bathyT(i,j) + htolerance)
then
730 eta(i,j,nz+1) = -g%bathyT(i,j)
731 contractions = contractions + 1
734 call sum_across_pes(contractions)
735 if ((contractions > 0) .and. (is_root_pe()))
then
736 write(mesg,
'("Thickness initial conditions were contracted ",'// &
737 '"to fit topography in ",I8," places.")') contractions
738 call mom_error(warning,
'adjustEtaToFitBathymetry: '//mesg)
743 do k=nz,1,-1 ;
do j=js,je ;
do i=is,ie
746 if (eta(i,j,k) < (eta(i,j,k+1) + gv%Angstrom_Z))
then
747 eta(i,j,k) = eta(i,j,k+1) + gv%Angstrom_Z
748 h(i,j,k) = gv%Angstrom_Z
750 h(i,j,k) = (eta(i,j,k) - eta(i,j,k+1))
752 enddo ;
enddo ;
enddo
755 do j=js,je ;
do i=is,ie
759 if (-eta(i,j,nz+1) < g%bathyT(i,j) - htolerance)
then
760 dilations = dilations + 1
761 if (eta(i,j,1) <= eta(i,j,nz+1))
then
762 do k=1,nz ; h(i,j,k) = (eta(i,j,1) + g%bathyT(i,j)) / real(nz) ;
enddo
764 dilate = (eta(i,j,1) + g%bathyT(i,j)) / (eta(i,j,1) - eta(i,j,nz+1))
765 do k=1,nz ; h(i,j,k) = h(i,j,k) * dilate ;
enddo
767 do k=nz,2,-1 ; eta(i,j,k) = eta(i,j,k+1) + h(i,j,k) ;
enddo
772 do k=1,nz ;
do j=js,je ;
do i=is,ie
773 h(i,j,k) = h(i,j,k)*gv%Z_to_H
774 enddo ;
enddo ;
enddo
776 call sum_across_pes(dilations)
777 if ((dilations > 0) .and. (is_root_pe()))
then
778 write(mesg,
'("Thickness initial conditions were dilated ",'// &
779 '"to fit topography in ",I8," places.")') dilations
780 call mom_error(warning,
'adjustEtaToFitBathymetry: '//mesg)
783 end subroutine adjustetatofitbathymetry
786 subroutine initialize_thickness_uniform(h, G, GV, param_file, just_read_params)
789 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
793 logical,
optional,
intent(in) :: just_read_params
796 character(len=40) :: mdl =
"initialize_thickness_uniform"
797 real :: e0(SZK_(G)+1)
799 real :: eta1D(SZK_(G)+1)
802 integer :: i, j, k, is, ie, js, je, nz
804 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
806 just_read = .false. ;
if (
present(just_read_params)) just_read = just_read_params
808 if (just_read)
return
810 call calltree_enter(trim(mdl)//
"(), MOM_state_initialization.F90")
812 if (g%max_depth<=0.)
call mom_error(fatal,
"initialize_thickness_uniform: "// &
813 "MAXIMUM_DEPTH has a non-sensical value! Was it set?")
816 e0(k) = -g%max_depth * real(k-1) / real(nz)
819 do j=js,je ;
do i=is,ie
825 eta1d(nz+1) = -g%bathyT(i,j)
828 if (eta1d(k) < (eta1d(k+1) + gv%Angstrom_Z))
then
829 eta1d(k) = eta1d(k+1) + gv%Angstrom_Z
830 h(i,j,k) = gv%Angstrom_H
832 h(i,j,k) = gv%Z_to_H * (eta1d(k) - eta1d(k+1))
837 call calltree_leave(trim(mdl)//
'()')
838 end subroutine initialize_thickness_uniform
841 subroutine initialize_thickness_list(h, G, GV, US, param_file, just_read_params)
845 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
849 logical,
optional,
intent(in) :: just_read_params
852 character(len=40) :: mdl =
"initialize_thickness_list"
853 real :: e0(SZK_(G)+1)
855 real :: eta1D(SZK_(G)+1)
858 character(len=200) :: filename, eta_file, inputdir
859 character(len=72) :: eta_var
860 integer :: i, j, k, is, ie, js, je, nz
862 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
864 just_read = .false. ;
if (
present(just_read_params)) just_read = just_read_params
866 call get_param(param_file, mdl,
"INTERFACE_IC_FILE", eta_file, &
867 "The file from which horizontal mean initial conditions "//&
868 "for interface depths can be read.", fail_if_missing=.true.)
869 call get_param(param_file, mdl,
"INTERFACE_IC_VAR", eta_var, &
870 "The variable name for horizontal mean initial conditions "//&
871 "for interface depths relative to mean sea level.", &
874 if (just_read)
return
876 call calltree_enter(trim(mdl)//
"(), MOM_state_initialization.F90")
878 call get_param(param_file, mdl,
"INPUTDIR", inputdir, default=
".")
879 filename = trim(slasher(inputdir))//trim(eta_file)
880 call log_param(param_file, mdl,
"INPUTDIR/INTERFACE_IC_FILE", filename)
883 call mom_read_data(filename, eta_var, e0(:), scale=us%m_to_Z)
885 if ((abs(e0(1)) - 0.0) > 0.001)
then
887 do k=nz+1,2,-1 ; e0(k) = e0(k-1) ;
enddo
891 if (e0(2) > e0(1))
then
892 do k=1,nz ; e0(k) = -e0(k) ;
enddo
895 do j=js,je ;
do i=is,ie
901 eta1d(nz+1) = -g%bathyT(i,j)
904 if (eta1d(k) < (eta1d(k+1) + gv%Angstrom_Z))
then
905 eta1d(k) = eta1d(k+1) + gv%Angstrom_Z
906 h(i,j,k) = gv%Angstrom_H
908 h(i,j,k) = gv%Z_to_H * (eta1d(k) - eta1d(k+1))
913 call calltree_leave(trim(mdl)//
'()')
914 end subroutine initialize_thickness_list
917 subroutine initialize_thickness_search
918 call mom_error(fatal,
" MOM_state_initialization.F90, initialize_thickness_search: NOT IMPLEMENTED")
919 end subroutine initialize_thickness_search
922 subroutine convert_thickness(h, G, GV, US, tv)
926 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
932 real,
dimension(SZI_(G),SZJ_(G)) :: &
934 real :: dz_geo(SZI_(G),SZJ_(G))
940 integer,
dimension(2) :: EOSdom
941 integer :: i, j, k, is, ie, js, je, Isq, Ieq, Jsq, Jeq, nz
942 integer :: itt, max_itt
944 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
945 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
948 if (gv%Boussinesq)
then
949 call mom_error(fatal,
"Not yet converting thickness with Boussinesq approx.")
951 i_gearth = gv%RZ_to_H / gv%g_Earth
952 hr_to_pres = gv%g_Earth * gv%H_to_Z
954 if (
associated(tv%eqn_of_state))
then
955 do j=jsq,jeq+1 ;
do i=isq,ieq+1
956 p_bot(i,j) = 0.0 ; p_top(i,j) = 0.0
958 eosdom(:) = eos_domain(g%HI)
961 do i=is,ie ; p_top(i,j) = p_bot(i,j) ;
enddo
963 tv%eqn_of_state, eosdom)
965 p_bot(i,j) = p_top(i,j) + hr_to_pres * (h(i,j,k) * rho(i))
970 call int_specific_vol_dp(tv%T(:,:,k), tv%S(:,:,k), p_top, p_bot, 0.0, g%HI, &
971 tv%eqn_of_state, us, dz_geo)
972 if (itt < max_itt)
then ;
do j=js,je
974 tv%eqn_of_state, eosdom)
978 p_bot(i,j) = p_bot(i,j) + rho(i) * (hr_to_pres*h(i,j,k) - dz_geo(i,j))
983 do j=js,je ;
do i=is,ie
984 h(i,j,k) = (p_bot(i,j) - p_top(i,j)) * i_gearth
988 do k=1,nz ;
do j=js,je ;
do i=is,ie
989 h(i,j,k) = h(i,j,k) * (gv%Rlay(k) / gv%Rho0)
990 enddo ;
enddo ;
enddo
994 end subroutine convert_thickness
997 subroutine depress_surface(h, G, GV, US, param_file, tv, just_read_params)
1001 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
1005 logical,
optional,
intent(in) :: just_read_params
1008 real,
dimension(SZI_(G),SZJ_(G)) :: &
1010 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1) :: &
1013 real :: scale_factor
1015 character(len=40) :: mdl =
"depress_surface"
1016 character(len=200) :: inputdir, eta_srf_file
1017 character(len=200) :: filename, eta_srf_var
1018 logical :: just_read
1019 integer :: i, j, k, is, ie, js, je, nz
1020 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
1022 just_read = .false. ;
if (
present(just_read_params)) just_read = just_read_params
1026 call get_param(param_file, mdl,
"INPUTDIR", inputdir, default=
".")
1027 inputdir = slasher(inputdir)
1028 call get_param(param_file, mdl,
"SURFACE_HEIGHT_IC_FILE", eta_srf_file,&
1029 "The initial condition file for the surface height.", &
1030 fail_if_missing=.not.just_read, do_not_log=just_read)
1031 call get_param(param_file, mdl,
"SURFACE_HEIGHT_IC_VAR", eta_srf_var, &
1032 "The initial condition variable for the surface height.",&
1033 default=
"SSH", do_not_log=just_read)
1034 filename = trim(inputdir)//trim(eta_srf_file)
1035 if (.not.just_read) &
1036 call log_param(param_file, mdl,
"INPUTDIR/SURFACE_HEIGHT_IC_FILE", filename)
1038 call get_param(param_file, mdl,
"SURFACE_HEIGHT_IC_SCALE", scale_factor, &
1039 "A scaling factor to convert SURFACE_HEIGHT_IC_VAR into units of m", &
1040 units=
"variable", default=1.0, scale=us%m_to_Z, do_not_log=just_read)
1042 if (just_read)
return
1044 call mom_read_data(filename, eta_srf_var, eta_sfc, g%Domain, scale=scale_factor)
1047 call find_eta(h, tv, g, gv, us, eta)
1049 do j=js,je ;
do i=is,ie ;
if (g%mask2dT(i,j) > 0.0)
then
1053 if (eta_sfc(i,j) > eta(i,j,1))
then
1055 if (eta_sfc(i,j) - eta(i,j,nz+1) > 10.0*(eta(i,j,1) - eta(i,j,nz+1)))
then
1057 call mom_error(warning,
"Free surface height dilation attempted "//&
1058 "to exceed 10-fold.", all_print=.true.)
1060 dilate = (eta_sfc(i,j) - eta(i,j,nz+1)) / (eta(i,j,1) - eta(i,j,nz+1))
1062 do k=1,nz ; h(i,j,k) = h(i,j,k) * dilate ;
enddo
1063 elseif (eta(i,j,1) > eta_sfc(i,j))
then
1066 if (eta(i,j,k) <= eta_sfc(i,j))
exit
1067 if (eta(i,j,k+1) >= eta_sfc(i,j))
then
1068 h(i,j,k) = gv%Angstrom_H
1070 h(i,j,k) = max(gv%Angstrom_H, h(i,j,k) * &
1071 (eta_sfc(i,j) - eta(i,j,k+1)) / (eta(i,j,k) - eta(i,j,k+1)) )
1075 endif ;
enddo ;
enddo
1077 end subroutine depress_surface
1081 subroutine trim_for_ice(PF, G, GV, US, ALE_CSp, tv, h, just_read_params)
1086 type(
ale_cs),
pointer :: ALE_CSp
1088 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
1090 logical,
optional,
intent(in) :: just_read_params
1093 character(len=200) :: mdl =
"trim_for_ice"
1094 real,
dimension(SZI_(G),SZJ_(G)) :: p_surf
1095 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: S_t, S_b
1096 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: T_t, T_b
1097 character(len=200) :: inputdir, filename, p_surf_file, p_surf_var
1098 real :: scale_factor
1099 real :: min_thickness
1101 logical :: default_2018_answers, remap_answers_2018
1102 logical :: just_read
1103 logical :: use_remapping
1106 just_read = .false. ;
if (
present(just_read_params)) just_read = just_read_params
1108 call get_param(pf, mdl,
"SURFACE_PRESSURE_FILE", p_surf_file, &
1109 "The initial condition file for the surface pressure exerted by ice.", &
1110 fail_if_missing=.not.just_read, do_not_log=just_read)
1111 call get_param(pf, mdl,
"SURFACE_PRESSURE_VAR", p_surf_var, &
1112 "The initial condition variable for the surface pressure exerted by ice.", &
1113 units=
"Pa", default=
"", do_not_log=just_read)
1114 call get_param(pf, mdl,
"INPUTDIR", inputdir, default=
".", do_not_log=.true.)
1115 filename = trim(slasher(inputdir))//trim(p_surf_file)
1116 if (.not.just_read)
call log_param(pf, mdl,
"!INPUTDIR/SURFACE_HEIGHT_IC_FILE", filename)
1118 call get_param(pf, mdl,
"SURFACE_PRESSURE_SCALE", scale_factor, &
1119 "A scaling factor to convert SURFACE_PRESSURE_VAR from "//&
1120 "file SURFACE_PRESSURE_FILE into a surface pressure.", &
1121 units=
"file dependent", default=1., do_not_log=just_read)
1122 call get_param(pf, mdl,
"MIN_THICKNESS", min_thickness,
'Minimum layer thickness', &
1123 units=
'm', default=1.e-3, do_not_log=just_read, scale=us%m_to_Z)
1124 call get_param(pf, mdl,
"TRIMMING_USES_REMAPPING", use_remapping, &
1125 'When trimming the column, also remap T and S.', &
1126 default=.false., do_not_log=just_read)
1127 remap_answers_2018 = .true.
1128 if (use_remapping)
then
1129 call get_param(pf, mdl,
"DEFAULT_2018_ANSWERS", default_2018_answers, &
1130 "This sets the default value for the various _2018_ANSWERS parameters.", &
1132 call get_param(pf, mdl,
"REMAPPING_2018_ANSWERS", remap_answers_2018, &
1133 "If true, use the order of arithmetic and expressions that recover the "//&
1134 "answers from the end of 2018. Otherwise, use updated and more robust "//&
1135 "forms of the same expressions.", default=default_2018_answers)
1138 if (just_read)
return
1140 call mom_read_data(filename, p_surf_var, p_surf, g%Domain, &
1141 scale=scale_factor*us%kg_m3_to_R*us%m_s_to_L_T**2)
1143 if (use_remapping)
then
1145 call initialize_remapping(remap_cs,
'PLM', boundary_extrapolation=.true.)
1149 if (
associated(ale_csp) )
then
1150 call ts_plm_edge_values(ale_csp, s_t, s_b, t_t, t_b, g, gv, tv, h, .true.)
1153 do k=1,g%ke ;
do j=g%jsc,g%jec ;
do i=g%isc,g%iec
1154 t_t(i,j,k) = tv%T(i,j,k) ; t_b(i,j,k) = tv%T(i,j,k)
1155 s_t(i,j,k) = tv%S(i,j,k) ; s_b(i,j,k) = tv%S(i,j,k)
1156 enddo ;
enddo ;
enddo
1159 do j=g%jsc,g%jec ;
do i=g%isc,g%iec
1160 call cut_off_column_top(gv%ke, tv, gv, us, gv%g_Earth, g%bathyT(i,j), &
1161 min_thickness, tv%T(i,j,:), t_t(i,j,:), t_b(i,j,:), &
1162 tv%S(i,j,:), s_t(i,j,:), s_b(i,j,:), p_surf(i,j), h(i,j,:), remap_cs, &
1163 z_tol=1.0e-5*us%m_to_Z, remap_answers_2018=remap_answers_2018)
1166 end subroutine trim_for_ice
1171 subroutine cut_off_column_top(nk, tv, GV, US, G_earth, depth, min_thickness, T, T_t, T_b, &
1172 S, S_t, S_b, p_surf, h, remap_CS, z_tol, remap_answers_2018)
1173 integer,
intent(in) :: nk
1177 real,
intent(in) :: G_earth
1178 real,
intent(in) :: depth
1179 real,
intent(in) :: min_thickness
1180 real,
dimension(nk),
intent(inout) :: T
1181 real,
dimension(nk),
intent(in) :: T_t
1182 real,
dimension(nk),
intent(in) :: T_b
1183 real,
dimension(nk),
intent(inout) :: S
1184 real,
dimension(nk),
intent(in) :: S_t
1185 real,
dimension(nk),
intent(in) :: S_b
1186 real,
intent(in) :: p_surf
1187 real,
dimension(nk),
intent(inout) :: h
1190 real,
optional,
intent(in) :: z_tol
1192 logical,
optional,
intent(in) :: remap_answers_2018
1198 real,
dimension(nk+1) :: e
1199 real,
dimension(nk) :: h0, S0, T0, h1, S1, T1
1201 real :: z_out, e_top
1202 logical :: answers_2018
1205 answers_2018 = .true. ;
if (
present(remap_answers_2018)) answers_2018 = remap_answers_2018
1210 e(k) = e(k+1) + gv%H_to_Z*h(k)
1217 call find_depth_of_pressure_in_cell(t_t(k), t_b(k), s_t(k), s_b(k), e(k), e(k+1), &
1218 p_t, p_surf, gv%Rho0, g_earth, tv%eqn_of_state, &
1219 us, p_b, z_out, z_tol=z_tol)
1220 if (z_out>=e(k))
then
1223 elseif (z_out<=e(k+1))
then
1233 if (e_top<e(1))
then
1236 if (e(k) > e_top)
then
1239 e_top = e_top - min_thickness
1242 h(k) = gv%Z_to_H * max( min_thickness, e(k) - e(k+1) )
1243 if (e(k) < e_top)
exit
1249 if (
associated(remap_cs))
then
1255 if (answers_2018)
then
1256 call remapping_core_h(remap_cs, nk, h0, t0, nk, h1, t1, 1.0e-30*gv%m_to_H, 1.0e-10*gv%m_to_H)
1257 call remapping_core_h(remap_cs, nk, h0, s0, nk, h1, s1, 1.0e-30*gv%m_to_H, 1.0e-10*gv%m_to_H)
1259 call remapping_core_h(remap_cs, nk, h0, t0, nk, h1, t1, gv%H_subroundoff, gv%H_subroundoff)
1260 call remapping_core_h(remap_cs, nk, h0, s0, nk, h1, s1, gv%H_subroundoff, gv%H_subroundoff)
1268 end subroutine cut_off_column_top
1271 subroutine initialize_velocity_from_file(u, v, G, US, param_file, just_read_params)
1273 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
1275 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
1280 logical,
optional,
intent(in) :: just_read_params
1283 character(len=40) :: mdl =
"initialize_velocity_from_file"
1284 character(len=200) :: filename,velocity_file,inputdir
1285 logical :: just_read
1287 just_read = .false. ;
if (
present(just_read_params)) just_read = just_read_params
1289 if (.not.just_read)
call calltree_enter(trim(mdl)//
"(), MOM_state_initialization.F90")
1291 call get_param(param_file, mdl,
"VELOCITY_FILE", velocity_file, &
1292 "The name of the velocity initial condition file.", &
1293 fail_if_missing=.not.just_read, do_not_log=just_read)
1294 call get_param(param_file, mdl,
"INPUTDIR", inputdir, default=
".")
1295 inputdir = slasher(inputdir)
1297 if (just_read)
return
1299 filename = trim(inputdir)//trim(velocity_file)
1300 call log_param(param_file, mdl,
"INPUTDIR/VELOCITY_FILE", filename)
1302 if (.not.
file_exists(filename, g%Domain))
call mom_error(fatal, &
1303 " initialize_velocity_from_file: Unable to open "//trim(filename))
1306 call mom_read_vector(filename,
"u",
"v", u(:,:,:), v(:,:,:), g%Domain, scale=us%m_s_to_L_T)
1308 call calltree_leave(trim(mdl)//
'()')
1309 end subroutine initialize_velocity_from_file
1312 subroutine initialize_velocity_zero(u, v, G, param_file, just_read_params)
1314 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
1316 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
1320 logical,
optional,
intent(in) :: just_read_params
1323 character(len=200) :: mdl =
"initialize_velocity_zero"
1324 logical :: just_read
1325 integer :: i, j, k, is, ie, js, je, Isq, Ieq, Jsq, Jeq, nz
1326 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
1327 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
1329 just_read = .false. ;
if (
present(just_read_params)) just_read = just_read_params
1331 if (.not.just_read)
call calltree_enter(trim(mdl)//
"(), MOM_state_initialization.F90")
1333 if (just_read)
return
1335 do k=1,nz ;
do j=js,je ;
do i=isq,ieq
1337 enddo ;
enddo ;
enddo
1338 do k=1,nz ;
do j=jsq,jeq ;
do i=is,ie
1340 enddo ;
enddo ;
enddo
1342 call calltree_leave(trim(mdl)//
'()')
1343 end subroutine initialize_velocity_zero
1346 subroutine initialize_velocity_uniform(u, v, G, US, param_file, just_read_params)
1348 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
1350 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
1355 logical,
optional,
intent(in) :: just_read_params
1358 integer :: i, j, k, is, ie, js, je, Isq, Ieq, Jsq, Jeq, nz
1359 real :: initial_u_const, initial_v_const
1360 logical :: just_read
1361 character(len=200) :: mdl =
"initialize_velocity_uniform"
1362 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
1363 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
1365 just_read = .false. ;
if (
present(just_read_params)) just_read = just_read_params
1367 call get_param(param_file, mdl,
"INITIAL_U_CONST", initial_u_const, &
1368 "A initial uniform value for the zonal flow.", &
1369 default=0.0, units=
"m s-1", scale=us%m_s_to_L_T, do_not_log=just_read)
1370 call get_param(param_file, mdl,
"INITIAL_V_CONST", initial_v_const, &
1371 "A initial uniform value for the meridional flow.", &
1372 default=0.0, units=
"m s-1", scale=us%m_s_to_L_T, do_not_log=just_read)
1374 if (just_read)
return
1376 do k=1,nz ;
do j=js,je ;
do i=isq,ieq
1377 u(i,j,k) = initial_u_const
1378 enddo ;
enddo ;
enddo
1379 do k=1,nz ;
do j=jsq,jeq ;
do i=is,ie
1380 v(i,j,k) = initial_v_const
1381 enddo ;
enddo ;
enddo
1383 end subroutine initialize_velocity_uniform
1387 subroutine initialize_velocity_circular(u, v, G, US, param_file, just_read_params)
1389 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
1391 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
1396 logical,
optional,
intent(in) :: just_read_params
1399 character(len=200) :: mdl =
"initialize_velocity_circular"
1400 real :: circular_max_u
1403 logical :: just_read
1404 integer :: i, j, k, is, ie, js, je, Isq, Ieq, Jsq, Jeq, nz
1405 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
1406 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
1408 just_read = .false. ;
if (
present(just_read_params)) just_read = just_read_params
1410 call get_param(param_file, mdl,
"CIRCULAR_MAX_U", circular_max_u, &
1411 "The amplitude of zonal flow from which to scale the "// &
1412 "circular stream function [m s-1].", &
1413 units=
"m s-1", default=0., scale=us%m_s_to_L_T, do_not_log=just_read)
1415 if (just_read)
return
1419 do k=1,nz ;
do j=js,je ;
do i=isq,ieq
1421 psi2 = my_psi(i,j-1)
1422 u(i,j,k) = (psi1 - psi2) / g%dy_Cu(i,j)
1423 enddo ;
enddo ;
enddo
1424 do k=1,nz ;
do j=jsq,jeq ;
do i=is,ie
1426 psi2 = my_psi(i-1,j)
1427 v(i,j,k) = (psi2 - psi1) / g%dx_Cv(i,j)
1428 enddo ;
enddo ;
enddo
1433 real function my_psi(ig,jg)
1439 x = 2.0*(g%geoLonBu(ig,jg)-g%west_lon) / g%len_lon - 1.0
1440 y = 2.0*(g%geoLatBu(ig,jg)-g%south_lat) / g%len_lat - 1.0
1441 r = sqrt( x**2 + y**2 )
1443 my_psi = 0.5*(1.0 - cos(dpi*r))
1444 my_psi = my_psi * (circular_max_u * g%US%m_to_L*g%len_lon*1e3 / dpi)
1447 end subroutine initialize_velocity_circular
1450 subroutine initialize_temp_salt_from_file(T, S, G, param_file, just_read_params)
1451 type(ocean_grid_type),
intent(in) :: G
1452 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(out) :: T
1454 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(out) :: S
1456 type(param_file_type),
intent(in) :: param_file
1457 logical,
optional,
intent(in) :: just_read_params
1460 logical :: just_read
1461 character(len=200) :: filename, salt_filename
1462 character(len=200) :: ts_file, salt_file, inputdir
1463 character(len=40) :: mdl =
"initialize_temp_salt_from_file"
1464 character(len=64) :: temp_var, salt_var
1466 just_read = .false. ;
if (
present(just_read_params)) just_read = just_read_params
1468 if (.not.just_read)
call calltree_enter(trim(mdl)//
"(), MOM_state_initialization.F90")
1470 call get_param(param_file, mdl,
"TS_FILE", ts_file, &
1471 "The initial condition file for temperature.", &
1472 fail_if_missing=.not.just_read, do_not_log=just_read)
1473 call get_param(param_file, mdl,
"INPUTDIR", inputdir, default=
".")
1474 inputdir = slasher(inputdir)
1476 filename = trim(inputdir)//trim(ts_file)
1477 if (.not.just_read)
call log_param(param_file, mdl,
"INPUTDIR/TS_FILE", filename)
1478 call get_param(param_file, mdl,
"TEMP_IC_VAR", temp_var, &
1479 "The initial condition variable for potential temperature.", &
1480 default=
"PTEMP", do_not_log=just_read)
1481 call get_param(param_file, mdl,
"SALT_IC_VAR", salt_var, &
1482 "The initial condition variable for salinity.", &
1483 default=
"SALT", do_not_log=just_read)
1484 call get_param(param_file, mdl,
"SALT_FILE", salt_file, &
1485 "The initial condition file for salinity.", &
1486 default=trim(ts_file), do_not_log=just_read)
1488 if (just_read)
return
1490 if (.not.file_exists(filename, g%Domain))
call mom_error(fatal, &
1491 " initialize_temp_salt_from_file: Unable to open "//trim(filename))
1494 call mom_read_data(filename, temp_var, t(:,:,:), g%Domain)
1496 salt_filename = trim(inputdir)//trim(salt_file)
1497 if (.not.file_exists(salt_filename, g%Domain))
call mom_error(fatal, &
1498 " initialize_temp_salt_from_file: Unable to open "//trim(salt_filename))
1500 call mom_read_data(salt_filename, salt_var, s(:,:,:), g%Domain)
1502 call calltree_leave(trim(mdl)//
'()')
1503 end subroutine initialize_temp_salt_from_file
1506 subroutine initialize_temp_salt_from_profile(T, S, G, param_file, just_read_params)
1507 type(ocean_grid_type),
intent(in) :: G
1508 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(out) :: T
1510 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(out) :: S
1512 type(param_file_type),
intent(in) :: param_file
1513 logical,
optional,
intent(in) :: just_read_params
1516 real,
dimension(SZK_(G)) :: T0, S0
1518 logical :: just_read
1519 character(len=200) :: filename, ts_file, inputdir
1520 character(len=40) :: mdl =
"initialize_temp_salt_from_profile"
1522 just_read = .false. ;
if (
present(just_read_params)) just_read = just_read_params
1524 if (.not.just_read)
call calltree_enter(trim(mdl)//
"(), MOM_state_initialization.F90")
1526 call get_param(param_file, mdl,
"TS_FILE", ts_file, &
1527 "The file with the reference profiles for temperature "//&
1528 "and salinity.", fail_if_missing=.not.just_read, do_not_log=just_read)
1530 if (just_read)
return
1532 call get_param(param_file, mdl,
"INPUTDIR", inputdir, default=
".")
1533 inputdir = slasher(inputdir)
1534 filename = trim(inputdir)//trim(ts_file)
1535 call log_param(param_file, mdl,
"INPUTDIR/TS_FILE", filename)
1536 if (.not.file_exists(filename))
call mom_error(fatal, &
1537 " initialize_temp_salt_from_profile: Unable to open "//trim(filename))
1540 call mom_read_data(filename,
"PTEMP", t0(:))
1541 call mom_read_data(filename,
"SALT", s0(:))
1543 do k=1,g%ke ;
do j=g%jsc,g%jec ;
do i=g%isc,g%iec
1544 t(i,j,k) = t0(k) ; s(i,j,k) = s0(k)
1545 enddo ;
enddo ;
enddo
1547 call calltree_leave(trim(mdl)//
'()')
1548 end subroutine initialize_temp_salt_from_profile
1551 subroutine initialize_temp_salt_fit(T, S, G, GV, US, param_file, eqn_of_state, P_Ref, just_read_params)
1552 type(ocean_grid_type),
intent(in) :: G
1553 type(verticalgrid_type),
intent(in) :: GV
1554 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(out) :: T
1556 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(out) :: S
1558 type(unit_scale_type),
intent(in) :: US
1559 type(param_file_type),
intent(in) :: param_file
1561 type(eos_type),
pointer :: eqn_of_state
1562 real,
intent(in) :: P_Ref
1564 logical,
optional,
intent(in) :: just_read_params
1571 real :: pres(SZK_(G))
1572 real :: drho_dT(SZK_(G))
1573 real :: drho_dS(SZK_(G))
1574 real :: rho_guess(SZK_(G))
1575 logical :: fit_salin
1576 logical :: just_read
1577 character(len=40) :: mdl =
"initialize_temp_salt_fit"
1578 integer :: i, j, k, itt, nz
1581 just_read = .false. ;
if (
present(just_read_params)) just_read = just_read_params
1583 if (.not.just_read)
call calltree_enter(trim(mdl)//
"(), MOM_state_initialization.F90")
1585 call get_param(param_file, mdl,
"T_REF", t_ref, &
1586 "A reference temperature used in initialization.", &
1587 units=
"degC", fail_if_missing=.not.just_read, do_not_log=just_read)
1588 call get_param(param_file, mdl,
"S_REF", s_ref, &
1589 "A reference salinity used in initialization.", units=
"PSU", &
1590 default=35.0, do_not_log=just_read)
1591 call get_param(param_file, mdl,
"FIT_SALINITY", fit_salin, &
1592 "If true, accept the prescribed temperature and fit the "//&
1593 "salinity; otherwise take salinity and fit temperature.", &
1594 default=.false., do_not_log=just_read)
1596 if (just_read)
return
1599 pres(k) = p_ref ; s0(k) = s_ref
1603 call calculate_density(t0(1), s0(1), pres(1), rho_guess(1), eqn_of_state)
1604 call calculate_density_derivs(t0, s0, pres, drho_dt, drho_ds, eqn_of_state, (/1,1/) )
1609 s0(k) = max(0.0, s0(1) + (gv%Rlay(k) - rho_guess(1)) / drho_ds(1))
1613 call calculate_density(t0, s0, pres, rho_guess, eqn_of_state)
1614 call calculate_density_derivs(t0, s0, pres, drho_dt, drho_ds, eqn_of_state)
1616 s0(k) = max(0.0, s0(k) + (gv%Rlay(k) - rho_guess(k)) / drho_ds(k))
1622 t0(k) = t0(1) + (gv%Rlay(k) - rho_guess(1)) / drho_dt(1)
1625 call calculate_density(t0, s0, pres, rho_guess, eqn_of_state)
1626 call calculate_density_derivs(t0, s0, pres, drho_dt, drho_ds, eqn_of_state)
1628 t0(k) = t0(k) + (gv%Rlay(k) - rho_guess(k)) / drho_dt(k)
1633 do k=1,nz ;
do j=g%jsd,g%jed ;
do i=g%isd,g%ied
1634 t(i,j,k) = t0(k) ; s(i,j,k) = s0(k)
1635 enddo ;
enddo ;
enddo
1637 call calltree_leave(trim(mdl)//
'()')
1638 end subroutine initialize_temp_salt_fit
1645 subroutine initialize_temp_salt_linear(T, S, G, param_file, just_read_params)
1646 type(ocean_grid_type),
intent(in) :: G
1647 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(out) :: T
1649 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(out) :: S
1651 type(param_file_type),
intent(in) :: param_file
1653 logical,
optional,
intent(in) :: just_read_params
1659 real :: delta_S, delta_T
1660 real :: S_top, T_top
1661 real :: S_range, T_range
1663 logical :: just_read
1664 character(len=40) :: mdl =
"initialize_temp_salt_linear"
1666 just_read = .false. ;
if (
present(just_read_params)) just_read = just_read_params
1668 if (.not.just_read)
call calltree_enter(trim(mdl)//
"(), MOM_state_initialization.F90")
1669 call get_param(param_file, mdl,
"T_TOP", t_top, &
1670 "Initial temperature of the top surface.", &
1671 units=
"degC", fail_if_missing=.not.just_read, do_not_log=just_read)
1672 call get_param(param_file, mdl,
"T_RANGE", t_range, &
1673 "Initial temperature difference (top-bottom).", &
1674 units=
"degC", fail_if_missing=.not.just_read, do_not_log=just_read)
1675 call get_param(param_file, mdl,
"S_TOP", s_top, &
1676 "Initial salinity of the top surface.", &
1677 units=
"PSU", fail_if_missing=.not.just_read, do_not_log=just_read)
1678 call get_param(param_file, mdl,
"S_RANGE", s_range, &
1679 "Initial salinity difference (top-bottom).", &
1680 units=
"PSU", fail_if_missing=.not.just_read, do_not_log=just_read)
1682 if (just_read)
return
1691 s(:,:,k) = s_top - s_range*((real(k)-0.5)/real(g%ke))
1692 t(:,:,k) = t_top - t_range*((real(k)-0.5)/real(g%ke))
1704 call calltree_leave(trim(mdl)//
'()')
1705 end subroutine initialize_temp_salt_linear
1712 subroutine initialize_sponges_file(G, GV, US, use_temperature, tv, param_file, Layer_CSp, ALE_CSp, Time)
1713 type(ocean_grid_type),
intent(in) :: G
1714 type(verticalgrid_type),
intent(in) :: GV
1715 type(unit_scale_type),
intent(in) :: US
1716 logical,
intent(in) :: use_temperature
1717 type(thermo_var_ptrs),
intent(in) :: tv
1719 type(param_file_type),
intent(in) :: param_file
1720 type(sponge_cs),
pointer :: Layer_CSp
1722 type(ale_sponge_cs),
pointer :: ALE_CSp
1724 type(time_type),
intent(in) :: Time
1727 real,
allocatable,
dimension(:,:,:) :: eta
1728 real,
allocatable,
dimension(:,:,:) :: h
1730 real,
dimension (SZI_(G),SZJ_(G),SZK_(G)) :: &
1732 real,
dimension (SZI_(G),SZJ_(G)) :: &
1734 real,
allocatable,
dimension(:,:,:) :: tmp_tr
1736 real :: Idamp(SZI_(G),SZJ_(G))
1737 real :: pres(SZI_(G))
1739 integer,
dimension(2) :: EOSdom
1740 integer :: i, j, k, is, ie, js, je, nz
1741 integer :: isd, ied, jsd, jed
1742 integer,
dimension(4) :: siz
1744 character(len=40) :: potemp_var, salin_var, Idamp_var, eta_var
1745 character(len=40) :: mdl =
"initialize_sponges_file"
1746 character(len=200) :: damping_file, state_file
1747 character(len=200) :: filename, inputdir
1750 logical :: time_space_interp_sponge
1755 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
1756 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
1758 pres(:) = 0.0 ; tmp(:,:,:) = 0.0 ; idamp(:,:) = 0.0
1760 call get_param(param_file, mdl,
"INPUTDIR", inputdir, default=
".")
1761 inputdir = slasher(inputdir)
1762 call get_param(param_file, mdl,
"SPONGE_DAMPING_FILE", damping_file, &
1763 "The name of the file with the sponge damping rates.", &
1764 fail_if_missing=.true.)
1765 call get_param(param_file, mdl,
"SPONGE_STATE_FILE", state_file, &
1766 "The name of the file with the state to damp toward.", &
1767 default=damping_file)
1768 call get_param(param_file, mdl,
"SPONGE_PTEMP_VAR", potemp_var, &
1769 "The name of the potential temperature variable in "//&
1770 "SPONGE_STATE_FILE.", default=
"PTEMP")
1771 call get_param(param_file, mdl,
"SPONGE_SALT_VAR", salin_var, &
1772 "The name of the salinity variable in "//&
1773 "SPONGE_STATE_FILE.", default=
"SALT")
1774 call get_param(param_file, mdl,
"SPONGE_ETA_VAR", eta_var, &
1775 "The name of the interface height variable in "//&
1776 "SPONGE_STATE_FILE.", default=
"ETA")
1777 call get_param(param_file, mdl,
"SPONGE_IDAMP_VAR", idamp_var, &
1778 "The name of the inverse damping rate variable in "//&
1779 "SPONGE_DAMPING_FILE.", default=
"Idamp")
1780 call get_param(param_file, mdl,
"USE_REGRIDDING", use_ale, do_not_log = .true.)
1781 time_space_interp_sponge = .false.
1782 call get_param(param_file, mdl,
"NEW_SPONGES", time_space_interp_sponge, &
1783 "Set True if using the newer sponging code which "//&
1784 "performs on-the-fly regridding in lat-lon-time.",&
1785 "of sponge restoring data.", default=.false.)
1786 if (time_space_interp_sponge)
then
1787 call mom_error(warning,
" initialize_sponges: NEW_SPONGES has been deprecated. "//&
1788 "Please use INTERPOLATE_SPONGE_TIME_SPACE instead. Setting "//&
1789 "INTERPOLATE_SPONGE_TIME_SPACE = True.")
1791 call get_param(param_file, mdl,
"INTERPOLATE_SPONGE_TIME_SPACE", time_space_interp_sponge, &
1792 "Set True if using the newer sponging code which "//&
1793 "performs on-the-fly regridding in lat-lon-time.",&
1794 "of sponge restoring data.", default=time_space_interp_sponge)
1796 filename = trim(inputdir)//trim(damping_file)
1797 call log_param(param_file, mdl,
"INPUTDIR/SPONGE_DAMPING_FILE", filename)
1798 if (.not.file_exists(filename, g%Domain)) &
1799 call mom_error(fatal,
" initialize_sponges: Unable to open "//trim(filename))
1801 if (time_space_interp_sponge .and. .not. use_ale) &
1802 call mom_error(fatal,
" initialize_sponges: Time-varying sponges are currently unavailable in layered mode ")
1804 call mom_read_data(filename, idamp_var, idamp(:,:), g%Domain, scale=us%T_to_s)
1810 filename = trim(inputdir)//trim(state_file)
1811 call log_param(param_file, mdl,
"INPUTDIR/SPONGE_STATE_FILE", filename)
1812 if (.not.file_exists(filename, g%Domain)) &
1813 call mom_error(fatal,
" initialize_sponges: Unable to open "//trim(filename))
1816 if (.not. use_ale)
then
1818 allocate(eta(isd:ied,jsd:jed,nz+1)); eta(:,:,:) = 0.0
1819 call mom_read_data(filename, eta_var, eta(:,:,:), g%Domain, scale=us%m_to_Z)
1821 do j=js,je ;
do i=is,ie
1822 eta(i,j,nz+1) = -g%bathyT(i,j)
1824 do k=nz,1,-1 ;
do j=js,je ;
do i=is,ie
1825 if (eta(i,j,k) < (eta(i,j,k+1) + gv%Angstrom_Z)) &
1826 eta(i,j,k) = eta(i,j,k+1) + gv%Angstrom_Z
1827 enddo ;
enddo ;
enddo
1830 call initialize_sponge(idamp, eta, g, param_file, layer_csp, gv)
1833 if ( gv%nkml>0)
then
1837 do i=is-1,ie ; pres(i) = tv%P_Ref ;
enddo
1838 eosdom(:) = eos_domain(g%HI)
1840 call mom_read_data(filename, potemp_var, tmp(:,:,:), g%Domain)
1841 call mom_read_data(filename, salin_var, tmp2(:,:,:), g%Domain)
1844 call calculate_density(tmp(:,j,1), tmp2(:,j,1), pres, tmp_2d(:,j), tv%eqn_of_state, eosdom)
1847 call set_up_sponge_ml_density(tmp_2d, g, layer_csp)
1856 if ( use_temperature)
then
1857 call mom_read_data(filename, potemp_var, tmp(:,:,:), g%Domain)
1858 call set_up_sponge_field(tmp, tv%T, g, nz, layer_csp)
1859 call mom_read_data(filename, salin_var, tmp(:,:,:), g%Domain)
1860 call set_up_sponge_field(tmp, tv%S, g, nz, layer_csp)
1867 if (.not. time_space_interp_sponge)
then
1868 call field_size(filename,eta_var,siz,no_domain=.true.)
1869 if (siz(1) /= g%ieg-g%isg+1 .or. siz(2) /= g%jeg-g%jsg+1) &
1870 call mom_error(fatal,
"initialize_sponge_file: Array size mismatch for sponge data.")
1872 allocate(eta(isd:ied,jsd:jed,nz_data+1))
1873 allocate(h(isd:ied,jsd:jed,nz_data))
1874 call mom_read_data(filename, eta_var, eta(:,:,:), g%Domain, scale=us%m_to_Z)
1875 do j=js,je ;
do i=is,ie
1876 eta(i,j,nz+1) = -g%bathyT(i,j)
1878 do k=nz,1,-1 ;
do j=js,je ;
do i=is,ie
1879 if (eta(i,j,k) < (eta(i,j,k+1) + gv%Angstrom_Z)) &
1880 eta(i,j,k) = eta(i,j,k+1) + gv%Angstrom_Z
1881 enddo ;
enddo ;
enddo
1882 do k=1,nz;
do j=js,je ;
do i=is,ie
1883 h(i,j,k) = gv%Z_to_H*(eta(i,j,k)-eta(i,j,k+1))
1884 enddo ;
enddo ;
enddo
1885 call initialize_ale_sponge(idamp, g, param_file, ale_csp, h, nz_data)
1888 if (use_temperature)
then
1889 allocate(tmp_tr(isd:ied,jsd:jed,nz_data))
1890 call mom_read_data(filename, potemp_var, tmp_tr(:,:,:), g%Domain)
1891 call set_up_ale_sponge_field(tmp_tr, g, tv%T, ale_csp)
1892 call mom_read_data(filename, salin_var, tmp_tr(:,:,:), g%Domain)
1893 call set_up_ale_sponge_field(tmp_tr, g, tv%S, ale_csp)
1898 call initialize_ale_sponge(idamp, g, param_file, ale_csp)
1900 if ( use_temperature)
then
1901 call set_up_ale_sponge_field(filename, potemp_var, time, g, gv, us, tv%T, ale_csp)
1902 call set_up_ale_sponge_field(filename, salin_var, time, g, gv, us, tv%S, ale_csp)
1911 end subroutine initialize_sponges_file
1915 subroutine set_velocity_depth_max(G)
1916 type(ocean_grid_type),
intent(inout) :: G
1920 do i=g%isd,g%ied-1 ;
do j=g%jsd,g%jed
1921 g%Dblock_u(i,j) = g%mask2dCu(i,j) * max(g%bathyT(i,j), g%bathyT(i+1,j))
1922 g%Dopen_u(i,j) = g%Dblock_u(i,j)
1924 do i=g%isd,g%ied ;
do j=g%jsd,g%jed-1
1925 g%Dblock_v(i,j) = g%mask2dCv(i,j) * max(g%bathyT(i,j), g%bathyT(i,j+1))
1926 g%Dopen_v(i,j) = g%Dblock_v(i,j)
1928 end subroutine set_velocity_depth_max
1932 subroutine compute_global_grid_integrals(G, US)
1933 type(ocean_grid_type),
intent(inout) :: G
1934 type(unit_scale_type),
intent(in) :: US
1936 real,
dimension(G%isc:G%iec, G%jsc:G%jec) :: tmpForSumming
1940 area_scale = us%L_to_m**2
1941 tmpforsumming(:,:) = 0.
1942 g%areaT_global = 0.0 ; g%IareaT_global = 0.0
1943 do j=g%jsc,g%jec ;
do i=g%isc,g%iec
1944 tmpforsumming(i,j) = area_scale*g%areaT(i,j) * g%mask2dT(i,j)
1946 g%areaT_global = reproducing_sum(tmpforsumming)
1947 g%IareaT_global = 1. / (g%areaT_global)
1948 end subroutine compute_global_grid_integrals
1952 subroutine set_velocity_depth_min(G)
1953 type(ocean_grid_type),
intent(inout) :: G
1957 do i=g%isd,g%ied-1 ;
do j=g%jsd,g%jed
1958 g%Dblock_u(i,j) = g%mask2dCu(i,j) * min(g%bathyT(i,j), g%bathyT(i+1,j))
1959 g%Dopen_u(i,j) = g%Dblock_u(i,j)
1961 do i=g%isd,g%ied ;
do j=g%jsd,g%jed-1
1962 g%Dblock_v(i,j) = g%mask2dCv(i,j) * min(g%bathyT(i,j), g%bathyT(i,j+1))
1963 g%Dopen_v(i,j) = g%Dblock_v(i,j)
1965 end subroutine set_velocity_depth_min
1970 subroutine mom_temp_salt_initialize_from_z(h, tv, G, GV, US, PF, just_read_params)
1971 type(ocean_grid_type),
intent(inout) :: G
1972 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
1974 type(thermo_var_ptrs),
intent(inout) :: tv
1976 type(verticalgrid_type),
intent(in) :: GV
1977 type(unit_scale_type),
intent(in) :: US
1978 type(param_file_type),
intent(in) :: PF
1980 logical,
optional,
intent(in) :: just_read_params
1984 character(len=200) :: filename
1986 character(len=200) :: tfilename
1988 character(len=200) :: sfilename
1990 character(len=200) :: shelf_file
1991 character(len=200) :: inputdir
1992 character(len=200) :: mesg, area_varname, ice_shelf_file
1994 type(eos_type),
pointer :: eos => null()
1995 type(thermo_var_ptrs) :: tv_loc
1996 type(verticalgrid_type) :: GV_loc
1998 # include "version_variable.h"
1999 character(len=40) :: mdl =
"MOM_initialize_layers_from_Z"
2001 integer,
dimension(2) :: EOSdom
2002 integer :: is, ie, js, je, nz
2003 integer :: isc,iec,jsc,jec
2004 integer :: isg, ieg, jsg, jeg
2005 integer :: isd, ied, jsd, jed
2007 integer :: i, j, k, ks, np, ni, nj
2010 integer :: kd, inconsistent
2016 real,
dimension(:,:),
pointer :: shelf_area => null()
2017 real :: Hmix_default
2020 real :: missing_value_temp, missing_value_salt
2021 logical :: correct_thickness
2022 character(len=40) :: potemp_var, salin_var
2023 character(len=8) :: laynum
2025 integer,
parameter :: niter=10
2026 logical :: just_read
2027 logical :: adjust_temperature = .true.
2028 real,
parameter :: missing_value = -1.e20
2029 real,
parameter :: temp_land_fill = 0.0, salt_land_fill = 35.0
2030 logical :: reentrant_x, tripolar_n,dbg
2031 logical :: debug = .false.
2034 real,
dimension(:),
allocatable :: z_edges_in, z_in
2035 real,
dimension(:),
allocatable :: Rb
2036 real,
dimension(:,:,:),
allocatable,
target :: temp_z, salt_z, mask_z
2037 real,
dimension(:,:,:),
allocatable :: rho_z
2038 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1) :: zi
2039 integer,
dimension(SZI_(G),SZJ_(G)) :: nlevs
2040 real,
dimension(SZI_(G)) :: press
2043 real,
dimension(:),
allocatable :: hTarget
2044 real,
dimension(:,:),
allocatable :: area_shelf_h
2045 real,
dimension(:,:),
allocatable,
target :: frac_shelf_h
2046 real,
dimension(:,:,:),
allocatable,
target :: tmpT1dIn, tmpS1dIn
2047 real,
dimension(:,:,:),
allocatable :: tmp_mask_in
2048 real,
dimension(:,:,:),
allocatable :: h1
2049 real,
dimension(:,:,:),
allocatable :: dz_interface
2050 real :: zTopOfCell, zBottomOfCell
2051 type(regridding_cs) :: regridCS
2052 type(remapping_cs) :: remapCS
2054 logical :: homogenize, useALEremapping, remap_full_column, remap_general, remap_old_alg
2055 logical :: answers_2018, default_2018_answers, hor_regrid_answers_2018
2056 logical :: use_ice_shelf
2057 logical :: pre_gridded
2058 logical :: separate_mixed_layer
2059 character(len=10) :: remappingScheme
2060 real :: tempAvg, saltAvg
2061 integer :: nPoints, ans
2062 integer :: id_clock_routine, id_clock_read, id_clock_interp, id_clock_fill, id_clock_ALE
2064 id_clock_routine = cpu_clock_id(
'(Initialize from Z)', grain=clock_routine)
2065 id_clock_ale = cpu_clock_id(
'(Initialize from Z) ALE', grain=clock_loop)
2067 call cpu_clock_begin(id_clock_routine)
2069 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
2070 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
2071 isg = g%isg ; ieg = g%ieg ; jsg = g%jsg ; jeg = g%jeg
2073 pi_180=atan(1.0)/45.
2075 just_read = .false. ;
if (
present(just_read_params)) just_read = just_read_params
2077 if (.not.just_read)
call calltree_enter(trim(mdl)//
"(), MOM_state_initialization.F90")
2078 if (.not.just_read)
call log_version(pf, mdl, version,
"")
2080 inputdir =
"." ;
call get_param(pf, mdl,
"INPUTDIR", inputdir)
2081 inputdir = slasher(inputdir)
2083 eos => tv%eqn_of_state
2085 reentrant_x = .false. ;
call get_param(pf, mdl,
"REENTRANT_X", reentrant_x, default=.true.)
2086 tripolar_n = .false. ;
call get_param(pf, mdl,
"TRIPOLAR_N", tripolar_n, default=.false.)
2088 call get_param(pf, mdl,
"TEMP_SALT_Z_INIT_FILE",filename, &
2089 "The name of the z-space input file used to initialize "//&
2090 "temperatures (T) and salinities (S). If T and S are not "//&
2091 "in the same file, TEMP_Z_INIT_FILE and SALT_Z_INIT_FILE "//&
2092 "must be set.",default=
"temp_salt_z.nc",do_not_log=just_read)
2093 call get_param(pf, mdl,
"TEMP_Z_INIT_FILE",tfilename, &
2094 "The name of the z-space input file used to initialize "//&
2095 "temperatures, only.", default=trim(filename),do_not_log=just_read)
2096 call get_param(pf, mdl,
"SALT_Z_INIT_FILE",sfilename, &
2097 "The name of the z-space input file used to initialize "//&
2098 "temperatures, only.", default=trim(filename),do_not_log=just_read)
2099 filename = trim(inputdir)//trim(filename)
2100 tfilename = trim(inputdir)//trim(tfilename)
2101 sfilename = trim(inputdir)//trim(sfilename)
2102 call get_param(pf, mdl,
"Z_INIT_FILE_PTEMP_VAR", potemp_var, &
2103 "The name of the potential temperature variable in "//&
2104 "TEMP_Z_INIT_FILE.", default=
"ptemp",do_not_log=just_read)
2105 call get_param(pf, mdl,
"Z_INIT_FILE_SALT_VAR", salin_var, &
2106 "The name of the salinity variable in "//&
2107 "SALT_Z_INIT_FILE.", default=
"salt",do_not_log=just_read)
2108 call get_param(pf, mdl,
"Z_INIT_HOMOGENIZE", homogenize, &
2109 "If True, then horizontally homogenize the interpolated "//&
2110 "initial conditions.", default=.false., do_not_log=just_read)
2111 call get_param(pf, mdl,
"Z_INIT_ALE_REMAPPING", usealeremapping, &
2112 "If True, then remap straight to model coordinate from file.", &
2113 default=.false., do_not_log=just_read)
2114 call get_param(pf, mdl,
"Z_INIT_REMAPPING_SCHEME", remappingscheme, &
2115 "The remapping scheme to use if using Z_INIT_ALE_REMAPPING "//&
2116 "is True.", default=
"PPM_IH4", do_not_log=just_read)
2117 call get_param(pf, mdl,
"Z_INIT_REMAP_GENERAL", remap_general, &
2118 "If false, only initializes to z* coordinates. "//&
2119 "If true, allows initialization directly to general coordinates.",&
2120 default=.false., do_not_log=just_read)
2121 call get_param(pf, mdl,
"Z_INIT_REMAP_FULL_COLUMN", remap_full_column, &
2122 "If false, only reconstructs profiles for valid data points. "//&
2123 "If true, inserts vanished layers below the valid data.", &
2124 default=remap_general, do_not_log=just_read)
2125 call get_param(pf, mdl,
"Z_INIT_REMAP_OLD_ALG", remap_old_alg, &
2126 "If false, uses the preferred remapping algorithm for initialization. "//&
2127 "If true, use an older, less robust algorithm for remapping.", &
2128 default=.false., do_not_log=just_read)
2129 call get_param(pf, mdl,
"DEFAULT_2018_ANSWERS", default_2018_answers, &
2130 "This sets the default value for the various _2018_ANSWERS parameters.", &
2132 call get_param(pf, mdl,
"TEMP_SALT_INIT_VERTICAL_REMAP_ONLY", pre_gridded, &
2133 "If true, initial conditions are on the model horizontal grid. " //&
2134 "Extrapolation over missing ocean values is done using an ICE-9 "//&
2135 "procedure with vertical ALE remapping .", &
2137 if (usealeremapping)
then
2138 call get_param(pf, mdl,
"REMAPPING_2018_ANSWERS", answers_2018, &
2139 "If true, use the order of arithmetic and expressions that recover the "//&
2140 "answers from the end of 2018. Otherwise, use updated and more robust "//&
2141 "forms of the same expressions.", default=default_2018_answers)
2143 call get_param(pf, mdl,
"HOR_REGRID_2018_ANSWERS", hor_regrid_answers_2018, &
2144 "If true, use the order of arithmetic for horizonal regridding that recovers "//&
2145 "the answers from the end of 2018. Otherwise, use rotationally symmetric "//&
2146 "forms of the same expressions.", default=default_2018_answers)
2147 call get_param(pf, mdl,
"ICE_SHELF", use_ice_shelf, default=.false.)
2148 if (use_ice_shelf)
then
2149 call get_param(pf, mdl,
"ICE_THICKNESS_FILE", ice_shelf_file, &
2150 "The file from which the ice bathymetry and area are read.", &
2151 fail_if_missing=.not.just_read, do_not_log=just_read)
2152 shelf_file = trim(inputdir)//trim(ice_shelf_file)
2153 if (.not.just_read)
call log_param(pf, mdl,
"INPUTDIR/THICKNESS_FILE", shelf_file)
2154 call get_param(pf, mdl,
"ICE_AREA_VARNAME", area_varname, &
2155 "The name of the area variable in ICE_THICKNESS_FILE.", &
2156 fail_if_missing=.not.just_read, do_not_log=just_read)
2158 if (.not.usealeremapping)
then
2159 call get_param(pf, mdl,
"ADJUST_THICKNESS", correct_thickness, &
2160 "If true, all mass below the bottom removed if the "//&
2161 "topography is shallower than the thickness input file "//&
2162 "would indicate.", default=.false., do_not_log=just_read)
2164 call get_param(pf, mdl,
"FIT_TO_TARGET_DENSITY_IC", adjust_temperature, &
2165 "If true, all the interior layers are adjusted to "//&
2166 "their target densities using mostly temperature "//&
2167 "This approach can be problematic, particularly in the "//&
2168 "high latitudes.", default=.true., do_not_log=just_read)
2169 call get_param(pf, mdl,
"Z_INIT_SEPARATE_MIXED_LAYER", separate_mixed_layer, &
2170 "If true, distribute the topmost Z_INIT_HMIX_DEPTH of water over NKML layers, "//&
2171 "and do not correct the density of the topmost NKML+NKBL layers. Otherwise "//&
2172 "all layers are initialized based on the depths of their target densities.", &
2173 default=.false., do_not_log=just_read.or.(gv%nkml==0))
2174 if (gv%nkml == 0) separate_mixed_layer = .false.
2175 call get_param(pf, mdl,
"MINIMUM_DEPTH", hmix_default, default=0.0)
2176 call get_param(pf, mdl,
"Z_INIT_HMIX_DEPTH", hmix_depth, &
2177 "The mixed layer depth in the initial conditions when Z_INIT_SEPARATE_MIXED_LAYER "//&
2178 "is set to true.", default=hmix_default, units=
"m", scale=us%m_to_Z, &
2179 do_not_log=(just_read .or. .not.separate_mixed_layer))
2184 call cpu_clock_end(id_clock_routine)
2188 eps_z = gv%Angstrom_Z
2189 eps_rho = 1.0e-10*us%kg_m3_to_R
2206 call horiz_interp_and_extrap_tracer(tfilename, potemp_var, 1.0, 1, &
2207 g, temp_z, mask_z, z_in, z_edges_in, missing_value_temp, reentrant_x, &
2208 tripolar_n, homogenize, m_to_z=us%m_to_Z, answers_2018=hor_regrid_answers_2018, ongrid=pre_gridded)
2210 call horiz_interp_and_extrap_tracer(sfilename, salin_var, 1.0, 1, &
2211 g, salt_z, mask_z, z_in, z_edges_in, missing_value_salt, reentrant_x, &
2212 tripolar_n, homogenize, m_to_z=us%m_to_Z, answers_2018=hor_regrid_answers_2018, ongrid=pre_gridded)
2217 do k=1,
size(z_edges_in,1) ; z_edges_in(k) = -z_edges_in(k) ;
enddo
2219 allocate(rho_z(isd:ied,jsd:jed,kd))
2220 allocate(area_shelf_h(isd:ied,jsd:jed))
2221 allocate(frac_shelf_h(isd:ied,jsd:jed))
2224 call convert_temp_salt_for_teos10(temp_z, salt_z, g%HI, kd, mask_z, eos)
2227 eosdom(:) = eos_domain(g%HI)
2228 do k=1,kd ;
do j=js,je
2229 call calculate_density(temp_z(:,j,k), salt_z(:,j,k), press, rho_z(:,j,k), eos, eosdom)
2232 call pass_var(temp_z,g%Domain)
2233 call pass_var(salt_z,g%Domain)
2234 call pass_var(mask_z,g%Domain)
2235 call pass_var(rho_z,g%Domain)
2238 if (use_ice_shelf)
then
2239 if (.not.file_exists(shelf_file, g%Domain))
call mom_error(fatal, &
2240 "MOM_temp_salt_initialize_from_Z: Unable to open shelf file "//trim(shelf_file))
2242 call mom_read_data(shelf_file, trim(area_varname), area_shelf_h, g%Domain, scale=us%m_to_L**2)
2245 frac_shelf_h(:,:) = 0.0
2247 do j=jsd,jed ;
do i=isd,ied
2248 if (g%areaT(i,j) > 0.0) &
2249 frac_shelf_h(i,j) = area_shelf_h(i,j) / g%areaT(i,j)
2252 shelf_area => frac_shelf_h
2258 if (usealeremapping)
then
2259 call cpu_clock_begin(id_clock_ale)
2260 nkd = max(gv%ke, kd)
2263 allocate( tmp_mask_in(isd:ied,jsd:jed,nkd) ) ; tmp_mask_in(:,:,:) = 0.
2264 allocate( h1(isd:ied,jsd:jed,nkd) ) ; h1(:,:,:) = 0.
2265 allocate( tmpt1din(isd:ied,jsd:jed,nkd) ) ; tmpt1din(:,:,:) = 0.
2266 allocate( tmps1din(isd:ied,jsd:jed,nkd) ) ; tmps1din(:,:,:) = 0.
2267 do j = js, je ;
do i = is, ie
2268 if (g%mask2dT(i,j)>0.)
then
2269 ztopofcell = 0. ; zbottomofcell = 0.
2270 tmp_mask_in(i,j,1:kd) = mask_z(i,j,:)
2272 if (tmp_mask_in(i,j,k)>0. .and. k<=kd)
then
2273 zbottomofcell = max( z_edges_in(k+1), -g%bathyT(i,j) )
2274 tmpt1din(i,j,k) = temp_z(i,j,k)
2275 tmps1din(i,j,k) = salt_z(i,j,k)
2277 zbottomofcell = -g%bathyT(i,j)
2278 tmpt1din(i,j,k) = tmpt1din(i,j,k-1)
2279 tmps1din(i,j,k) = tmps1din(i,j,k-1)
2281 tmpt1din(i,j,k) = -99.9
2282 tmps1din(i,j,k) = -99.9
2284 h1(i,j,k) = gv%Z_to_H * (ztopofcell - zbottomofcell)
2285 ztopofcell = zbottomofcell
2287 h1(i,j,kd) = h1(i,j,kd) + gv%Z_to_H * max(0., ztopofcell + g%bathyT(i,j) )
2291 deallocate( tmp_mask_in )
2292 call pass_var(h1, g%Domain)
2293 call pass_var(tmpt1din, g%Domain)
2294 call pass_var(tmps1din, g%Domain)
2298 call ale_initregridding( gv, us, g%max_depth, pf, mdl, regridcs )
2300 if (.not. remap_general)
then
2302 allocate( htarget(nz) )
2303 htarget = getcoordinateresolution( regridcs )
2304 do j = js, je ;
do i = is, ie
2306 if (g%mask2dT(i,j)>0.)
then
2308 ztopofcell = 0. ; zbottomofcell = 0.
2310 zbottomofcell = max( ztopofcell - htarget(k), -g%bathyT(i,j) )
2311 h(i,j,k) = gv%Z_to_H * (ztopofcell - zbottomofcell)
2312 ztopofcell = zbottomofcell
2318 call pass_var(h, g%Domain)
2319 deallocate( htarget )
2323 call initialize_remapping( remapcs, remappingscheme, boundary_extrapolation=.false., answers_2018=answers_2018 )
2324 if (remap_general)
then
2325 call set_regrid_params( regridcs, min_thickness=0. )
2327 tv_loc%T => tmpt1din
2328 tv_loc%S => tmps1din
2331 allocate( dz_interface(isd:ied,jsd:jed,nkd+1) )
2332 if (use_ice_shelf)
then
2333 call regridding_main( remapcs, regridcs, g, gv_loc, h1, tv_loc, h, dz_interface, shelf_area )
2335 call regridding_main( remapcs, regridcs, g, gv_loc, h1, tv_loc, h, dz_interface )
2337 deallocate( dz_interface )
2339 call ale_remap_scalar(remapcs, g, gv, nkd, h1, tmpt1din, h, tv%T, all_cells=remap_full_column, &
2340 old_remap=remap_old_alg, answers_2018=answers_2018 )
2341 call ale_remap_scalar(remapcs, g, gv, nkd, h1, tmps1din, h, tv%S, all_cells=remap_full_column, &
2342 old_remap=remap_old_alg, answers_2018=answers_2018 )
2344 deallocate( tmpt1din )
2345 deallocate( tmps1din )
2347 call cpu_clock_end(id_clock_ale)
2353 nlevs = int(sum(mask_z,dim=3))
2357 do k=2,nz ; rb(k) = 0.5*(gv%Rlay(k-1)+gv%Rlay(k)) ;
enddo
2360 rb(nz+1) = 2.0*gv%Rlay(nz) - gv%Rlay(nz-1)
2362 rb(nz+1) = 2.0 * gv%Rlay(1)
2365 nkml = 0 ;
if (separate_mixed_layer) nkml = gv%nkml
2367 call find_interfaces(rho_z, z_in, kd, rb, g%bathyT, zi, g, us, &
2368 nlevs, nkml, hml=hmix_depth, eps_z=eps_z, eps_rho=eps_rho)
2370 if (correct_thickness)
then
2371 call adjustetatofitbathymetry(g, gv, us, zi, h)
2373 do k=nz,1,-1 ;
do j=js,je ;
do i=is,ie
2374 if (zi(i,j,k) < (zi(i,j,k+1) + gv%Angstrom_Z))
then
2375 zi(i,j,k) = zi(i,j,k+1) + gv%Angstrom_Z
2376 h(i,j,k) = gv%Angstrom_H
2378 h(i,j,k) = gv%Z_to_H * (zi(i,j,k) - zi(i,j,k+1))
2380 enddo ;
enddo ;
enddo
2382 do j=js,je ;
do i=is,ie
2383 if (abs(zi(i,j,nz+1) + g%bathyT(i,j)) > 1.0*us%m_to_Z) &
2384 inconsistent = inconsistent + 1
2386 call sum_across_pes(inconsistent)
2388 if ((inconsistent > 0) .and. (is_root_pe()))
then
2389 write(mesg,
'("Thickness initial conditions are inconsistent ",'// &
2390 '"with topography in ",I5," places.")') inconsistent
2391 call mom_error(warning, mesg)
2395 call tracer_z_init_array(temp_z, z_edges_in, kd, zi, missing_value, g, nz, nlevs, eps_z, tv%T)
2396 call tracer_z_init_array(salt_z, z_edges_in, kd, zi, missing_value, g, nz, nlevs, eps_z, tv%S)
2399 npoints = 0 ; tempavg = 0. ; saltavg = 0.
2400 do j=js,je ;
do i=is,ie ;
if (g%mask2dT(i,j) >= 1.0)
then
2401 npoints = npoints + 1
2402 tempavg = tempavg + tv%T(i,j,k)
2403 saltavg = saltavg + tv%S(i,j,k)
2404 endif ;
enddo ;
enddo
2407 if (homogenize)
then
2408 call sum_across_pes(npoints)
2409 call sum_across_pes(tempavg)
2410 call sum_across_pes(saltavg)
2412 tempavg = tempavg / real(npoints)
2413 saltavg = saltavg / real(npoints)
2415 tv%T(:,:,k) = tempavg
2416 tv%S(:,:,k) = saltavg
2423 do k=1,nz ;
do j=js,je ;
do i=is,ie
2424 if (tv%T(i,j,k) == missing_value)
then
2425 tv%T(i,j,k) = temp_land_fill
2426 tv%S(i,j,k) = salt_land_fill
2428 enddo ;
enddo ;
enddo
2431 if (adjust_temperature .and. .not. usealeremapping)
then
2433 ks = 1 ;
if (separate_mixed_layer) ks = gv%nk_rho_varies + 1
2434 call determine_temperature(tv%T, tv%S, gv%Rlay(1:nz), tv%P_Ref, niter, &
2435 missing_value, h, ks, g, us, eos)
2438 deallocate(z_in, z_edges_in, temp_z, salt_z, mask_z)
2439 deallocate(rho_z, area_shelf_h, frac_shelf_h)
2441 call pass_var(h, g%Domain)
2442 call pass_var(tv%T, g%Domain)
2443 call pass_var(tv%S, g%Domain)
2445 call calltree_leave(trim(mdl)//
'()')
2446 call cpu_clock_end(id_clock_routine)
2448 end subroutine mom_temp_salt_initialize_from_z
2452 subroutine find_interfaces(rho, zin, nk_data, Rb, depth, zi, G, US, nlevs, nkml, hml, &
2454 type(ocean_grid_type),
intent(in) :: G
2455 integer,
intent(in) :: nk_data
2456 real,
dimension(SZI_(G),SZJ_(G),nk_data), &
2458 real,
dimension(nk_data),
intent(in) :: zin
2459 real,
dimension(SZK_(G)+1),
intent(in) :: Rb
2460 real,
dimension(SZI_(G),SZJ_(G)), &
2462 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1), &
2464 type(unit_scale_type),
intent(in) :: US
2465 integer,
dimension(SZI_(G),SZJ_(G)), &
2467 integer,
intent(in) :: nkml
2469 real,
intent(in) :: hml
2470 real,
intent(in) :: eps_z
2471 real,
intent(in) :: eps_rho
2474 real,
dimension(nk_data) :: rho_
2475 real,
dimension(SZK_(G)+1) :: zi_
2478 real,
parameter :: zoff=0.999
2480 integer :: nlevs_data
2481 logical :: work_down
2482 integer :: k_int, lo_int, hi_int, mid
2483 integer :: i, j, k, is, ie, js, je, nz
2485 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
2489 do j=js,je ;
do i=is,ie
2490 nlevs_data = nlevs(i,j)
2491 do k=1,nlevs_data ; rho_(k) = rho(i,j,k) ;
enddo
2499 do k=2,nlevs_data-1 ;
if (rho_(k) - rho_(k-1) < 0.0 )
then
2501 rho_(k-1) = rho_(k) - eps_rho
2503 drhodz = (rho_(k+1)-rho_(k-1)) / (zin(k+1)-zin(k-1))
2504 if (drhodz < 0.0) unstable=.true.
2505 rho_(k) = rho_(k-1) + drhodz*zoff*(zin(k)-zin(k-1))
2510 do k=nlevs_data-1,2,-1 ;
if (rho_(k+1) - rho_(k) < 0.0)
then
2511 if (k == nlevs_data-1)
then
2512 rho_(k+1) = rho_(k-1) + eps_rho
2514 drhodz = (rho_(k+1)-rho_(k-1)) / (zin(k+1)-zin(k-1))
2515 if (drhodz < 0.0) unstable=.true.
2516 rho_(k) = rho_(k+1) - drhodz*(zin(k+1)-zin(k))
2529 lo_int = 1 ; hi_int = nlevs_data
2530 do while (lo_int < hi_int)
2531 mid = (lo_int+hi_int) / 2
2532 if (rb(k) < rho_(mid))
then ; hi_int = mid
2533 else ; lo_int = mid+1 ;
endif
2535 k_int = max(1, lo_int-1)
2538 slope = (zin(k_int+1) - zin(k_int)) / max(rho_(k_int+1) - rho_(k_int), eps_rho)
2539 zi_(k) = -1.0*(zin(k_int) + slope*(rb(k)-rho_(k_int)))
2540 zi_(k) = min(max(zi_(k), -depth(i,j)), -1.0*hml)
2542 zi_(nz+1) = -depth(i,j)
2543 if (nkml > 0)
then ;
do k=2,nkml+1
2544 zi_(k) = max(hml*((1.0-real(k))/real(nkml)), -depth(i,j))
2546 do k=nz,max(nkml+2,2),-1
2547 if (zi_(k) < zi_(k+1) + eps_z) zi_(k) = zi_(k+1) + eps_z
2548 if (zi_(k) > -1.0*hml) zi_(k) = max(-1.0*hml, -depth(i,j))
2556 end subroutine find_interfaces
2559 subroutine mom_state_init_tests(G, GV, US, tv)
2560 type(ocean_grid_type),
intent(inout) :: G
2561 type(verticalgrid_type),
intent(in) :: GV
2562 type(unit_scale_type),
intent(in) :: US
2563 type(thermo_var_ptrs),
intent(in) :: tv
2566 integer,
parameter :: nk=5
2567 real,
dimension(nk) :: T, T_t, T_b
2568 real,
dimension(nk) :: S, S_t, S_b
2569 real,
dimension(nk) :: rho
2570 real,
dimension(nk) :: h
2571 real,
dimension(nk) :: z
2572 real,
dimension(nk+1) :: e
2574 real :: P_tot, P_t, P_b
2577 type(remapping_cs),
pointer :: remap_CS => null()
2579 i_z_scale = 1.0 / (500.0*us%m_to_Z)
2581 h(k) = 100.0*gv%m_to_H
2585 e(k+1) = e(k) - gv%H_to_Z * h(k)
2589 z(k) = 0.5 * ( e(k) + e(k+1) )
2590 t_t(k) = 20. + (0. * i_z_scale) * e(k)
2591 t(k) = 20. + (0. * i_z_scale)*z(k)
2592 t_b(k) = 20. + (0. * i_z_scale)*e(k+1)
2593 s_t(k) = 35. - (0. * i_z_scale)*e(k)
2594 s(k) = 35. + (0. * i_z_scale)*z(k)
2595 s_b(k) = 35. - (0. * i_z_scale)*e(k+1)
2596 call calculate_density(0.5*(t_t(k)+t_b(k)), 0.5*(s_t(k)+s_b(k)), -gv%Rho0*gv%g_Earth*us%m_to_Z*z(k), &
2597 rho(k), tv%eqn_of_state)
2598 p_tot = p_tot + gv%g_Earth * rho(k) * gv%H_to_Z*h(k)
2603 call find_depth_of_pressure_in_cell(t_t(k), t_b(k), s_t(k), s_b(k), e(k), e(k+1), p_t, 0.5*p_tot, &
2604 gv%Rho0, gv%g_Earth, tv%eqn_of_state, us, p_b, z_out)
2605 write(0,*) k, us%RL2_T2_to_Pa*p_t, us%RL2_T2_to_Pa*p_b, 0.5*us%RL2_T2_to_Pa*p_tot, &
2606 us%Z_to_m*e(k), us%Z_to_m*e(k+1), us%Z_to_m*z_out
2609 write(0,*) us%RL2_T2_to_Pa*p_b, us%RL2_T2_to_Pa*p_tot
2612 write(0,*)
' ==================================================================== '
2614 write(0,*) gv%H_to_m*h
2615 call cut_off_column_top(nk, tv, gv, us, gv%g_Earth, -e(nk+1), gv%Angstrom_Z, &
2616 t, t_t, t_b, s, s_t, s_b, 0.5*p_tot, h, remap_cs)
2617 write(0,*) gv%H_to_m*h
2619 end subroutine mom_state_init_tests