RK2 splitting for time stepping MOM adiabatic dynamics.
245 type(ocean_grid_type),
intent(inout) :: g
246 type(verticalgrid_type),
intent(in) :: gv
247 type(unit_scale_type),
intent(in) :: us
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)), &
254 type(thermo_var_ptrs),
intent(in) :: tv
255 type(vertvisc_type),
intent(inout) :: visc
256 type(time_type),
intent(in) :: time_local
257 real,
intent(in) :: dt
258 type(mech_forcing),
intent(in) :: forces
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
277 type(mom_dyn_split_rk2_cs),
pointer :: cs
278 logical,
intent(in) :: calc_dtbt
279 type(varmix_cs),
pointer :: varmix
280 type(meke_type),
pointer :: meke
281 type(thickness_diffuse_cs),
pointer :: thickness_diffuse_csp
283 type(wave_parameters_cs),
optional,
pointer :: waves
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)
372 call check_redundant(
"Start predictor u ", u, v, g)
373 call check_redundant(
"Start predictor uh ", uh, vh, g)
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)
411 call create_group_pass(cs%pass_eta, eta, g%Domain, halo=1)
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))
414 call create_group_pass(cs%pass_uvp, up, vp, g%Domain, halo=max(1,cont_stencil))
415 call create_group_pass(cs%pass_hp_uv, hp, g%Domain, halo=2)
416 call create_group_pass(cs%pass_hp_uv, u_av, v_av, g%Domain, halo=2)
417 call create_group_pass(cs%pass_hp_uv, uh(:,:,:), vh(:,:,:), g%Domain, halo=2)
419 call create_group_pass(cs%pass_uv, u, v, g%Domain, halo=max(2,cont_stencil))
420 call create_group_pass(cs%pass_h, h, g%Domain, halo=max(2,cont_stencil))
421 call create_group_pass(cs%pass_av_uvh, u_av, v_av, g%Domain, halo=2)
422 call create_group_pass(cs%pass_av_uvh, uh(:,:,:), vh(:,:,:), g%Domain, halo=2)
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, &
480 call check_redundant(
"pre-btstep CS%Ca ", cs%Cau, cs%Cav, g)
481 call check_redundant(
"pre-btstep CS%PF ", cs%PFu, cs%PFv, g)
482 call check_redundant(
"pre-btstep CS%diff ", cs%diffu, cs%diffv, g)
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, &
588 call check_redundant(
"Predictor 1 up", up, vp, g)
589 call check_redundant(
"Predictor 1 uh", uh, vh, g)
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)
696 call check_redundant(
"Predictor up ", up, vp, g)
697 call check_redundant(
"Predictor uh ", uh, vh, g)
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, &
738 call check_redundant(
"corr pre-btstep CS%Ca ", cs%Cau, cs%Cav, g)
739 call check_redundant(
"corr pre-btstep CS%PF ", cs%PFu, cs%PFv, g)
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()")
766 call check_redundant(
"u_accel_bt ", cs%u_accel_bt, cs%v_accel_bt, g)
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()")