15 use mom,
only : extract_surface_state, allocate_surface_state, finish_mom_initialization
16 use mom,
only : get_mom_state_elements, mom_state_is_synchronized
17 use mom,
only : get_ocean_stocks, step_offline
36 use mom_surface_forcing_gfdl,
only : surface_forcing_init, convert_iob_to_fluxes
37 use mom_surface_forcing_gfdl,
only : convert_iob_to_forces, ice_ocn_bnd_type_chksum
39 use mom_surface_forcing_gfdl,
only : forcing_save_restart
40 use mom_time_manager,
only : time_type,
operator(>),
operator(+),
operator(-)
50 use mom_ice_shelf,
only : add_shelf_forces, ice_shelf_end, ice_shelf_save_restart
51 use coupler_types_mod,
only : coupler_1d_bc_type, coupler_2d_bc_type
52 use coupler_types_mod,
only : coupler_type_spawn, coupler_type_write_chksums
53 use coupler_types_mod,
only : coupler_type_initialized, coupler_type_copy_data
54 use coupler_types_mod,
only : coupler_type_set_diags, coupler_type_send_data
55 use mpp_domains_mod,
only : domain2d, mpp_get_layout, mpp_get_global_domain
56 use mpp_domains_mod,
only : mpp_define_domains, mpp_get_compute_domain, mpp_get_data_domain
58 use fms_mod,
only : stdout
59 use mpp_mod,
only : mpp_chksum
60 use mom_eos,
only : gsw_sp_from_sr, gsw_pt_from_ct
64 #include <MOM_memory.h>
66 #ifdef _USE_GENERIC_TRACER
70 implicit none ;
private
72 public ocean_model_init, ocean_model_end, update_ocean_model
73 public ocean_model_save_restart, ocean_stock_pe
75 public ocean_model_init_sfc, ocean_model_flux_init
76 public ocean_model_restart
77 public ice_ocn_bnd_type_chksum
78 public ocean_public_type_chksum
81 public ocean_model_get_uv_surf
85 module procedure ocean_model_data1d_get
86 module procedure ocean_model_data2d_get
94 type(domain2d) :: domain
95 logical :: is_ocean_pe
96 character(len=32) :: instance_name =
''
99 integer,
pointer,
dimension(:) :: pelist => null()
100 logical,
pointer,
dimension(:,:) :: maskmap =>null()
106 integer :: stagger = -999
115 real,
pointer,
dimension(:,:) :: &
125 type(coupler_2d_bc_type) :: fields
130 integer,
dimension(2) :: axes = 0
140 logical :: is_ocean_pe = .false.
141 type(time_type) :: time
142 type(time_type) :: time_dyn
144 integer :: restart_control
152 integer :: nstep_thermo = 0
153 logical :: use_ice_shelf
156 logical :: icebergs_alter_ocean
161 logical :: offline_tracer_mode = .false.
169 logical :: single_step_call
176 logical :: thermo_spans_coupling
178 logical :: diabatic_first
202 ice_shelf_csp => null()
206 marine_ice_csp => null()
211 forcing_csp => null()
213 restart_csp => null()
227 subroutine ocean_model_init(Ocean_sfc, OS, Time_init, Time_in, wind_stagger, gas_fields_ocn)
229 intent(inout) :: ocean_sfc
235 type(time_type),
intent(in) :: time_init
236 type(time_type),
intent(in) :: time_in
237 integer,
optional,
intent(in) :: wind_stagger
239 type(coupler_1d_bc_type), &
240 optional,
intent(in) :: gas_fields_ocn
252 logical :: use_melt_pot
255 #include "version_variable.h"
256 character(len=40) :: mdl =
"ocean_model_init"
257 character(len=48) :: stagger
260 logical :: use_temperature
262 call calltree_enter(
"ocean_model_init(), ocean_model_MOM.F90")
263 if (
associated(os))
then
264 call mom_error(warning,
"ocean_model_init called with an associated "// &
265 "ocean_state_type structure. Model is already initialized.")
270 os%is_ocean_pe = ocean_sfc%is_ocean_pe
271 if (.not.os%is_ocean_pe)
return
273 os%Time = time_in ; os%Time_dyn = time_in
274 call initialize_mom(os%Time, time_init, param_file, os%dirs, os%MOM_CSp, &
275 os%restart_CSp, time_in, offline_tracer_mode=os%offline_tracer_mode, &
276 diag_ptr=os%diag, count_calls=.true.)
277 call get_mom_state_elements(os%MOM_CSp, g=os%grid, gv=os%GV, us=os%US, c_p=os%C_p, &
278 c_p_scaled=os%fluxes%C_p, use_temp=use_temperature)
283 call get_param(param_file, mdl,
"SINGLE_STEPPING_CALL", os%single_step_call, &
284 "If true, advance the state of MOM with a single step "//&
285 "including both dynamics and thermodynamics. If false, "//&
286 "the two phases are advanced with separate calls.", default=.true.)
287 call get_param(param_file, mdl,
"DT", os%dt, &
288 "The (baroclinic) dynamics time step. The time-step that "//&
289 "is actually used will be an integer fraction of the "//&
290 "forcing time-step.", units=
"s", fail_if_missing=.true.)
291 call get_param(param_file, mdl,
"DT_THERM", os%dt_therm, &
292 "The thermodynamic and tracer advection time step. "//&
293 "Ideally DT_THERM should be an integer multiple of DT "//&
294 "and less than the forcing or coupling time-step, unless "//&
295 "THERMO_SPANS_COUPLING is true, in which case DT_THERM "//&
296 "can be an integer multiple of the coupling timestep. By "//&
297 "default DT_THERM is set to DT.", units=
"s", default=os%dt)
298 call get_param(param_file,
"MOM",
"THERMO_SPANS_COUPLING", os%thermo_spans_coupling, &
299 "If true, the MOM will take thermodynamic and tracer "//&
300 "timesteps that can be longer than the coupling timestep. "//&
301 "The actual thermodynamic timestep that is used in this "//&
302 "case is the largest integer multiple of the coupling "//&
303 "timestep that is less than or equal to DT_THERM.", default=.false.)
304 call get_param(param_file, mdl,
"DIABATIC_FIRST", os%diabatic_first, &
305 "If true, apply diabatic and thermodynamic processes, "//&
306 "including buoyancy forcing and mass gain or loss, "//&
307 "before stepping the dynamics forward.", default=.false.)
309 call get_param(param_file, mdl,
"RESTART_CONTROL", os%Restart_control, &
310 "An integer whose bits encode which restart files are "//&
311 "written. Add 2 (bit 1) for a time-stamped file, and odd "//&
312 "(bit 0) for a non-time-stamped file. A restart file "//&
313 "will be saved at the end of the run segment for any "//&
314 "non-negative value.", default=1)
315 call get_param(param_file, mdl,
"OCEAN_SURFACE_STAGGER", stagger, &
316 "A case-insensitive character string to indicate the "//&
317 "staggering of the surface velocity field that is "//&
318 "returned to the coupler. Valid values include "//&
319 "'A', 'B', or 'C'.", default=
"C")
320 if (uppercase(stagger(1:1)) ==
'A')
then ; ocean_sfc%stagger = agrid
321 elseif (uppercase(stagger(1:1)) ==
'B')
then ; ocean_sfc%stagger = bgrid_ne
322 elseif (uppercase(stagger(1:1)) ==
'C')
then ; ocean_sfc%stagger = cgrid_ne
323 else ;
call mom_error(fatal,
"ocean_model_init: OCEAN_SURFACE_STAGGER = "// &
324 trim(stagger)//
" is invalid.") ;
endif
326 call get_param(param_file, mdl,
"RHO_0", rho0, &
327 "The mean ocean density used with BOUSSINESQ true to "//&
328 "calculate accelerations and the mass for conservation "//&
329 "properties, or with BOUSSINSEQ false to convert some "//&
330 "parameters from vertical units of m to kg m-2.", &
331 units=
"kg m-3", default=1035.0)
332 call get_param(param_file, mdl,
"G_EARTH", g_earth, &
333 "The gravitational acceleration of the Earth.", &
334 units=
"m s-2", default = 9.80)
336 call get_param(param_file, mdl,
"ICE_SHELF", os%use_ice_shelf, &
337 "If true, enables the ice shelf model.", default=.false.)
339 call get_param(param_file, mdl,
"ICEBERGS_APPLY_RIGID_BOUNDARY", os%icebergs_alter_ocean, &
340 "If true, allows icebergs to change boundary condition felt by ocean", default=.false.)
342 os%press_to_z = 1.0/(rho0*g_earth)
346 call get_param(param_file, mdl,
"HFREEZE", hfrz, &
347 "If HFREEZE > 0, melt potential will be computed. The actual depth "//&
348 "over which melt potential is computed will be min(HFREEZE, OBLD), "//&
349 "where OBLD is the boundary layer depth. If HFREEZE <= 0 (default), "//&
350 "melt potential will not be computed.", units=
"m", default=-1.0, do_not_log=.true.)
352 if (hfrz .gt. 0.0)
then
358 call allocate_surface_state(os%sfc_state, os%grid, use_temperature, do_integrals=.true., &
359 gas_fields_ocn=gas_fields_ocn, use_meltpot=use_melt_pot)
361 if (
present(wind_stagger))
then
362 call surface_forcing_init(time_in, os%grid, os%US, param_file, os%diag, &
363 os%forcing_CSp, wind_stagger)
365 call surface_forcing_init(time_in, os%grid, os%US, param_file, os%diag, &
369 if (os%use_ice_shelf)
then
370 call initialize_ice_shelf(param_file, os%grid, os%Time, os%ice_shelf_CSp, &
371 os%diag, os%forces, os%fluxes)
373 if (os%icebergs_alter_ocean)
then
374 call marine_ice_init(os%Time, os%grid, param_file, os%diag, os%marine_ice_CSp)
375 if (.not. os%use_ice_shelf) &
379 call get_param(param_file, mdl,
"USE_WAVES", os%Use_Waves, &
380 "If true, enables surface wave modules.", default=.false.)
381 if (os%use_waves)
then
382 call mom_wave_interface_init(os%Time, os%grid, os%GV, os%US, param_file, os%Waves, os%diag)
384 call mom_wave_interface_init_lite(param_file)
387 if (
associated(os%grid%Domain%maskmap))
then
388 call initialize_ocean_public_type(os%grid%Domain%mpp_domain, ocean_sfc, &
389 os%diag, maskmap=os%grid%Domain%maskmap, &
390 gas_fields_ocn=gas_fields_ocn)
392 call initialize_ocean_public_type(os%grid%Domain%mpp_domain, ocean_sfc, &
393 os%diag, gas_fields_ocn=gas_fields_ocn)
398 if (
present(gas_fields_ocn))
then
399 call coupler_type_set_diags(ocean_sfc%fields,
"ocean_sfc", &
400 ocean_sfc%axes(1:2), time_in)
402 call extract_surface_state(os%MOM_CSp, os%sfc_state)
404 call convert_state_to_ocean_type(os%sfc_state, ocean_sfc, os%grid, os%US)
408 call close_param_file(param_file)
409 call diag_mediator_close_registration(os%diag)
411 call mom_mesg(
'==== Completed MOM6 Coupled Initialization ====', 2)
413 call calltree_leave(
"ocean_model_init(")
414 end subroutine ocean_model_init
421 subroutine update_ocean_model(Ice_ocean_boundary, OS, Ocean_sfc, time_start_update, &
422 Ocean_coupling_time_step, update_dyn, update_thermo, &
423 Ocn_fluxes_used, start_cycle, end_cycle, cycle_length)
425 intent(in) :: ice_ocean_boundary
431 intent(inout) :: ocean_sfc
434 type(time_type),
intent(in) :: time_start_update
435 type(time_type),
intent(in) :: ocean_coupling_time_step
437 logical,
optional,
intent(in) :: update_dyn
439 logical,
optional,
intent(in) :: update_thermo
441 logical,
optional,
intent(in) :: ocn_fluxes_used
444 logical,
optional,
intent(in) :: start_cycle
447 logical,
optional,
intent(in) :: end_cycle
450 real,
optional,
intent(in) :: cycle_length
453 type(time_type) :: time_seg_start
456 type(time_type) :: time_thermo_start
459 type(time_type) :: time1
460 integer :: index_bnds(4)
466 real :: t_elapsed_seg
470 integer :: n_last_thermo
471 logical :: thermo_does_span_coupling
474 logical :: step_thermo
475 integer :: is, ie, js, je
477 call calltree_enter(
"update_ocean_model(), ocean_model_MOM.F90")
478 dt_coupling = time_type_to_real(ocean_coupling_time_step)
480 if (.not.
associated(os))
then
481 call mom_error(fatal,
"update_ocean_model called with an unassociated "// &
482 "ocean_state_type structure. ocean_model_init must be "// &
483 "called first to allocate this structure.")
487 do_dyn = .true. ;
if (
present(update_dyn)) do_dyn = update_dyn
488 do_thermo = .true. ;
if (
present(update_thermo)) do_thermo = update_thermo
490 if (do_thermo .and. (time_start_update /= os%Time)) &
491 call mom_error(warning,
"update_ocean_model: internal clock does not "//&
492 "agree with time_start_update argument.")
493 if (do_dyn .and. (time_start_update /= os%Time_dyn)) &
494 call mom_error(warning,
"update_ocean_model: internal dynamics clock does not "//&
495 "agree with time_start_update argument.")
497 if (.not.(do_dyn .or. do_thermo))
call mom_error(fatal, &
498 "update_ocean_model called without updating either dynamics or thermodynamics.")
499 if (do_dyn .and. do_thermo .and. (os%Time /= os%Time_dyn))
call mom_error(fatal, &
500 "update_ocean_model called to update both dynamics and thermodynamics with inconsistent clocks.")
504 is = os%grid%isc ; ie = os%grid%iec ; js = os%grid%jsc ; je = os%grid%jec
505 call coupler_type_spawn(ocean_sfc%fields, os%sfc_state%tr_fields, &
506 (/is,is,ie,ie/), (/js,js,je,je/), as_needed=.true.)
509 call mpp_get_compute_domain(ocean_sfc%Domain, index_bnds(1), index_bnds(2), &
510 index_bnds(3), index_bnds(4))
513 call convert_iob_to_forces(ice_ocean_boundary, os%forces, index_bnds, os%Time_dyn, os%grid, os%US, &
514 os%forcing_CSp, dt_forcing=dt_coupling, reset_avg=os%fluxes%fluxes_used)
515 if (os%use_ice_shelf) &
516 call add_shelf_forces(os%grid, os%US, os%Ice_shelf_CSp, os%forces)
517 if (os%icebergs_alter_ocean) &
518 call iceberg_forces(os%grid, os%forces, os%use_ice_shelf, &
519 os%sfc_state, dt_coupling, os%marine_ice_CSp)
523 if (os%fluxes%fluxes_used)
then
524 call convert_iob_to_fluxes(ice_ocean_boundary, os%fluxes, index_bnds, os%Time, dt_coupling, &
525 os%grid, os%US, os%forcing_CSp, os%sfc_state)
528 if (os%use_ice_shelf) &
529 call shelf_calc_flux(os%sfc_state, os%fluxes, os%Time, dt_coupling, os%Ice_shelf_CSp)
530 if (os%icebergs_alter_ocean) &
531 call iceberg_fluxes(os%grid, os%US, os%fluxes, os%use_ice_shelf, &
532 os%sfc_state, dt_coupling, os%marine_ice_CSp)
534 #ifdef _USE_GENERIC_TRACER
535 call enable_averaging(dt_coupling, os%Time + ocean_coupling_time_step, os%diag)
536 call mom_generic_tracer_fluxes_accumulate(os%fluxes, 1.0)
537 call disable_averaging(os%diag)
542 os%flux_tmp%C_p = os%fluxes%C_p
543 call convert_iob_to_fluxes(ice_ocean_boundary, os%flux_tmp, index_bnds, os%Time, dt_coupling, &
544 os%grid, os%US, os%forcing_CSp, os%sfc_state)
546 if (os%use_ice_shelf) &
547 call shelf_calc_flux(os%sfc_state, os%flux_tmp, os%Time, dt_coupling, os%Ice_shelf_CSp)
548 if (os%icebergs_alter_ocean) &
549 call iceberg_fluxes(os%grid, os%US, os%flux_tmp, os%use_ice_shelf, &
550 os%sfc_state, dt_coupling, os%marine_ice_CSp)
552 call fluxes_accumulate(os%flux_tmp, os%fluxes, os%grid, weight)
553 #ifdef _USE_GENERIC_TRACER
555 call mom_generic_tracer_fluxes_accumulate(os%flux_tmp, weight)
561 if (do_dyn .and.
associated(os%forces%net_mass_src) .and. .not.os%forces%net_mass_src_set) &
562 call get_net_mass_forcing(os%fluxes, os%grid, os%US, os%forces%net_mass_src)
564 if (os%use_waves .and. do_thermo)
then
568 call update_surface_waves(os%grid, os%GV, os%US, os%time, ocean_coupling_time_step, os%waves)
571 if ((os%nstep==0) .and. (os%nstep_thermo==0))
then
572 call finish_mom_initialization(os%Time, os%dirs, os%MOM_CSp, os%restart_CSp)
575 time_thermo_start = os%Time
576 time_seg_start = os%Time ;
if (do_dyn) time_seg_start = os%Time_dyn
577 time1 = time_seg_start
579 if (os%offline_tracer_mode .and. do_thermo)
then
580 call step_offline(os%forces, os%fluxes, os%sfc_state, time1, dt_coupling, os%MOM_CSp)
581 elseif ((.not.do_thermo) .or. (.not.do_dyn))
then
583 call step_mom(os%forces, os%fluxes, os%sfc_state, time1, dt_coupling, os%MOM_CSp, &
584 waves=os%Waves, do_dynamics=do_dyn, do_thermodynamics=do_thermo, &
585 start_cycle=start_cycle, end_cycle=end_cycle, cycle_length=cycle_length, &
586 reset_therm=ocn_fluxes_used)
587 elseif (os%single_step_call)
then
588 call step_mom(os%forces, os%fluxes, os%sfc_state, time1, dt_coupling, os%MOM_CSp, waves=os%Waves)
590 n_max = 1 ;
if (dt_coupling > os%dt) n_max = ceiling(dt_coupling/os%dt - 0.001)
591 dt_dyn = dt_coupling / real(n_max)
592 thermo_does_span_coupling = (os%thermo_spans_coupling .and. (os%dt_therm > 1.5*dt_coupling))
594 if (thermo_does_span_coupling)
then
595 dt_therm = dt_coupling * floor(os%dt_therm / dt_coupling + 0.001)
596 nts = floor(dt_therm/dt_dyn + 0.001)
598 nts = max(1,min(n_max,floor(os%dt_therm/dt_dyn + 0.001)))
602 time1 = time_seg_start ; t_elapsed_seg = 0.0
604 if (os%diabatic_first)
then
605 if (thermo_does_span_coupling)
call mom_error(fatal, &
606 "MOM is not yet set up to have restarts that work with "//&
607 "THERMO_SPANS_COUPLING and DIABATIC_FIRST.")
608 if (modulo(n-1,nts)==0)
then
609 dtdia = dt_dyn*min(nts,n_max-(n-1))
610 call step_mom(os%forces, os%fluxes, os%sfc_state, time1, dtdia, os%MOM_CSp, &
611 waves=os%Waves, do_dynamics=.false., do_thermodynamics=.true., &
612 start_cycle=(n==1), end_cycle=.false., cycle_length=dt_coupling)
615 call step_mom(os%forces, os%fluxes, os%sfc_state, time1, dt_dyn, os%MOM_CSp, &
616 waves=os%Waves, do_dynamics=.true., do_thermodynamics=.false., &
617 start_cycle=.false., end_cycle=(n==n_max), cycle_length=dt_coupling)
619 call step_mom(os%forces, os%fluxes, os%sfc_state, time1, dt_dyn, os%MOM_CSp, &
620 waves=os%Waves, do_dynamics=.true., do_thermodynamics=.false., &
621 start_cycle=(n==1), end_cycle=.false., cycle_length=dt_coupling)
623 step_thermo = .false.
624 if (thermo_does_span_coupling)
then
626 step_thermo = mom_state_is_synchronized(os%MOM_CSp, adv_dyn=.true.)
627 elseif ((modulo(n,nts)==0) .or. (n==n_max))
then
628 dtdia = dt_dyn*(n - n_last_thermo)
633 if (step_thermo)
then
635 time1 = time1 - real_to_time(dtdia - dt_dyn)
636 call step_mom(os%forces, os%fluxes, os%sfc_state, time1, dtdia, os%MOM_CSp, &
637 waves=os%Waves, do_dynamics=.false., do_thermodynamics=.true., &
638 start_cycle=.false., end_cycle=(n==n_max), cycle_length=dt_coupling)
642 t_elapsed_seg = t_elapsed_seg + dt_dyn
643 time1 = time_seg_start + real_to_time(t_elapsed_seg)
647 if (do_dyn) os%Time_dyn = time_seg_start + ocean_coupling_time_step
648 if (do_dyn) os%nstep = os%nstep + 1
649 os%Time = time_thermo_start
650 if (do_thermo) os%Time = os%Time + ocean_coupling_time_step
651 if (do_thermo) os%nstep_thermo = os%nstep_thermo + 1
654 call mech_forcing_diags(os%forces, dt_coupling, os%grid, os%Time_dyn, os%diag, os%forcing_CSp%handles)
657 if (os%fluxes%fluxes_used .and. do_thermo)
then
658 call forcing_diagnostics(os%fluxes, os%sfc_state, os%grid, os%US, os%Time, os%diag, os%forcing_CSp%handles)
664 call convert_state_to_ocean_type(os%sfc_state, ocean_sfc, os%grid, os%US)
665 time1 = os%Time ;
if (do_dyn) time1 = os%Time_dyn
666 call coupler_type_send_data(ocean_sfc%fields, time1)
668 call calltree_leave(
"update_ocean_model()")
669 end subroutine update_ocean_model
672 subroutine ocean_model_restart(OS, timestamp)
675 character(len=*),
optional,
intent(in) :: timestamp
678 if (.not.mom_state_is_synchronized(os%MOM_CSp)) &
679 call mom_error(warning,
"End of MOM_main reached with inconsistent "//&
680 "dynamics and advective times. Additional restart fields "//&
681 "that have not been coded yet would be required for reproducibility.")
682 if (.not.os%fluxes%fluxes_used)
call mom_error(fatal,
"ocean_model_restart "//&
683 "was called with unused buoyancy fluxes. For conservation, the ocean "//&
684 "restart files can only be created after the buoyancy forcing is applied.")
686 if (btest(os%Restart_control,1))
then
687 call save_restart(os%dirs%restart_output_dir, os%Time, os%grid, &
688 os%restart_CSp, .true., gv=os%GV)
689 call forcing_save_restart(os%forcing_CSp, os%grid, os%Time, &
690 os%dirs%restart_output_dir, .true.)
691 if (os%use_ice_shelf)
then
692 call ice_shelf_save_restart(os%Ice_shelf_CSp, os%Time, os%dirs%restart_output_dir, .true.)
695 if (btest(os%Restart_control,0))
then
696 call save_restart(os%dirs%restart_output_dir, os%Time, os%grid, &
697 os%restart_CSp, gv=os%GV)
698 call forcing_save_restart(os%forcing_CSp, os%grid, os%Time, &
699 os%dirs%restart_output_dir)
700 if (os%use_ice_shelf)
then
701 call ice_shelf_save_restart(os%Ice_shelf_CSp, os%Time, os%dirs%restart_output_dir)
705 end subroutine ocean_model_restart
710 subroutine ocean_model_end(Ocean_sfc, Ocean_state, Time)
716 type(time_type),
intent(in) :: time
718 call ocean_model_save_restart(ocean_state, time)
719 call diag_mediator_end(time, ocean_state%diag)
720 call mom_end(ocean_state%MOM_CSp)
721 if (ocean_state%use_ice_shelf)
call ice_shelf_end(ocean_state%Ice_shelf_CSp)
722 end subroutine ocean_model_end
726 subroutine ocean_model_save_restart(OS, Time, directory, filename_suffix)
729 type(time_type),
intent(in) :: time
730 character(len=*),
optional,
intent(in) :: directory
732 character(len=*),
optional,
intent(in) :: filename_suffix
737 character(len=200) :: restart_dir
739 if (.not.mom_state_is_synchronized(os%MOM_CSp)) &
740 call mom_error(warning,
"ocean_model_save_restart called with inconsistent "//&
741 "dynamics and advective times. Additional restart fields "//&
742 "that have not been coded yet would be required for reproducibility.")
743 if (.not.os%fluxes%fluxes_used)
call mom_error(fatal,
"ocean_model_save_restart "//&
744 "was called with unused buoyancy fluxes. For conservation, the ocean "//&
745 "restart files can only be created after the buoyancy forcing is applied.")
747 if (
present(directory))
then ; restart_dir = directory
748 else ; restart_dir = os%dirs%restart_output_dir ;
endif
750 call save_restart(restart_dir, time, os%grid, os%restart_CSp, gv=os%GV)
752 call forcing_save_restart(os%forcing_CSp, os%grid, time, restart_dir)
754 if (os%use_ice_shelf)
then
755 call ice_shelf_save_restart(os%Ice_shelf_CSp, os%Time, os%dirs%restart_output_dir)
758 end subroutine ocean_model_save_restart
761 subroutine initialize_ocean_public_type(input_domain, Ocean_sfc, diag, maskmap, &
763 type(domain2d),
intent(in) :: input_domain
768 logical,
dimension(:,:), &
769 optional,
intent(in) :: maskmap
771 type(coupler_1d_bc_type), &
772 optional,
intent(in) :: gas_fields_ocn
777 integer :: xsz, ysz, layout(2)
780 integer :: isc, iec, jsc, jec
782 call mpp_get_layout(input_domain,layout)
783 call mpp_get_global_domain(input_domain, xsize=xsz, ysize=ysz)
784 if (
PRESENT(maskmap))
then
785 call mpp_define_domains((/1,xsz,1,ysz/),layout,ocean_sfc%Domain, maskmap=maskmap)
787 call mpp_define_domains((/1,xsz,1,ysz/),layout,ocean_sfc%Domain)
789 call mpp_get_compute_domain(ocean_sfc%Domain, isc, iec, jsc, jec)
791 allocate ( ocean_sfc%t_surf (isc:iec,jsc:jec), &
792 ocean_sfc%s_surf (isc:iec,jsc:jec), &
793 ocean_sfc%u_surf (isc:iec,jsc:jec), &
794 ocean_sfc%v_surf (isc:iec,jsc:jec), &
795 ocean_sfc%sea_lev(isc:iec,jsc:jec), &
796 ocean_sfc%area (isc:iec,jsc:jec), &
797 ocean_sfc%frazil (isc:iec,jsc:jec))
799 ocean_sfc%t_surf(:,:) = 0.0
800 ocean_sfc%s_surf(:,:) = 0.0
801 ocean_sfc%u_surf(:,:) = 0.0
802 ocean_sfc%v_surf(:,:) = 0.0
803 ocean_sfc%sea_lev(:,:) = 0.0
804 ocean_sfc%frazil(:,:) = 0.0
805 ocean_sfc%area(:,:) = 0.0
806 ocean_sfc%axes = diag%axesT1%handles
808 if (
present(gas_fields_ocn))
then
809 call coupler_type_spawn(gas_fields_ocn, ocean_sfc%fields, (/isc,isc,iec,iec/), &
810 (/jsc,jsc,jec,jec/), suffix =
'_ocn', as_needed=.true.)
813 end subroutine initialize_ocean_public_type
820 subroutine convert_state_to_ocean_type(sfc_state, Ocean_sfc, G, US, patm, press_to_z)
821 type(
surface),
intent(inout) :: sfc_state
824 target,
intent(inout) :: Ocean_sfc
829 real,
optional,
intent(in) :: patm(:,:)
830 real,
optional,
intent(in) :: press_to_z
834 character(len=48) :: val_str
835 integer :: isc_bnd, iec_bnd, jsc_bnd, jec_bnd
836 integer :: i, j, i0, j0, is, ie, js, je
838 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
839 call pass_vector(sfc_state%u, sfc_state%v, g%Domain)
841 call mpp_get_compute_domain(ocean_sfc%Domain, isc_bnd, iec_bnd, &
843 if (
present(patm))
then
845 if (.not.
present(press_to_z))
call mom_error(fatal, &
846 'convert_state_to_ocean_type: press_to_z must be present if patm is.')
849 i0 = is - isc_bnd ; j0 = js - jsc_bnd
850 if (sfc_state%T_is_conT)
then
852 do j=jsc_bnd,jec_bnd ;
do i=isc_bnd,iec_bnd
853 ocean_sfc%t_surf(i,j) = gsw_pt_from_ct(sfc_state%SSS(i+i0,j+j0), &
854 sfc_state%SST(i+i0,j+j0)) + celsius_kelvin_offset
857 do j=jsc_bnd,jec_bnd ;
do i=isc_bnd,iec_bnd
858 ocean_sfc%t_surf(i,j) = sfc_state%SST(i+i0,j+j0) + celsius_kelvin_offset
861 if (sfc_state%S_is_absS)
then
863 do j=jsc_bnd,jec_bnd ;
do i=isc_bnd,iec_bnd
864 ocean_sfc%s_surf(i,j) = gsw_sp_from_sr(sfc_state%SSS(i+i0,j+j0))
867 do j=jsc_bnd,jec_bnd ;
do i=isc_bnd,iec_bnd
868 ocean_sfc%s_surf(i,j) = sfc_state%SSS(i+i0,j+j0)
872 if (
present(patm))
then
873 do j=jsc_bnd,jec_bnd ;
do i=isc_bnd,iec_bnd
874 ocean_sfc%sea_lev(i,j) = us%Z_to_m * sfc_state%sea_lev(i+i0,j+j0) + patm(i,j) * press_to_z
875 ocean_sfc%area(i,j) = us%L_to_m**2 * g%areaT(i+i0,j+j0)
878 do j=jsc_bnd,jec_bnd ;
do i=isc_bnd,iec_bnd
879 ocean_sfc%sea_lev(i,j) = us%Z_to_m * sfc_state%sea_lev(i+i0,j+j0)
880 ocean_sfc%area(i,j) = us%L_to_m**2 * g%areaT(i+i0,j+j0)
884 if (
allocated(sfc_state%frazil))
then
885 do j=jsc_bnd,jec_bnd ;
do i=isc_bnd,iec_bnd
886 ocean_sfc%frazil(i,j) = us%Q_to_J_kg*us%RZ_to_kg_m2 * sfc_state%frazil(i+i0,j+j0)
890 if (ocean_sfc%stagger == agrid)
then
891 do j=jsc_bnd,jec_bnd ;
do i=isc_bnd,iec_bnd
892 ocean_sfc%u_surf(i,j) = g%mask2dT(i+i0,j+j0) * us%L_T_to_m_s * &
893 0.5*(sfc_state%u(i+i0,j+j0)+sfc_state%u(i-1+i0,j+j0))
894 ocean_sfc%v_surf(i,j) = g%mask2dT(i+i0,j+j0) * us%L_T_to_m_s * &
895 0.5*(sfc_state%v(i+i0,j+j0)+sfc_state%v(i+i0,j-1+j0))
897 elseif (ocean_sfc%stagger == bgrid_ne)
then
898 do j=jsc_bnd,jec_bnd ;
do i=isc_bnd,iec_bnd
899 ocean_sfc%u_surf(i,j) = g%mask2dBu(i+i0,j+j0) * us%L_T_to_m_s * &
900 0.5*(sfc_state%u(i+i0,j+j0)+sfc_state%u(i+i0,j+j0+1))
901 ocean_sfc%v_surf(i,j) = g%mask2dBu(i+i0,j+j0) * us%L_T_to_m_s * &
902 0.5*(sfc_state%v(i+i0,j+j0)+sfc_state%v(i+i0+1,j+j0))
904 elseif (ocean_sfc%stagger == cgrid_ne)
then
905 do j=jsc_bnd,jec_bnd ;
do i=isc_bnd,iec_bnd
906 ocean_sfc%u_surf(i,j) = g%mask2dCu(i+i0,j+j0)*us%L_T_to_m_s * sfc_state%u(i+i0,j+j0)
907 ocean_sfc%v_surf(i,j) = g%mask2dCv(i+i0,j+j0)*us%L_T_to_m_s * sfc_state%v(i+i0,j+j0)
910 write(val_str,
'(I8)') ocean_sfc%stagger
911 call mom_error(fatal,
"convert_state_to_ocean_type: "//&
912 "Ocean_sfc%stagger has the unrecognized value of "//trim(val_str))
915 if (coupler_type_initialized(sfc_state%tr_fields))
then
916 if (.not.coupler_type_initialized(ocean_sfc%fields))
then
917 call mom_error(fatal,
"convert_state_to_ocean_type: "//&
918 "Ocean_sfc%fields has not been initialized.")
920 call coupler_type_copy_data(sfc_state%tr_fields, ocean_sfc%fields)
923 end subroutine convert_state_to_ocean_type
929 subroutine ocean_model_init_sfc(OS, Ocean_sfc)
934 integer :: is, ie, js, je
936 is = os%grid%isc ; ie = os%grid%iec ; js = os%grid%jsc ; je = os%grid%jec
937 call coupler_type_spawn(ocean_sfc%fields, os%sfc_state%tr_fields, &
938 (/is,is,ie,ie/), (/js,js,je,je/), as_needed=.true.)
940 call extract_surface_state(os%MOM_CSp, os%sfc_state)
942 call convert_state_to_ocean_type(os%sfc_state, ocean_sfc, os%grid, os%US)
944 end subroutine ocean_model_init_sfc
950 subroutine ocean_model_flux_init(OS, verbosity)
954 integer,
optional,
intent(in) :: verbosity
959 os_is_set = .false. ;
if (
present(os)) os_is_set =
associated(os)
962 verbose = 5 ;
if (os_is_set) verbose = 3
963 if (
present(verbosity)) verbose = verbosity
965 call call_tracer_flux_init(verbosity=verbose)
967 end subroutine ocean_model_flux_init
972 subroutine ocean_stock_pe(OS, index, value, time_index)
973 use stock_constants_mod,
only : istock_water, istock_heat,istock_salt
976 integer,
intent(in) :: index
977 real,
intent(out) ::
value
978 integer,
optional,
intent(in) :: time_index
988 if (.not.
associated(os))
return
989 if (.not.os%is_ocean_pe)
return
993 if (os%GV%Boussinesq)
then
994 call get_ocean_stocks(os%MOM_CSp, mass=
value, on_pe_only=.true.)
996 call get_ocean_stocks(os%MOM_CSp, mass=
value, salt=salt, on_pe_only=.true.)
1000 call get_ocean_stocks(os%MOM_CSp, heat=
value, on_pe_only=.true.)
1002 call get_ocean_stocks(os%MOM_CSp, salt=
value, on_pe_only=.true.)
1003 case default ;
value = 0.0
1009 end subroutine ocean_stock_pe
1012 subroutine ocean_model_data2d_get(OS, Ocean, name, array2D, isc, jsc)
1018 character(len=*) ,
intent(in) :: name
1019 real,
dimension(isc:,jsc:),
intent(out):: array2D
1021 integer ,
intent(in) :: isc
1022 integer ,
intent(in) :: jsc
1024 integer :: g_isc, g_iec, g_jsc, g_jec,g_isd, g_ied, g_jsd, g_jed, i, j
1026 if (.not.
associated(os))
return
1027 if (.not.os%is_ocean_pe)
return
1032 call mpp_get_compute_domain(os%grid%Domain%mpp_domain, g_isc, g_iec, g_jsc, g_jec)
1033 call mpp_get_data_domain (os%grid%Domain%mpp_domain, g_isd, g_ied, g_jsd, g_jed)
1035 g_isc = g_isc-g_isd+1 ; g_iec = g_iec-g_isd+1 ; g_jsc = g_jsc-g_jsd+1 ; g_jec = g_jec-g_jsd+1
1040 array2d(isc:,jsc:) = os%US%L_to_m**2*os%grid%areaT(g_isc:g_iec,g_jsc:g_jec)
1042 array2d(isc:,jsc:) = os%grid%mask2dT(g_isc:g_iec,g_jsc:g_jec)
1048 array2d(isc:,jsc:) = ocean%t_surf(isc:,jsc:)-celsius_kelvin_offset
1050 array2d(isc:,jsc:) = ocean%t_surf(isc:,jsc:)-celsius_kelvin_offset
1052 array2d(isc:,jsc:) = ocean%t_surf(isc:,jsc:)-celsius_kelvin_offset
1054 array2d(isc:,jsc:) = ocean%t_surf(isc:,jsc:)-celsius_kelvin_offset
1056 array2d(isc:,jsc:) = 0
1058 array2d(isc:,jsc:) = os%grid%cos_rot(g_isc:g_iec,g_jsc:g_jec)
1060 array2d(isc:,jsc:) = os%grid%sin_rot(g_isc:g_iec,g_jsc:g_jec)
1062 array2d(isc:,jsc:) = ocean%s_surf(isc:,jsc:)
1064 array2d(isc:,jsc:) = ocean%sea_lev(isc:,jsc:)
1066 array2d(isc:,jsc:) = ocean%frazil(isc:,jsc:)
1068 call mom_error(fatal,
'get_ocean_grid_data2D: unknown argument name='//name)
1071 end subroutine ocean_model_data2d_get
1074 subroutine ocean_model_data1d_get(OS, Ocean, name, value)
1079 character(len=*) ,
intent(in) :: name
1080 real ,
intent(out)::
value
1082 if (.not.
associated(os))
return
1083 if (.not.os%is_ocean_pe)
return
1089 call mom_error(fatal,
'get_ocean_grid_data1D: unknown argument name='//name)
1092 end subroutine ocean_model_data1d_get
1095 subroutine ocean_public_type_chksum(id, timestep, ocn)
1097 character(len=*),
intent(in) :: id
1098 integer,
intent(in) :: timestep
1101 integer :: n, m, outunit
1105 write(outunit,*)
"BEGIN CHECKSUM(ocean_type):: ", id, timestep
1106 write(outunit,100)
'ocean%t_surf ',mpp_chksum(ocn%t_surf )
1107 write(outunit,100)
'ocean%s_surf ',mpp_chksum(ocn%s_surf )
1108 write(outunit,100)
'ocean%u_surf ',mpp_chksum(ocn%u_surf )
1109 write(outunit,100)
'ocean%v_surf ',mpp_chksum(ocn%v_surf )
1110 write(outunit,100)
'ocean%sea_lev ',mpp_chksum(ocn%sea_lev)
1111 write(outunit,100)
'ocean%frazil ',mpp_chksum(ocn%frazil )
1113 call coupler_type_write_chksums(ocn%fields, outunit,
'ocean%')
1114 100
FORMAT(
" CHECKSUM::",a20,
" = ",z20)
1116 end subroutine ocean_public_type_chksum
1119 subroutine get_ocean_grid(OS, Gridp)
1127 end subroutine get_ocean_grid
1130 subroutine ocean_model_get_uv_surf(OS, Ocean, name, array2D, isc, jsc)
1136 character(len=*) ,
intent(in) :: name
1137 real,
dimension(isc:,jsc:),
intent(out):: array2d
1139 integer ,
intent(in) :: isc
1140 integer ,
intent(in) :: jsc
1143 type(
surface),
pointer :: sfc_state
1146 integer :: isc_bnd, iec_bnd, jsc_bnd, jec_bnd
1147 integer :: i, j, i0, j0
1148 integer :: is, ie, js, je
1150 if (.not.
associated(os))
return
1151 if (.not.os%is_ocean_pe)
return
1154 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
1156 call mpp_get_compute_domain(ocean%Domain, isc_bnd, iec_bnd, &
1159 i0 = is - isc_bnd ; j0 = js - jsc_bnd
1161 sfc_state => os%sfc_state
1165 do j=jsc_bnd,jec_bnd ;
do i=isc_bnd,iec_bnd
1166 array2d(i,j) = g%mask2dT(i+i0,j+j0) * &
1167 0.5*(sfc_state%u(i+i0,j+j0)+sfc_state%u(i-1+i0,j+j0))
1170 do j=jsc_bnd,jec_bnd ;
do i=isc_bnd,iec_bnd
1171 array2d(i,j) = g%mask2dT(i+i0,j+j0) * &
1172 0.5*(sfc_state%v(i+i0,j+j0)+sfc_state%v(i+i0,j-1+j0))
1175 do j=jsc_bnd,jec_bnd ;
do i=isc_bnd,iec_bnd
1176 array2d(i,j) = g%mask2dBu(i+i0,j+j0) * &
1177 0.5*(sfc_state%u(i+i0,j+j0)+sfc_state%u(i+i0,j+j0+1))
1180 do j=jsc_bnd,jec_bnd ;
do i=isc_bnd,iec_bnd
1181 array2d(i,j) = g%mask2dBu(i+i0,j+j0) * &
1182 0.5*(sfc_state%v(i+i0,j+j0)+sfc_state%v(i+i0+1,j+j0))
1185 call mom_error(fatal,
'ocean_model_get_UV_surf: unknown argument name='//name)
1188 end subroutine ocean_model_get_uv_surf