12 use mom_cpu_clock,
only : cpu_clock_id, cpu_clock_begin, cpu_clock_end
14 use mom_cpu_clock,
only : clock_module_driver, clock_module, clock_routine
20 use mom_domains,
only : to_south, to_west, to_all, cgrid_ne, scalar_pair
21 use mom_domains,
only : to_north, to_east, omit_corners
35 use mom_time_manager,
only :
operator(-),
operator(>),
operator(*),
operator(/)
38 use mom_barotropic,
only : barotropic_init, btstep, btcalc, bt_mass_source
66 implicit none ;
private
68 #include <MOM_memory.h>
72 real allocable_,
dimension(NIMEMB_PTR_,NJMEM_,NKMEM_) :: &
73 cau, & !< CAu = f*v - u.grad(u) [L T-2 ~> m s-2]
74 pfu, & !< PFu = -dM/dx [L T-2 ~> m s-2]
77 real allocable_,
dimension(NIMEM_,NJMEMB_PTR_,NKMEM_) :: &
78 cav, & !< CAv = -f*u - u.grad(v) [L T-2 ~> m s-2]
79 pfv, & !< PFv = -dM/dy [L T-2 ~> m s-2]
82 real allocable_,
dimension(NIMEMB_PTR_,NJMEM_,NKMEM_) :: visc_rem_u
88 real allocable_,
dimension(NIMEMB_PTR_,NJMEM_,NKMEM_) :: u_accel_bt
92 real allocable_,
dimension(NIMEM_,NJMEMB_PTR_,NKMEM_) :: visc_rem_v
98 real allocable_,
dimension(NIMEM_,NJMEMB_PTR_,NKMEM_) :: v_accel_bt
104 real allocable_,
dimension(NIMEM_,NJMEM_) :: eta
107 real allocable_,
dimension(NIMEMB_PTR_,NJMEM_,NKMEM_) :: u_av
110 real allocable_,
dimension(NIMEM_,NJMEMB_PTR_,NKMEM_) :: v_av
113 real allocable_,
dimension(NIMEM_,NJMEM_,NKMEM_) :: h_av
115 real allocable_,
dimension(NIMEM_,NJMEM_) :: eta_pf
117 real allocable_,
dimension(NIMEMB_PTR_,NJMEM_) :: uhbt
120 real allocable_,
dimension(NIMEM_,NJMEMB_PTR_) :: vhbt
123 real allocable_,
dimension(NIMEM_,NJMEM_,NKMEM_) :: pbce
127 real,
pointer,
dimension(:,:) :: taux_bot => null()
129 real,
pointer,
dimension(:,:) :: tauy_bot => null()
137 logical :: bt_use_layer_fluxes
140 logical :: split_bottom_stress
155 logical :: module_is_initialized = .false.
158 integer :: id_uh = -1, id_vh = -1
159 integer :: id_umo = -1, id_vmo = -1
160 integer :: id_umo_2d = -1, id_vmo_2d = -1
161 integer :: id_pfu = -1, id_pfv = -1
162 integer :: id_cau = -1, id_cav = -1
164 integer :: id_hf_pfu_2d = -1, id_hf_pfv_2d = -1
166 integer :: id_hf_cau_2d = -1, id_hf_cav_2d = -1
169 integer :: id_uav = -1, id_vav = -1
170 integer :: id_u_bt_accel = -1, id_v_bt_accel = -1
172 integer :: id_hf_u_bt_accel_2d = -1, id_hf_v_bt_accel_2d = -1
206 type(
ale_cs),
pointer :: ale_csp => null()
215 type(group_pass_type) :: pass_eta
216 type(group_pass_type) :: pass_visc_rem
217 type(group_pass_type) :: pass_uvp
218 type(group_pass_type) :: pass_hp_uv
219 type(group_pass_type) :: pass_uv
220 type(group_pass_type) :: pass_h
221 type(group_pass_type) :: pass_av_uvh
226 public step_mom_dyn_split_rk2
227 public register_restarts_dyn_split_rk2
228 public initialize_dyn_split_rk2
229 public end_dyn_split_rk2
232 integer :: id_clock_cor, id_clock_pres, id_clock_vertvisc
233 integer :: id_clock_horvisc, id_clock_mom_update
234 integer :: id_clock_continuity, id_clock_thick_diff
235 integer :: id_clock_btstep, id_clock_btcalc, id_clock_btforce
236 integer :: id_clock_pass, id_clock_pass_init
242 subroutine step_mom_dyn_split_rk2(u, v, h, tv, visc, Time_local, dt, forces, p_surf_begin, p_surf_end, &
243 uh, vh, uhtr, vhtr, eta_av, G, GV, US, CS, calc_dtbt, VarMix, &
244 MEKE, thickness_diffuse_CSp, Waves)
248 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
249 target,
intent(inout) :: u
250 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
251 target,
intent(inout) :: v
252 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
256 type(time_type),
intent(in) :: time_local
257 real,
intent(in) :: dt
259 real,
dimension(:,:),
pointer :: p_surf_begin
261 real,
dimension(:,:),
pointer :: p_surf_end
263 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
264 target,
intent(inout) :: uh
266 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
267 target,
intent(inout) :: vh
269 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
270 intent(inout) :: uhtr
272 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
273 intent(inout) :: vhtr
275 real,
dimension(SZI_(G),SZJ_(G)),
intent(out) :: eta_av
278 logical,
intent(in) :: calc_dtbt
287 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)) :: up
288 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)) :: vp
289 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: hp
291 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)) :: u_bc_accel
292 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)) :: v_bc_accel
296 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)),
target :: uh_in
297 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)),
target :: vh_in
301 real,
dimension(SZIB_(G),SZJ_(G)) :: uhbt_out
302 real,
dimension(SZI_(G),SZJB_(G)) :: vhbt_out
306 real,
dimension(SZI_(G),SZJ_(G)) :: eta_pred
310 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)) :: u_old_rad_obc
311 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)) :: v_old_rad_obc
317 real,
pointer,
dimension(:,:) :: &
318 p_surf => null(), eta_pf_start => null(), &
319 taux_bot => null(), tauy_bot => null(), &
322 real,
pointer,
dimension(:,:,:) :: &
323 uh_ptr => null(), u_ptr => null(), vh_ptr => null(), v_ptr => null(), &
335 real,
allocatable,
dimension(:,:) :: &
336 hf_pfu_2d, hf_pfv_2d, &
337 hf_cau_2d, hf_cav_2d, &
338 hf_u_bt_accel_2d, hf_v_bt_accel_2d
342 logical :: dyn_p_surf
343 logical :: bt_cont_bt_thick
347 logical :: showcalltree, sym
349 integer :: i, j, k, is, ie, js, je, isq, ieq, jsq, jeq, nz
350 integer :: cont_stencil
351 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
352 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
353 u_av => cs%u_av ; v_av => cs%v_av ; h_av => cs%h_av ; eta => cs%eta
355 sym=.false.;
if (g%Domain%symmetric) sym=.true.
357 showcalltree = calltree_showquery()
358 if (showcalltree)
call calltree_enter(
"step_MOM_dyn_split_RK2(), MOM_dynamics_split_RK2.F90")
362 do j=g%jsd,g%jed ;
do i=g%isdB,g%iedB ; up(i,j,k) = 0.0 ;
enddo ;
enddo
363 do j=g%jsdB,g%jedB ;
do i=g%isd,g%ied ; vp(i,j,k) = 0.0 ;
enddo ;
enddo
364 do j=g%jsd,g%jed ;
do i=g%isd,g%ied ; hp(i,j,k) = h(i,j,k) ;
enddo ;
enddo
368 call updatecfltruncationvalue(time_local, cs%vertvisc_CSp)
371 call mom_state_chksum(
"Start predictor ", u, v, h, uh, vh, g, gv, us, symmetric=sym)
376 dyn_p_surf =
associated(p_surf_begin) .and.
associated(p_surf_end)
379 call safe_alloc_ptr(eta_pf_start,g%isd,g%ied,g%jsd,g%jed)
380 eta_pf_start(:,:) = 0.0
382 p_surf => forces%p_surf
385 if (
associated(cs%OBC))
then
386 if (cs%debug_OBC)
call open_boundary_test_extern_h(g, gv, cs%OBC, h)
389 call update_obc_ramp(time_local, cs%OBC)
391 do k=1,nz ;
do j=g%jsd,g%jed ;
do i=g%IsdB,g%IedB
392 u_old_rad_obc(i,j,k) = u_av(i,j,k)
393 enddo ;
enddo ;
enddo
394 do k=1,nz ;
do j=g%JsdB,g%JedB ;
do i=g%isd,g%ied
395 v_old_rad_obc(i,j,k) = v_av(i,j,k)
396 enddo ;
enddo ;
enddo
399 bt_cont_bt_thick = .false.
400 if (
associated(cs%BT_cont)) bt_cont_bt_thick = &
401 (
allocated(cs%BT_cont%h_u) .and.
allocated(cs%BT_cont%h_v))
403 if (cs%split_bottom_stress)
then
404 taux_bot => cs%taux_bot ; tauy_bot => cs%tauy_bot
409 cont_stencil = continuity_stencil(cs%continuity_CSp)
410 call cpu_clock_begin(id_clock_pass)
412 call create_group_pass(cs%pass_visc_rem, cs%visc_rem_u, cs%visc_rem_v, g%Domain, &
413 to_all+scalar_pair, cgrid_ne, halo=max(1,cont_stencil))
423 call cpu_clock_end(id_clock_pass)
429 if (cs%begw == 0.0)
call enable_averages(dt, time_local, cs%diag)
430 call cpu_clock_begin(id_clock_pres)
431 call pressureforce(h, tv, cs%PFu, cs%PFv, g, gv, us, cs%PressureForce_CSp, &
432 cs%ALE_CSp, p_surf, cs%pbce, cs%eta_PF)
434 pres_to_eta = 1.0 / (gv%g_Earth * gv%H_to_RZ)
436 do j=jsq,jeq+1 ;
do i=isq,ieq+1
437 eta_pf_start(i,j) = cs%eta_PF(i,j) - pres_to_eta * (p_surf_begin(i,j) - p_surf_end(i,j))
440 call cpu_clock_end(id_clock_pres)
441 call disable_averaging(cs%diag)
442 if (showcalltree)
call calltree_waypoint(
"done with PressureForce (step_MOM_dyn_split_RK2)")
444 if (
associated(cs%OBC)) then;
if (cs%OBC%update_OBC)
then
445 call update_obc_data(cs%OBC, g, gv, us, tv, h, cs%update_OBC_CSp, time_local)
447 if (
associated(cs%OBC) .and. cs%debug_OBC) &
448 call open_boundary_zero_normal_flow(cs%OBC, g, cs%PFu, cs%PFv)
450 if (g%nonblocking_updates) &
451 call start_group_pass(cs%pass_eta, g%Domain, clock=id_clock_pass)
454 call cpu_clock_begin(id_clock_cor)
455 call coradcalc(u_av, v_av, h_av, uh, vh, cs%CAu, cs%CAv, cs%OBC, cs%ADp, &
456 g, gv, us, cs%CoriolisAdv_CSp)
457 call cpu_clock_end(id_clock_cor)
458 if (showcalltree)
call calltree_waypoint(
"done with CorAdCalc (step_MOM_dyn_split_RK2)")
461 call cpu_clock_begin(id_clock_btforce)
464 do j=js,je ;
do i=isq,ieq
465 u_bc_accel(i,j,k) = (cs%CAu(i,j,k) + cs%PFu(i,j,k)) + cs%diffu(i,j,k)
467 do j=jsq,jeq ;
do i=is,ie
468 v_bc_accel(i,j,k) = (cs%CAv(i,j,k) + cs%PFv(i,j,k)) + cs%diffv(i,j,k)
471 if (
associated(cs%OBC))
then
472 call open_boundary_zero_normal_flow(cs%OBC, g, u_bc_accel, v_bc_accel)
474 call cpu_clock_end(id_clock_btforce)
477 call mom_accel_chksum(
"pre-btstep accel", cs%CAu, cs%CAv, cs%PFu, cs%PFv, &
478 cs%diffu, cs%diffv, g, gv, us, cs%pbce, u_bc_accel, v_bc_accel, &
483 call check_redundant(
"pre-btstep u_bc_accel ", u_bc_accel, v_bc_accel, g)
486 call cpu_clock_begin(id_clock_vertvisc)
489 do j=js,je ;
do i=isq,ieq
490 up(i,j,k) = g%mask2dCu(i,j) * (u(i,j,k) + dt * u_bc_accel(i,j,k))
492 do j=jsq,jeq ;
do i=is,ie
493 vp(i,j,k) = g%mask2dCv(i,j) * (v(i,j,k) + dt * v_bc_accel(i,j,k))
497 call enable_averages(dt, time_local, cs%diag)
498 call set_viscous_ml(u, v, h, tv, forces, visc, dt, g, gv, us, &
500 call disable_averaging(cs%diag)
503 call uvchksum(
"before vertvisc: up", up, vp, g%HI, haloshift=0, symmetric=sym, scale=us%L_T_to_m_s)
505 call vertvisc_coef(up, vp, h, forces, visc, dt, g, gv, us, cs%vertvisc_CSp, cs%OBC)
506 call vertvisc_remnant(visc, cs%visc_rem_u, cs%visc_rem_v, dt, g, gv, us, cs%vertvisc_CSp)
507 call cpu_clock_end(id_clock_vertvisc)
508 if (showcalltree)
call calltree_waypoint(
"done with vertvisc_coef (step_MOM_dyn_split_RK2)")
511 call cpu_clock_begin(id_clock_pass)
512 if (g%nonblocking_updates)
then
513 call complete_group_pass(cs%pass_eta, g%Domain)
514 call start_group_pass(cs%pass_visc_rem, g%Domain)
516 call do_group_pass(cs%pass_eta, g%Domain)
517 call do_group_pass(cs%pass_visc_rem, g%Domain)
519 call cpu_clock_end(id_clock_pass)
521 call cpu_clock_begin(id_clock_btcalc)
523 if (.not.bt_cont_bt_thick) &
524 call btcalc(h, g, gv, cs%barotropic_CSp, obc=cs%OBC)
525 call bt_mass_source(h, eta, .true., g, gv, cs%barotropic_CSp)
526 call cpu_clock_end(id_clock_btcalc)
528 if (g%nonblocking_updates) &
529 call complete_group_pass(cs%pass_visc_rem, g%Domain, clock=id_clock_pass)
532 if (
associated(cs%BT_cont) .or. cs%BT_use_layer_fluxes)
then
533 call cpu_clock_begin(id_clock_continuity)
534 call continuity(u, v, h, hp, uh_in, vh_in, dt, g, gv, us, cs%continuity_CSp, &
535 obc=cs%OBC, visc_rem_u=cs%visc_rem_u, visc_rem_v=cs%visc_rem_v, bt_cont=cs%BT_cont)
536 call cpu_clock_end(id_clock_continuity)
537 if (bt_cont_bt_thick)
then
538 call btcalc(h, g, gv, cs%barotropic_CSp, cs%BT_cont%h_u, cs%BT_cont%h_v, &
541 if (showcalltree)
call calltree_waypoint(
"done with continuity[BT_cont] (step_MOM_dyn_split_RK2)")
544 if (cs%BT_use_layer_fluxes)
then
545 uh_ptr => uh_in ; vh_ptr => vh_in; u_ptr => u ; v_ptr => v
548 call cpu_clock_begin(id_clock_btstep)
549 if (calc_dtbt)
call set_dtbt(g, gv, us, cs%barotropic_CSp, eta, cs%pbce)
550 if (showcalltree)
call calltree_enter(
"btstep(), MOM_barotropic.F90")
552 call btstep(u, v, eta, dt, u_bc_accel, v_bc_accel, forces, cs%pbce, cs%eta_PF, &
553 u_av, v_av, cs%u_accel_bt, cs%v_accel_bt, eta_pred, cs%uhbt, cs%vhbt, &
554 g, gv, us, cs%barotropic_CSp, cs%visc_rem_u, cs%visc_rem_v, adp=cs%ADp, &
555 obc=cs%OBC, bt_cont=cs%BT_cont, eta_pf_start=eta_pf_start, &
556 taux_bot=taux_bot, tauy_bot=tauy_bot, &
557 uh0=uh_ptr, vh0=vh_ptr, u_uh0=u_ptr, v_vh0=v_ptr)
558 if (showcalltree)
call calltree_leave(
"btstep()")
559 call cpu_clock_end(id_clock_btstep)
563 call cpu_clock_begin(id_clock_mom_update)
567 do j=jsq,jeq ;
do i=is,ie
568 vp(i,j,k) = g%mask2dCv(i,j) * (v(i,j,k) + dt_pred * &
569 (v_bc_accel(i,j,k) + cs%v_accel_bt(i,j,k)))
571 do j=js,je ;
do i=isq,ieq
572 up(i,j,k) = g%mask2dCu(i,j) * (u(i,j,k) + dt_pred * &
573 (u_bc_accel(i,j,k) + cs%u_accel_bt(i,j,k)))
576 call cpu_clock_end(id_clock_mom_update)
579 call uvchksum(
"Predictor 1 [uv]", up, vp, g%HI, haloshift=0, symmetric=sym, scale=us%L_T_to_m_s)
580 call hchksum(h,
"Predictor 1 h", g%HI, haloshift=1, scale=gv%H_to_m)
581 call uvchksum(
"Predictor 1 [uv]h", uh, vh, g%HI,haloshift=2, &
582 symmetric=sym, scale=gv%H_to_m*us%L_to_m**2*us%s_to_T)
584 call mom_accel_chksum(
"Predictor accel", cs%CAu, cs%CAv, cs%PFu, cs%PFv, &
585 cs%diffu, cs%diffv, g, gv, us, cs%pbce, cs%u_accel_bt, cs%v_accel_bt, symmetric=sym)
586 call mom_state_chksum(
"Predictor 1 init", u, v, h, uh, vh, g, gv, us, haloshift=2, &
594 call cpu_clock_begin(id_clock_vertvisc)
596 call uvchksum(
"0 before vertvisc: [uv]p", up, vp, g%HI,haloshift=0, symmetric=sym, scale=us%L_T_to_m_s)
598 call vertvisc_coef(up, vp, h, forces, visc, dt_pred, g, gv, us, cs%vertvisc_CSp, &
600 call vertvisc(up, vp, h, forces, visc, dt_pred, cs%OBC, cs%ADp, cs%CDp, g, &
601 gv, us, cs%vertvisc_CSp, cs%taux_bot, cs%tauy_bot, waves=waves)
602 if (showcalltree)
call calltree_waypoint(
"done with vertvisc (step_MOM_dyn_split_RK2)")
603 if (g%nonblocking_updates)
then
604 call cpu_clock_end(id_clock_vertvisc)
605 call start_group_pass(cs%pass_uvp, g%Domain, clock=id_clock_pass)
606 call cpu_clock_begin(id_clock_vertvisc)
608 call vertvisc_remnant(visc, cs%visc_rem_u, cs%visc_rem_v, dt_pred, g, gv, us, cs%vertvisc_CSp)
609 call cpu_clock_end(id_clock_vertvisc)
611 call do_group_pass(cs%pass_visc_rem, g%Domain, clock=id_clock_pass)
612 if (g%nonblocking_updates)
then
613 call complete_group_pass(cs%pass_uvp, g%Domain, clock=id_clock_pass)
615 call do_group_pass(cs%pass_uvp, g%Domain, clock=id_clock_pass)
620 call cpu_clock_begin(id_clock_continuity)
621 call continuity(up, vp, h, hp, uh, vh, dt, g, gv, us, cs%continuity_CSp, &
622 cs%uhbt, cs%vhbt, cs%OBC, cs%visc_rem_u, cs%visc_rem_v, &
623 u_av, v_av, bt_cont=cs%BT_cont)
624 call cpu_clock_end(id_clock_continuity)
625 if (showcalltree)
call calltree_waypoint(
"done with continuity (step_MOM_dyn_split_RK2)")
627 call do_group_pass(cs%pass_hp_uv, g%Domain, clock=id_clock_pass)
629 if (
associated(cs%OBC))
then
632 call uvchksum(
"Pre OBC avg [uv]", u_av, v_av, g%HI, haloshift=1, symmetric=sym, scale=us%L_T_to_m_s)
634 call radiation_open_bdry_conds(cs%OBC, u_av, u_old_rad_obc, v_av, v_old_rad_obc, g, us, dt_pred)
637 call uvchksum(
"Post OBC avg [uv]", u_av, v_av, g%HI, haloshift=1, symmetric=sym, scale=us%L_T_to_m_s)
643 if (g%nonblocking_updates)
then
644 call start_group_pass(cs%pass_av_uvh, g%Domain, clock=id_clock_pass)
649 do k=1,nz ;
do j=js-2,je+2 ;
do i=is-2,ie+2
650 h_av(i,j,k) = 0.5*(h(i,j,k) + hp(i,j,k))
651 enddo ;
enddo ;
enddo
654 call enable_averages(dt, time_local, cs%diag)
660 call cpu_clock_begin(id_clock_btcalc)
661 call bt_mass_source(hp, eta_pred, .false., g, gv, cs%barotropic_CSp)
662 call cpu_clock_end(id_clock_btcalc)
664 if (cs%begw /= 0.0)
then
669 do k=1,nz ;
do j=js-1,je+1 ;
do i=is-1,ie+1
670 hp(i,j,k) = (1.0-cs%begw)*h(i,j,k) + cs%begw*hp(i,j,k)
671 enddo ;
enddo ;
enddo
675 call cpu_clock_begin(id_clock_pres)
676 call pressureforce(hp, tv, cs%PFu, cs%PFv, g, gv, us, cs%PressureForce_CSp, &
677 cs%ALE_CSp, p_surf, cs%pbce, cs%eta_PF)
678 call cpu_clock_end(id_clock_pres)
679 if (showcalltree)
call calltree_waypoint(
"done with PressureForce[hp=(1-b).h+b.h] (step_MOM_dyn_split_RK2)")
682 if (g%nonblocking_updates) &
683 call complete_group_pass(cs%pass_av_uvh, g%Domain, clock=id_clock_pass)
685 if (bt_cont_bt_thick)
then
686 call btcalc(h, g, gv, cs%barotropic_CSp, cs%BT_cont%h_u, cs%BT_cont%h_v, &
688 if (showcalltree)
call calltree_waypoint(
"done with btcalc[BT_cont_BT_thick] (step_MOM_dyn_split_RK2)")
692 call mom_state_chksum(
"Predictor ", up, vp, hp, uh, vh, g, gv, us, symmetric=sym)
693 call uvchksum(
"Predictor avg [uv]", u_av, v_av, g%HI, haloshift=1, symmetric=sym, scale=us%L_T_to_m_s)
694 call hchksum(h_av,
"Predictor avg h", g%HI, haloshift=0, scale=gv%H_to_m)
701 call cpu_clock_begin(id_clock_horvisc)
702 call horizontal_viscosity(u_av, v_av, h_av, cs%diffu, cs%diffv, &
703 meke, varmix, g, gv, us, cs%hor_visc_CSp, &
704 obc=cs%OBC, bt=cs%barotropic_CSp, td=thickness_diffuse_csp, &
706 call cpu_clock_end(id_clock_horvisc)
707 if (showcalltree)
call calltree_waypoint(
"done with horizontal_viscosity (step_MOM_dyn_split_RK2)")
710 call cpu_clock_begin(id_clock_cor)
711 call coradcalc(u_av, v_av, h_av, uh, vh, cs%CAu, cs%CAv, cs%OBC, cs%ADp, &
712 g, gv, us, cs%CoriolisAdv_CSp)
713 call cpu_clock_end(id_clock_cor)
714 if (showcalltree)
call calltree_waypoint(
"done with CorAdCalc (step_MOM_dyn_split_RK2)")
719 call cpu_clock_begin(id_clock_btforce)
722 do j=js,je ;
do i=isq,ieq
723 u_bc_accel(i,j,k) = (cs%Cau(i,j,k) + cs%PFu(i,j,k)) + cs%diffu(i,j,k)
725 do j=jsq,jeq ;
do i=is,ie
726 v_bc_accel(i,j,k) = (cs%Cav(i,j,k) + cs%PFv(i,j,k)) + cs%diffv(i,j,k)
729 if (
associated(cs%OBC))
then
730 call open_boundary_zero_normal_flow(cs%OBC, g, u_bc_accel, v_bc_accel)
732 call cpu_clock_end(id_clock_btforce)
735 call mom_accel_chksum(
"corr pre-btstep accel", cs%CAu, cs%CAv, cs%PFu, cs%PFv, &
736 cs%diffu, cs%diffv, g, gv, us, cs%pbce, u_bc_accel, v_bc_accel, &
740 call check_redundant(
"corr pre-btstep CS%diff ", cs%diffu, cs%diffv, g)
741 call check_redundant(
"corr pre-btstep u_bc_accel ", u_bc_accel, v_bc_accel, g)
746 call cpu_clock_begin(id_clock_btstep)
747 if (cs%BT_use_layer_fluxes)
then
748 uh_ptr => uh ; vh_ptr => vh ; u_ptr => u_av ; v_ptr => v_av
751 if (showcalltree)
call calltree_enter(
"btstep(), MOM_barotropic.F90")
753 call btstep(u, v, eta, dt, u_bc_accel, v_bc_accel, forces, cs%pbce, &
754 cs%eta_PF, u_av, v_av, cs%u_accel_bt, cs%v_accel_bt, &
755 eta_pred, cs%uhbt, cs%vhbt, g, gv, us, cs%barotropic_CSp, &
756 cs%visc_rem_u, cs%visc_rem_v, etaav=eta_av, adp=cs%ADp, &
757 obc=cs%OBC, bt_cont = cs%BT_cont, eta_pf_start=eta_pf_start, &
758 taux_bot=taux_bot, tauy_bot=tauy_bot, &
759 uh0=uh_ptr, vh0=vh_ptr, u_uh0=u_ptr, v_vh0=v_ptr)
760 do j=js,je ;
do i=is,ie ; eta(i,j) = eta_pred(i,j) ;
enddo ;
enddo
762 call cpu_clock_end(id_clock_btstep)
763 if (showcalltree)
call calltree_leave(
"btstep()")
770 call cpu_clock_begin(id_clock_mom_update)
773 do j=js,je ;
do i=isq,ieq
774 u(i,j,k) = g%mask2dCu(i,j) * (u(i,j,k) + dt * &
775 (u_bc_accel(i,j,k) + cs%u_accel_bt(i,j,k)))
777 do j=jsq,jeq ;
do i=is,ie
778 v(i,j,k) = g%mask2dCv(i,j) * (v(i,j,k) + dt * &
779 (v_bc_accel(i,j,k) + cs%v_accel_bt(i,j,k)))
782 call cpu_clock_end(id_clock_mom_update)
785 call uvchksum(
"Corrector 1 [uv]", u, v, g%HI,haloshift=0, symmetric=sym, scale=us%L_T_to_m_s)
786 call hchksum(h,
"Corrector 1 h", g%HI, haloshift=2, scale=gv%H_to_m)
787 call uvchksum(
"Corrector 1 [uv]h", uh, vh, g%HI, haloshift=2, &
788 symmetric=sym, scale=gv%H_to_m*us%L_to_m**2*us%s_to_T)
790 call mom_accel_chksum(
"Corrector accel", cs%CAu, cs%CAv, cs%PFu, cs%PFv, &
791 cs%diffu, cs%diffv, g, gv, us, cs%pbce, cs%u_accel_bt, cs%v_accel_bt, &
797 call cpu_clock_begin(id_clock_vertvisc)
798 call vertvisc_coef(u, v, h, forces, visc, dt, g, gv, us, cs%vertvisc_CSp, cs%OBC)
799 call vertvisc(u, v, h, forces, visc, dt, cs%OBC, cs%ADp, cs%CDp, g, gv, us, &
800 cs%vertvisc_CSp, cs%taux_bot, cs%tauy_bot,waves=waves)
801 if (g%nonblocking_updates)
then
802 call cpu_clock_end(id_clock_vertvisc)
803 call start_group_pass(cs%pass_uv, g%Domain, clock=id_clock_pass)
804 call cpu_clock_begin(id_clock_vertvisc)
806 call vertvisc_remnant(visc, cs%visc_rem_u, cs%visc_rem_v, dt, g, gv, us, cs%vertvisc_CSp)
807 call cpu_clock_end(id_clock_vertvisc)
808 if (showcalltree)
call calltree_waypoint(
"done with vertvisc (step_MOM_dyn_split_RK2)")
812 do k=1,nz ;
do j=js-2,je+2 ;
do i=is-2,ie+2
813 h_av(i,j,k) = h(i,j,k)
814 enddo ;
enddo ;
enddo
816 call do_group_pass(cs%pass_visc_rem, g%Domain, clock=id_clock_pass)
817 if (g%nonblocking_updates)
then
818 call complete_group_pass(cs%pass_uv, g%Domain, clock=id_clock_pass)
820 call do_group_pass(cs%pass_uv, g%Domain, clock=id_clock_pass)
826 call cpu_clock_begin(id_clock_continuity)
827 call continuity(u, v, h, h, uh, vh, dt, g, gv, us, cs%continuity_CSp, &
828 cs%uhbt, cs%vhbt, cs%OBC, cs%visc_rem_u, cs%visc_rem_v, u_av, v_av)
829 call cpu_clock_end(id_clock_continuity)
830 call do_group_pass(cs%pass_h, g%Domain, clock=id_clock_pass)
833 call diag_update_remap_grids(cs%diag)
834 if (showcalltree)
call calltree_waypoint(
"done with continuity (step_MOM_dyn_split_RK2)")
836 if (g%nonblocking_updates)
then
837 call start_group_pass(cs%pass_av_uvh, g%Domain, clock=id_clock_pass)
839 call do_group_pass(cs%pass_av_uvh, g%domain, clock=id_clock_pass)
842 if (
associated(cs%OBC))
then
843 call radiation_open_bdry_conds(cs%OBC, u, u_old_rad_obc, v, v_old_rad_obc, g, us, dt)
848 do k=1,nz ;
do j=js-2,je+2 ;
do i=is-2,ie+2
849 h_av(i,j,k) = 0.5*(h_av(i,j,k) + h(i,j,k))
850 enddo ;
enddo ;
enddo
852 if (g%nonblocking_updates) &
853 call complete_group_pass(cs%pass_av_uvh, g%Domain, clock=id_clock_pass)
857 do j=js-2,je+2 ;
do i=isq-2,ieq+2
858 uhtr(i,j,k) = uhtr(i,j,k) + uh(i,j,k)*dt
860 do j=jsq-2,jeq+2 ;
do i=is-2,ie+2
861 vhtr(i,j,k) = vhtr(i,j,k) + vh(i,j,k)*dt
870 if (cs%id_PFu > 0)
call post_data(cs%id_PFu, cs%PFu, cs%diag)
871 if (cs%id_PFv > 0)
call post_data(cs%id_PFv, cs%PFv, cs%diag)
872 if (cs%id_CAu > 0)
call post_data(cs%id_CAu, cs%CAu, cs%diag)
873 if (cs%id_CAv > 0)
call post_data(cs%id_CAv, cs%CAv, cs%diag)
876 if (cs%id_uh > 0)
call post_data(cs%id_uh , uh, cs%diag)
877 if (cs%id_vh > 0)
call post_data(cs%id_vh , vh, cs%diag)
878 if (cs%id_uav > 0)
call post_data(cs%id_uav, u_av, cs%diag)
879 if (cs%id_vav > 0)
call post_data(cs%id_vav, v_av, cs%diag)
880 if (cs%id_u_BT_accel > 0)
call post_data(cs%id_u_BT_accel, cs%u_accel_bt, cs%diag)
881 if (cs%id_v_BT_accel > 0)
call post_data(cs%id_v_BT_accel, cs%v_accel_bt, cs%diag)
901 if (cs%id_hf_PFu_2d > 0)
then
902 allocate(hf_pfu_2d(g%IsdB:g%IedB,g%jsd:g%jed))
904 do k=1,nz ;
do j=js,je ;
do i=isq,ieq
905 hf_pfu_2d(i,j) = hf_pfu_2d(i,j) + cs%PFu(i,j,k) * cs%ADp%diag_hfrac_u(i,j,k)
906 enddo ;
enddo ;
enddo
907 call post_data(cs%id_hf_PFu_2d, hf_pfu_2d, cs%diag)
908 deallocate(hf_pfu_2d)
910 if (cs%id_hf_PFv_2d > 0)
then
911 allocate(hf_pfv_2d(g%isd:g%ied,g%JsdB:g%JedB))
913 do k=1,nz ;
do j=jsq,jeq ;
do i=is,ie
914 hf_pfv_2d(i,j) = hf_pfv_2d(i,j) + cs%PFv(i,j,k) * cs%ADp%diag_hfrac_v(i,j,k)
915 enddo ;
enddo ;
enddo
916 call post_data(cs%id_hf_PFv_2d, hf_pfv_2d, cs%diag)
917 deallocate(hf_pfv_2d)
934 if (cs%id_hf_CAu_2d > 0)
then
935 allocate(hf_cau_2d(g%IsdB:g%IedB,g%jsd:g%jed))
937 do k=1,nz ;
do j=js,je ;
do i=isq,ieq
938 hf_cau_2d(i,j) = hf_cau_2d(i,j) + cs%CAu(i,j,k) * cs%ADp%diag_hfrac_u(i,j,k)
939 enddo ;
enddo ;
enddo
940 call post_data(cs%id_hf_CAu_2d, hf_cau_2d, cs%diag)
941 deallocate(hf_cau_2d)
943 if (cs%id_hf_CAv_2d > 0)
then
944 allocate(hf_cav_2d(g%isd:g%ied,g%JsdB:g%JedB))
946 do k=1,nz ;
do j=jsq,jeq ;
do i=is,ie
947 hf_cav_2d(i,j) = hf_cav_2d(i,j) + cs%CAv(i,j,k) * cs%ADp%diag_hfrac_v(i,j,k)
948 enddo ;
enddo ;
enddo
949 call post_data(cs%id_hf_CAv_2d, hf_cav_2d, cs%diag)
950 deallocate(hf_cav_2d)
967 if (cs%id_hf_u_BT_accel_2d > 0)
then
968 allocate(hf_u_bt_accel_2d(g%IsdB:g%IedB,g%jsd:g%jed))
969 hf_u_bt_accel_2d(:,:) = 0.0
970 do k=1,nz ;
do j=js,je ;
do i=isq,ieq
971 hf_u_bt_accel_2d(i,j) = hf_u_bt_accel_2d(i,j) + cs%u_accel_bt(i,j,k) * cs%ADp%diag_hfrac_u(i,j,k)
972 enddo ;
enddo ;
enddo
973 call post_data(cs%id_hf_u_BT_accel_2d, hf_u_bt_accel_2d, cs%diag)
974 deallocate(hf_u_bt_accel_2d)
976 if (cs%id_hf_v_BT_accel_2d > 0)
then
977 allocate(hf_v_bt_accel_2d(g%isd:g%ied,g%JsdB:g%JedB))
978 hf_v_bt_accel_2d(:,:) = 0.0
979 do k=1,nz ;
do j=jsq,jeq ;
do i=is,ie
980 hf_v_bt_accel_2d(i,j) = hf_v_bt_accel_2d(i,j) + cs%v_accel_bt(i,j,k) * cs%ADp%diag_hfrac_v(i,j,k)
981 enddo ;
enddo ;
enddo
982 call post_data(cs%id_hf_v_BT_accel_2d, hf_v_bt_accel_2d, cs%diag)
983 deallocate(hf_v_bt_accel_2d)
987 call mom_state_chksum(
"Corrector ", u, v, h, uh, vh, g, gv, us, symmetric=sym)
988 call uvchksum(
"Corrector avg [uv]", u_av, v_av, g%HI, haloshift=1, symmetric=sym, scale=us%L_T_to_m_s)
989 call hchksum(h_av,
"Corrector avg h", g%HI, haloshift=1, scale=gv%H_to_m)
993 if (showcalltree)
call calltree_leave(
"step_MOM_dyn_split_RK2()")
995 end subroutine step_mom_dyn_split_rk2
1000 subroutine register_restarts_dyn_split_rk2(HI, GV, param_file, CS, restart_CS, uh, vh)
1006 real,
dimension(SZIB_(HI),SZJ_(HI),SZK_(GV)), &
1007 target,
intent(inout) :: uh
1008 real,
dimension(SZI_(HI),SZJB_(HI),SZK_(GV)), &
1009 target,
intent(inout) :: vh
1012 character(len=40) :: mdl =
"MOM_dynamics_split_RK2"
1013 character(len=48) :: thickness_units, flux_units
1015 integer :: isd, ied, jsd, jed, nz, isdb, iedb, jsdb, jedb
1017 isd = hi%isd ; ied = hi%ied ; jsd = hi%jsd ; jed = hi%jed ; nz = gv%ke
1018 isdb = hi%IsdB ; iedb = hi%IedB ; jsdb = hi%JsdB ; jedb = hi%JedB
1021 if (
associated(cs))
then
1022 call mom_error(warning,
"register_restarts_dyn_split_RK2 called with an associated "// &
1023 "control structure.")
1028 alloc_(cs%diffu(isdb:iedb,jsd:jed,nz)) ; cs%diffu(:,:,:) = 0.0
1029 alloc_(cs%diffv(isd:ied,jsdb:jedb,nz)) ; cs%diffv(:,:,:) = 0.0
1030 alloc_(cs%CAu(isdb:iedb,jsd:jed,nz)) ; cs%CAu(:,:,:) = 0.0
1031 alloc_(cs%CAv(isd:ied,jsdb:jedb,nz)) ; cs%CAv(:,:,:) = 0.0
1032 alloc_(cs%PFu(isdb:iedb,jsd:jed,nz)) ; cs%PFu(:,:,:) = 0.0
1033 alloc_(cs%PFv(isd:ied,jsdb:jedb,nz)) ; cs%PFv(:,:,:) = 0.0
1035 alloc_(cs%eta(isd:ied,jsd:jed)) ; cs%eta(:,:) = 0.0
1036 alloc_(cs%u_av(isdb:iedb,jsd:jed,nz)) ; cs%u_av(:,:,:) = 0.0
1037 alloc_(cs%v_av(isd:ied,jsdb:jedb,nz)) ; cs%v_av(:,:,:) = 0.0
1038 alloc_(cs%h_av(isd:ied,jsd:jed,nz)) ; cs%h_av(:,:,:) = gv%Angstrom_H
1040 thickness_units = get_thickness_units(gv)
1041 flux_units = get_flux_units(gv)
1043 if (gv%Boussinesq)
then
1044 vd(1) = var_desc(
"sfc",thickness_units,
"Free surface Height",
'h',
'1')
1046 vd(1) = var_desc(
"p_bot",thickness_units,
"Bottom Pressure",
'h',
'1')
1050 vd(1) = var_desc(
"u2",
"m s-1",
"Auxiliary Zonal velocity",
'u',
'L')
1051 vd(2) = var_desc(
"v2",
"m s-1",
"Auxiliary Meridional velocity",
'v',
'L')
1054 vd(1) = var_desc(
"h2",thickness_units,
"Auxiliary Layer Thickness",
'h',
'L')
1057 vd(1) = var_desc(
"uh",flux_units,
"Zonal thickness flux",
'u',
'L')
1058 vd(2) = var_desc(
"vh",flux_units,
"Meridional thickness flux",
'v',
'L')
1061 vd(1) = var_desc(
"diffu",
"m s-2",
"Zonal horizontal viscous acceleration",
'u',
'L')
1062 vd(2) = var_desc(
"diffv",
"m s-2",
"Meridional horizontal viscous acceleration",
'v',
'L')
1065 call register_barotropic_restarts(hi, gv, param_file, cs%barotropic_CSp, &
1068 end subroutine register_restarts_dyn_split_rk2
1072 subroutine initialize_dyn_split_rk2(u, v, h, uh, vh, eta, Time, G, GV, US, param_file, &
1073 diag, CS, restart_CS, dt, Accel_diag, Cont_diag, MIS, &
1074 VarMix, MEKE, thickness_diffuse_CSp, &
1075 OBC, update_OBC_CSp, ALE_CSp, setVisc_CSp, &
1076 visc, dirs, ntrunc, calc_dtbt, cont_stencil)
1080 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
1082 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
1084 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) ,
intent(inout) :: h
1085 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
1086 target,
intent(inout) :: uh
1087 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
1088 target,
intent(inout) :: vh
1089 real,
dimension(SZI_(G),SZJ_(G)),
intent(inout) :: eta
1090 type(time_type),
target,
intent(in) :: time
1092 type(
diag_ctrl),
target,
intent(inout) :: diag
1095 real,
intent(in) :: dt
1109 type(
ale_cs),
pointer :: ale_csp
1113 integer,
target,
intent(inout) :: ntrunc
1116 logical,
intent(out) :: calc_dtbt
1117 integer,
optional,
intent(out) :: cont_stencil
1121 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: h_tmp
1122 character(len=40) :: mdl =
"MOM_dynamics_split_RK2"
1124 # include "version_variable.h"
1125 character(len=48) :: thickness_units, flux_units, eta_rest_name
1132 real :: accel_rescale
1135 type(group_pass_type) :: pass_av_h_uvh
1136 logical :: use_tides, debug_truncations
1138 integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz
1139 integer :: isdb, iedb, jsdb, jedb
1140 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
1141 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
1142 isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
1144 if (.not.
associated(cs))
call mom_error(fatal, &
1145 "initialize_dyn_split_RK2 called with an unassociated control structure.")
1146 if (cs%module_is_initialized)
then
1147 call mom_error(warning,
"initialize_dyn_split_RK2 called with a control "// &
1148 "structure that has already been initialized.")
1151 cs%module_is_initialized = .true.
1156 call get_param(param_file, mdl,
"TIDES", use_tides, &
1157 "If true, apply tidal momentum forcing.", default=.false.)
1158 call get_param(param_file, mdl,
"BE", cs%be, &
1159 "If SPLIT is true, BE determines the relative weighting "//&
1160 "of a 2nd-order Runga-Kutta baroclinic time stepping "//&
1161 "scheme (0.5) and a backward Euler scheme (1) that is "//&
1162 "used for the Coriolis and inertial terms. BE may be "//&
1163 "from 0.5 to 1, but instability may occur near 0.5. "//&
1164 "BE is also applicable if SPLIT is false and USE_RK2 "//&
1165 "is true.", units=
"nondim", default=0.6)
1166 call get_param(param_file, mdl,
"BEGW", cs%begw, &
1167 "If SPLIT is true, BEGW is a number from 0 to 1 that "//&
1168 "controls the extent to which the treatment of gravity "//&
1169 "waves is forward-backward (0) or simulated backward "//&
1170 "Euler (1). 0 is almost always used. "//&
1171 "If SPLIT is false and USE_RK2 is true, BEGW can be "//&
1172 "between 0 and 0.5 to damp gravity waves.", &
1173 units=
"nondim", default=0.0)
1175 call get_param(param_file, mdl,
"SPLIT_BOTTOM_STRESS", cs%split_bottom_stress, &
1176 "If true, provide the bottom stress calculated by the "//&
1177 "vertical viscosity to the barotropic solver.", default=.false.)
1178 call get_param(param_file, mdl,
"BT_USE_LAYER_FLUXES", cs%BT_use_layer_fluxes, &
1179 "If true, use the summed layered fluxes plus an "//&
1180 "adjustment due to the change in the barotropic velocity "//&
1181 "in the barotropic continuity equation.", default=.true.)
1182 call get_param(param_file, mdl,
"DEBUG", cs%debug, &
1183 "If true, write out verbose debugging data.", &
1184 default=.false., debuggingparam=.true.)
1185 call get_param(param_file, mdl,
"DEBUG_OBC", cs%debug_OBC, default=.false.)
1186 call get_param(param_file, mdl,
"DEBUG_TRUNCATIONS", debug_truncations, &
1189 allocate(cs%taux_bot(isdb:iedb,jsd:jed)) ; cs%taux_bot(:,:) = 0.0
1190 allocate(cs%tauy_bot(isd:ied,jsdb:jedb)) ; cs%tauy_bot(:,:) = 0.0
1192 alloc_(cs%uhbt(isdb:iedb,jsd:jed)) ; cs%uhbt(:,:) = 0.0
1193 alloc_(cs%vhbt(isd:ied,jsdb:jedb)) ; cs%vhbt(:,:) = 0.0
1194 alloc_(cs%visc_rem_u(isdb:iedb,jsd:jed,nz)) ; cs%visc_rem_u(:,:,:) = 0.0
1195 alloc_(cs%visc_rem_v(isd:ied,jsdb:jedb,nz)) ; cs%visc_rem_v(:,:,:) = 0.0
1196 alloc_(cs%eta_PF(isd:ied,jsd:jed)) ; cs%eta_PF(:,:) = 0.0
1197 alloc_(cs%pbce(isd:ied,jsd:jed,nz)) ; cs%pbce(:,:,:) = 0.0
1199 alloc_(cs%u_accel_bt(isdb:iedb,jsd:jed,nz)) ; cs%u_accel_bt(:,:,:) = 0.0
1200 alloc_(cs%v_accel_bt(isd:ied,jsdb:jedb,nz)) ; cs%v_accel_bt(:,:,:) = 0.0
1202 mis%diffu => cs%diffu
1203 mis%diffv => cs%diffv
1209 mis%u_accel_bt => cs%u_accel_bt
1210 mis%v_accel_bt => cs%v_accel_bt
1214 cs%ADp => accel_diag
1216 accel_diag%diffu => cs%diffu
1217 accel_diag%diffv => cs%diffv
1218 accel_diag%PFu => cs%PFu
1219 accel_diag%PFv => cs%PFv
1220 accel_diag%CAu => cs%CAu
1221 accel_diag%CAv => cs%CAv
1227 id_clock_pass_init = cpu_clock_id(
'(Ocean init message passing)', &
1228 grain=clock_routine)
1230 call continuity_init(time, g, gv, us, param_file, diag, cs%continuity_CSp)
1231 if (
present(cont_stencil)) cont_stencil = continuity_stencil(cs%continuity_CSp)
1232 call coriolisadv_init(time, g, gv, us, param_file, diag, cs%ADp, cs%CoriolisAdv_CSp)
1233 if (use_tides)
call tidal_forcing_init(time, g, param_file, cs%tides_CSp)
1234 call pressureforce_init(time, g, gv, us, param_file, diag, cs%PressureForce_CSp, &
1236 call hor_visc_init(time, g, us, param_file, diag, cs%hor_visc_CSp, meke, adp=cs%ADp)
1237 call vertvisc_init(mis, time, g, gv, us, param_file, diag, cs%ADp, dirs, &
1238 ntrunc, cs%vertvisc_CSp)
1239 if (.not.
associated(setvisc_csp))
call mom_error(fatal, &
1240 "initialize_dyn_split_RK2 called with setVisc_CSp unassociated.")
1241 cs%set_visc_CSp => setvisc_csp
1242 call updatecfltruncationvalue(time, cs%vertvisc_CSp, &
1243 activate=is_new_run(restart_cs) )
1245 if (
associated(ale_csp)) cs%ALE_CSp => ale_csp
1246 if (
associated(obc))
then
1248 if (obc%ramp)
call update_obc_ramp(time, cs%OBC, &
1249 activate=is_new_run(restart_cs) )
1251 if (
associated(update_obc_csp)) cs%update_OBC_CSp => update_obc_csp
1253 eta_rest_name =
"sfc" ;
if (.not.gv%Boussinesq) eta_rest_name =
"p_bot"
1260 if (gv%Boussinesq)
then
1261 do j=js,je ;
do i=is,ie ; cs%eta(i,j) = -gv%Z_to_H * g%bathyT(i,j) ;
enddo ;
enddo
1263 do k=1,nz ;
do j=js,je ;
do i=is,ie
1264 cs%eta(i,j) = cs%eta(i,j) + h(i,j,k)
1265 enddo ;
enddo ;
enddo
1266 elseif ((gv%m_to_H_restart /= 0.0) .and. (gv%m_to_H_restart /= gv%m_to_H))
then
1267 h_rescale = gv%m_to_H / gv%m_to_H_restart
1268 do j=js,je ;
do i=is,ie ; cs%eta(i,j) = h_rescale * cs%eta(i,j) ;
enddo ;
enddo
1271 do j=js,je ;
do i=is,ie ; eta(i,j) = cs%eta(i,j) ;
enddo ;
enddo
1273 call barotropic_init(u, v, h, cs%eta, time, g, gv, us, param_file, diag, &
1274 cs%barotropic_CSp, restart_cs, calc_dtbt, cs%BT_cont, &
1279 call horizontal_viscosity(u, v, h, cs%diffu, cs%diffv, meke, varmix, &
1280 g, gv, us, cs%hor_visc_CSp, &
1281 obc=cs%OBC, bt=cs%barotropic_CSp, &
1282 td=thickness_diffuse_csp)
1284 if ( (us%s_to_T_restart * us%m_to_L_restart /= 0.0) .and. &
1285 (us%m_to_L * us%s_to_T_restart**2 /= us%m_to_L_restart * us%s_to_T**2) )
then
1286 accel_rescale = (us%m_to_L * us%s_to_T_restart**2) / (us%m_to_L_restart * us%s_to_T**2)
1287 do k=1,nz ;
do j=js,je ;
do i=g%IscB,g%IecB
1288 cs%diffu(i,j,k) = accel_rescale * cs%diffu(i,j,k)
1289 enddo ;
enddo ;
enddo
1290 do k=1,nz ;
do j=g%JscB,g%JecB ;
do i=is,ie
1291 cs%diffv(i,j,k) = accel_rescale * cs%diffv(i,j,k)
1292 enddo ;
enddo ;
enddo
1298 do k=1,nz ;
do j=jsd,jed ;
do i=isdb,iedb ; cs%u_av(i,j,k) = u(i,j,k) ;
enddo ;
enddo ;
enddo
1299 do k=1,nz ;
do j=jsdb,jedb ;
do i=isd,ied ; cs%v_av(i,j,k) = v(i,j,k) ;
enddo ;
enddo ;
enddo
1300 elseif ( (us%s_to_T_restart * us%m_to_L_restart /= 0.0) .and. &
1301 ((us%m_to_L * us%s_to_T_restart) /= (us%m_to_L_restart * us%s_to_T)) )
then
1302 vel_rescale = (us%m_to_L * us%s_to_T_restart) / (us%m_to_L_restart * us%s_to_T)
1303 do k=1,nz ;
do j=jsd,jed ;
do i=isdb,iedb ; cs%u_av(i,j,k) = vel_rescale * cs%u_av(i,j,k) ;
enddo ;
enddo ;
enddo
1304 do k=1,nz ;
do j=jsdb,jedb ;
do i=isd,ied ; cs%v_av(i,j,k) = vel_rescale * cs%v_av(i,j,k) ;
enddo ;
enddo ;
enddo
1310 do k=1,nz ;
do j=jsd,jed ;
do i=isd,ied ; h_tmp(i,j,k) = h(i,j,k) ;
enddo ;
enddo ;
enddo
1311 call continuity(u, v, h, h_tmp, uh, vh, dt, g, gv, us, cs%continuity_CSp, obc=cs%OBC)
1312 call pass_var(h_tmp, g%Domain, clock=id_clock_pass_init)
1313 do k=1,nz ;
do j=jsd,jed ;
do i=isd,ied
1314 cs%h_av(i,j,k) = 0.5*(h(i,j,k) + h_tmp(i,j,k))
1315 enddo ;
enddo ;
enddo
1318 cs%h_av(:,:,:) = h(:,:,:)
1319 elseif ((gv%m_to_H_restart /= 0.0) .and. (gv%m_to_H_restart /= gv%m_to_H))
then
1320 h_rescale = gv%m_to_H / gv%m_to_H_restart
1321 do k=1,nz ;
do j=js,je ;
do i=is,ie ; cs%h_av(i,j,k) = h_rescale * cs%h_av(i,j,k) ;
enddo ;
enddo ;
enddo
1323 if ( (gv%m_to_H_restart * us%s_to_T_restart * us%m_to_L_restart /= 0.0) .and. &
1324 ((gv%m_to_H * us%m_to_L**2 * us%s_to_T_restart) /= &
1325 (gv%m_to_H_restart * us%m_to_L_restart**2 * us%s_to_T)) )
then
1326 uh_rescale = (gv%m_to_H * us%m_to_L**2 * us%s_to_T_restart) / &
1327 (gv%m_to_H_restart * us%m_to_L_restart**2 * us%s_to_T)
1328 do k=1,nz ;
do j=js,je ;
do i=g%IscB,g%IecB ; uh(i,j,k) = uh_rescale * uh(i,j,k) ;
enddo ;
enddo ;
enddo
1329 do k=1,nz ;
do j=g%JscB,g%JecB ;
do i=is,ie ; vh(i,j,k) = uh_rescale * vh(i,j,k) ;
enddo ;
enddo ;
enddo
1333 call cpu_clock_begin(id_clock_pass_init)
1337 call do_group_pass(pass_av_h_uvh, g%Domain)
1338 call cpu_clock_end(id_clock_pass_init)
1340 flux_units = get_flux_units(gv)
1341 h_convert = gv%H_to_m ;
if (.not.gv%Boussinesq) h_convert = gv%H_to_kg_m2
1342 cs%id_uh = register_diag_field(
'ocean_model',
'uh', diag%axesCuL, time, &
1343 'Zonal Thickness Flux', flux_units, y_cell_method=
'sum', v_extensive=.true., &
1344 conversion=h_convert*us%L_to_m**2*us%s_to_T)
1345 cs%id_vh = register_diag_field(
'ocean_model',
'vh', diag%axesCvL, time, &
1346 'Meridional Thickness Flux', flux_units, x_cell_method=
'sum', v_extensive=.true., &
1347 conversion=h_convert*us%L_to_m**2*us%s_to_T)
1349 cs%id_CAu = register_diag_field(
'ocean_model',
'CAu', diag%axesCuL, time, &
1350 'Zonal Coriolis and Advective Acceleration',
'm s-2', conversion=us%L_T2_to_m_s2)
1351 cs%id_CAv = register_diag_field(
'ocean_model',
'CAv', diag%axesCvL, time, &
1352 'Meridional Coriolis and Advective Acceleration',
'm s-2', conversion=us%L_T2_to_m_s2)
1353 cs%id_PFu = register_diag_field(
'ocean_model',
'PFu', diag%axesCuL, time, &
1354 'Zonal Pressure Force Acceleration',
'm s-2', conversion=us%L_T2_to_m_s2)
1355 cs%id_PFv = register_diag_field(
'ocean_model',
'PFv', diag%axesCvL, time, &
1356 'Meridional Pressure Force Acceleration',
'm s-2', conversion=us%L_T2_to_m_s2)
1378 cs%id_hf_PFu_2d = register_diag_field(
'ocean_model',
'hf_PFu_2d', diag%axesCu1, time, &
1379 'Depth-sum Fractional Thickness-weighted Zonal Pressure Force Acceleration',
'm s-2', &
1380 conversion=us%L_T2_to_m_s2)
1381 if(cs%id_hf_PFu_2d > 0)
call safe_alloc_ptr(cs%ADp%diag_hfrac_u,isdb,iedb,jsd,jed,nz)
1383 cs%id_hf_PFv_2d = register_diag_field(
'ocean_model',
'hf_PFv_2d', diag%axesCv1, time, &
1384 'Depth-sum Fractional Thickness-weighted Meridional Pressure Force Acceleration',
'm s-2', &
1385 conversion=us%L_T2_to_m_s2)
1386 if(cs%id_hf_PFv_2d > 0)
call safe_alloc_ptr(cs%ADp%diag_hfrac_v,isd,ied,jsd,jedb,nz)
1388 cs%id_hf_CAu_2d = register_diag_field(
'ocean_model',
'hf_CAu_2d', diag%axesCu1, time, &
1389 'Depth-sum Fractional Thickness-weighted Zonal Coriolis and Advective Acceleration',
'm s-2', &
1390 conversion=us%L_T2_to_m_s2)
1391 if(cs%id_hf_CAu_2d > 0)
call safe_alloc_ptr(cs%ADp%diag_hfrac_u,isdb,iedb,jsd,jed,nz)
1393 cs%id_hf_CAv_2d = register_diag_field(
'ocean_model',
'hf_CAv_2d', diag%axesCv1, time, &
1394 'Depth-sum Fractional Thickness-weighted Meridional Coriolis and Advective Acceleration',
'm s-2', &
1395 conversion=us%L_T2_to_m_s2)
1396 if(cs%id_hf_CAv_2d > 0)
call safe_alloc_ptr(cs%ADp%diag_hfrac_v,isd,ied,jsd,jedb,nz)
1398 cs%id_uav = register_diag_field(
'ocean_model',
'uav', diag%axesCuL, time, &
1399 'Barotropic-step Averaged Zonal Velocity',
'm s-1', conversion=us%L_T_to_m_s)
1400 cs%id_vav = register_diag_field(
'ocean_model',
'vav', diag%axesCvL, time, &
1401 'Barotropic-step Averaged Meridional Velocity',
'm s-1', conversion=us%L_T_to_m_s)
1403 cs%id_u_BT_accel = register_diag_field(
'ocean_model',
'u_BT_accel', diag%axesCuL, time, &
1404 'Barotropic Anomaly Zonal Acceleration',
'm s-2', conversion=us%L_T2_to_m_s2)
1405 cs%id_v_BT_accel = register_diag_field(
'ocean_model',
'v_BT_accel', diag%axesCvL, time, &
1406 'Barotropic Anomaly Meridional Acceleration',
'm s-2', conversion=us%L_T2_to_m_s2)
1418 cs%id_hf_u_BT_accel_2d = register_diag_field(
'ocean_model',
'hf_u_BT_accel_2d', diag%axesCu1, time, &
1419 'Depth-sum Fractional Thickness-weighted Barotropic Anomaly Zonal Acceleration',
'm s-2', &
1420 conversion=us%L_T2_to_m_s2)
1421 if(cs%id_hf_u_BT_accel_2d > 0)
call safe_alloc_ptr(cs%ADp%diag_hfrac_u,isdb,iedb,jsd,jed,nz)
1423 cs%id_hf_v_BT_accel_2d = register_diag_field(
'ocean_model',
'hf_v_BT_accel_2d', diag%axesCv1, time, &
1424 'Depth-sum Fractional Thickness-weighted Barotropic Anomaly Meridional Acceleration',
'm s-2', &
1425 conversion=us%L_T2_to_m_s2)
1426 if(cs%id_hf_v_BT_accel_2d > 0)
call safe_alloc_ptr(cs%ADp%diag_hfrac_v,isd,ied,jsd,jedb,nz)
1428 id_clock_cor = cpu_clock_id(
'(Ocean Coriolis & mom advection)', grain=clock_module)
1429 id_clock_continuity = cpu_clock_id(
'(Ocean continuity equation)', grain=clock_module)
1430 id_clock_pres = cpu_clock_id(
'(Ocean pressure force)', grain=clock_module)
1431 id_clock_vertvisc = cpu_clock_id(
'(Ocean vertical viscosity)', grain=clock_module)
1432 id_clock_horvisc = cpu_clock_id(
'(Ocean horizontal viscosity)', grain=clock_module)
1433 id_clock_mom_update = cpu_clock_id(
'(Ocean momentum increments)', grain=clock_module)
1434 id_clock_pass = cpu_clock_id(
'(Ocean message passing)', grain=clock_module)
1435 id_clock_btcalc = cpu_clock_id(
'(Ocean barotropic mode calc)', grain=clock_module)
1436 id_clock_btstep = cpu_clock_id(
'(Ocean barotropic mode stepping)', grain=clock_module)
1437 id_clock_btforce = cpu_clock_id(
'(Ocean barotropic forcing calc)', grain=clock_module)
1439 end subroutine initialize_dyn_split_rk2
1443 subroutine end_dyn_split_rk2(CS)
1446 dealloc_(cs%diffu) ; dealloc_(cs%diffv)
1447 dealloc_(cs%CAu) ; dealloc_(cs%CAv)
1448 dealloc_(cs%PFu) ; dealloc_(cs%PFv)
1450 if (
associated(cs%taux_bot))
deallocate(cs%taux_bot)
1451 if (
associated(cs%tauy_bot))
deallocate(cs%tauy_bot)
1452 dealloc_(cs%uhbt) ; dealloc_(cs%vhbt)
1453 dealloc_(cs%u_accel_bt) ; dealloc_(cs%v_accel_bt)
1454 dealloc_(cs%visc_rem_u) ; dealloc_(cs%visc_rem_v)
1456 dealloc_(cs%eta) ; dealloc_(cs%eta_PF) ; dealloc_(cs%pbce)
1457 dealloc_(cs%h_av) ; dealloc_(cs%u_av) ; dealloc_(cs%v_av)
1459 call dealloc_bt_cont_type(cs%BT_cont)
1462 end subroutine end_dyn_split_rk2