124 use constants_mod
, only: grav, pi=>pi_8, radius, hlv, rdgas, omega, rvgas, cp_vapor
131 use fv_mp_mod, only: group_halo_update_type
132 use fv_mp_mod, only: start_group_halo_update, complete_group_halo_update
134 use diag_manager_mod
, only: send_data
136 use mpp_domains_mod
, only: dgrid_ne, cgrid_ne, mpp_update_domains, domain2d
137 use mpp_mod
, only: mpp_pe
138 use field_manager_mod
, only: model_atmos
139 use tracer_manager_mod
, only: get_tracer_index
157 real,
allocatable ::
rf(:)
163 real,
allocatable:: u00(:,:,:), v00(:,:,:)
174 subroutine fv_dynamics(npx, npy, npz, nq_tot, ng, bdt, consv_te, fill, &
175 reproduce_sum, kappa, cp_air, zvir, ptop, ks, ncnst, n_split, &
176 q_split, u, v, w, delz, hydrostatic, pt, delp, q, &
177 ps, pe, pk, peln, pkz, phis, q_con, omga, ua, va, uc, vc, &
178 ak, bk, mfx, mfy, cx, cy, ze0, hybrid_z, &
179 gridstruct, flagstruct, neststruct, idiag, bd, &
180 parent_grid, domain, diss_est, time_total)
183 use mpp_mod
, only: fatal, mpp_error
184 use ccpp_data
, only: ccpp_interstitial
187 real,
intent(IN) :: bdt
188 real,
intent(IN) :: consv_te
189 real,
intent(IN) :: kappa, cp_air
190 real,
intent(IN) :: zvir, ptop
191 real,
intent(IN),
optional :: time_total
193 integer,
intent(IN) :: npx
194 integer,
intent(IN) :: npy
195 integer,
intent(IN) :: npz
196 integer,
intent(IN) :: nq_tot
197 integer,
intent(IN) :: ng
198 integer,
intent(IN) :: ks
199 integer,
intent(IN) :: ncnst
200 integer,
intent(IN) :: n_split
201 integer,
intent(IN) :: q_split
202 logical,
intent(IN) :: fill
203 logical,
intent(IN) :: reproduce_sum
204 logical,
intent(IN) :: hydrostatic
205 logical,
intent(IN) :: hybrid_z
208 real,
intent(inout),
dimension(bd%isd:bd%ied ,bd%jsd:bd%jed+1,npz) :: u
209 real,
intent(inout),
dimension(bd%isd:bd%ied+1,bd%jsd:bd%jed ,npz) :: v
210 real,
intent(inout) :: w( bd%isd: ,bd%jsd: ,1:)
211 real,
intent(inout) :: pt( bd%isd:bd%ied ,bd%jsd:bd%jed ,npz)
212 real,
intent(inout) :: delp(bd%isd:bd%ied ,bd%jsd:bd%jed ,npz)
213 real,
intent(inout) :: q( bd%isd:bd%ied ,bd%jsd:bd%jed ,npz, ncnst)
214 real,
intent(inout) :: delz(bd%isd:,bd%jsd:,1:)
215 real,
intent(inout) :: ze0(bd%is:, bd%js: ,1:)
216 real,
intent(inout),
dimension(bd%isd:bd%ied ,bd%jsd:bd%jed,npz) :: diss_est
224 real,
intent(inout) :: ps (bd%isd:bd%ied ,bd%jsd:bd%jed)
225 real,
intent(inout) :: pe (bd%is-1:bd%ie+1, npz+1,bd%js-1:bd%je+1)
226 real,
intent(inout) :: pk (bd%is:bd%ie,bd%js:bd%je, npz+1)
227 real,
intent(inout) :: peln(bd%is:bd%ie,npz+1,bd%js:bd%je)
228 real,
intent(inout) :: pkz (bd%is:bd%ie,bd%js:bd%je,npz)
229 real,
intent(inout):: q_con(bd%isd:, bd%jsd:, 1:)
234 real,
intent(inout) :: phis(bd%isd:bd%ied,bd%jsd:bd%jed)
235 real,
intent(inout) :: omga(bd%isd:bd%ied,bd%jsd:bd%jed,npz)
236 real,
intent(inout) :: uc(bd%isd:bd%ied+1,bd%jsd:bd%jed ,npz)
237 real,
intent(inout) :: vc(bd%isd:bd%ied ,bd%jsd:bd%jed+1,npz)
239 real,
intent(inout),
dimension(bd%isd:bd%ied ,bd%jsd:bd%jed ,npz):: ua, va
240 real,
intent(in),
dimension(npz+1):: ak, bk
243 real,
intent(inout) :: mfx(bd%is:bd%ie+1, bd%js:bd%je, npz)
244 real,
intent(inout) :: mfy(bd%is:bd%ie , bd%js:bd%je+1, npz)
246 real,
intent(inout) :: cx(bd%is:bd%ie+1, bd%jsd:bd%jed, npz)
247 real,
intent(inout) :: cy(bd%isd:bd%ied ,bd%js:bd%je+1, npz)
252 type(domain2d),
intent(INOUT) :: domain
257 real :: ws(bd%is:bd%ie,bd%js:bd%je)
259 real :: te_2d(bd%is:bd%ie,bd%js:bd%je)
261 real :: teq(bd%is:bd%ie,bd%js:bd%je)
262 real :: ps2(bd%isd:bd%ied,bd%jsd:bd%jed)
263 real :: m_fac(bd%is:bd%ie,bd%js:bd%je)
265 real,
dimension(bd%is:bd%ie):: cvm
267 real,
allocatable :: dp1(:,:,:), dtdt_m(:,:,:), cappa(:,:,:)
270 real,
allocatable :: kapad(:,:,:)
272 real:: akap, rdg, ph1, ph2, mdt, gam, amdt, u0
273 real:: recip_k_split,reg_bc_update_time
274 integer :: kord_tracer(ncnst)
275 integer :: i,j,k, n, iq, n_map, nq, nwat, k_split
276 integer :: sphum, liq_wat = -999, ice_wat = -999
277 integer :: rainwat = -999, snowwat = -999, graupel = -999, cld_amt = -999
278 integer :: theta_d = -999
280 logical used, do_omega
282 logical used, last_step, do_omega
284 integer,
parameter :: max_packs=12
285 type(group_halo_update_type),
save :: i_pack(max_packs)
286 integer :: is, ie, js, je
287 integer :: isd, ied, jsd, jed
294 ccpp_associate: associate( cappa => ccpp_interstitial%cappa, &
295 dp1 => ccpp_interstitial%te0, &
296 dtdt_m => ccpp_interstitial%dtdt, &
297 last_step => ccpp_interstitial%last_step, &
298 te_2d => ccpp_interstitial%te0_2d )
314 k_split = flagstruct%k_split
315 recip_k_split=1./
real(k_split)
316 nwat = flagstruct%nwat
317 nq = nq_tot - flagstruct%dnats
324 call ccpp_interstitial%reset()
325 if (flagstruct%do_sat_adj)
then 326 ccpp_interstitial%out_dt = (idiag%id_mdt > 0)
332 allocate ( dp1(isd:ied, jsd:jed, 1:npz) )
336 allocate ( cappa(isd:ied,jsd:jed,npz) )
339 allocate ( cappa(isd:isd,jsd:jsd,1) )
345 allocate ( kapad(isd:ied, jsd:jed, npz) )
351 if (gridstruct%nested .or. any(neststruct%child_grids))
then 354 u, v, w, pt, delp, delz, q, uc, vc, pkz, &
355 neststruct%nested, flagstruct%inline_q, flagstruct%make_nh, ng, &
356 gridstruct, flagstruct, neststruct, &
357 neststruct%nest_timestep, neststruct%tracer_nest_timestep, &
361 if (gridstruct%nested)
then 364 0, 0, npx, npy, npz, bd, 1., 1., &
365 neststruct%pt_BC, bctype=neststruct%nestbctype )
368 0, 0, npx, npy, npz, bd, 1., 1., &
369 neststruct%q_con_BC, bctype=neststruct%nestbctype )
372 0, 0, npx, npy, npz, bd, 1., 1., &
373 neststruct%cappa_BC, bctype=neststruct%nestbctype )
385 if (flagstruct%regional)
then 389 call set_regional_bcs &
397 ,q,u,v,uc,vc, bd, npz, reg_bc_update_time )
402 if ( flagstruct%no_dycore )
then 403 if ( nwat.eq.2 .and. (.not.hydrostatic) )
then 404 sphum = get_tracer_index(model_atmos,
'sphum')
413 sphum = get_tracer_index(model_atmos,
'sphum')
414 liq_wat = get_tracer_index(model_atmos,
'liq_wat')
415 ice_wat = get_tracer_index(model_atmos,
'ice_wat')
416 rainwat = get_tracer_index(model_atmos,
'rainwat')
417 snowwat = get_tracer_index(model_atmos,
'snowwat')
418 graupel = get_tracer_index(model_atmos,
'graupel')
419 cld_amt = get_tracer_index(model_atmos,
'cld_amt')
422 theta_d = get_tracer_index(model_atmos,
'theta_d')
426 pfull(1) = 0.5*flagstruct%p_ref
433 ph1 = ak(k ) + bk(k )*flagstruct%p_ref
434 ph2 = ak(k+1) + bk(k+1)*flagstruct%p_ref
435 pfull(k) = (ph2 - ph1) / log(ph2/ph1)
438 if ( hydrostatic )
then 439 #if defined(CCPP) && defined(__GFORTRAN__) 448 call moist_cp(is,ie,isd,ied,jsd,jed, npz, j, k, nwat, sphum, liq_wat, rainwat, &
449 ice_wat, snowwat, graupel, q, q_con(is:ie,j,k), cvm)
452 dp1(i,j,k) = zvir*q(i,j,k,sphum)
457 #if defined(CCPP) && defined(__GFORTRAN__) 466 #if defined(CCPP) && defined(__GFORTRAN__) 473 if ( flagstruct%moist_phys )
then 476 call moist_cv(is,ie,isd,ied,jsd,jed, npz, j, k, nwat, sphum, liq_wat, rainwat, &
477 ice_wat, snowwat, graupel, q, q_con(is:ie,j,k), cvm)
481 dp1(i,j,k) =
virq(q(i,j,k,:))-1.
482 kapad(i,j,k)= kappa * (
virqd(q(i,j,k,:))/
vicpqd(q(i,j,k,:)))
484 dp1(i,j,k) = zvir*q(i,j,k,sphum)
488 cappa(i,j,k) = rdgas/(rdgas + cvm(i)/(1.+dp1(i,j,k)))
489 pkz(i,j,k) = exp(cappa(i,j,k)*log(rdg*delp(i,j,k)*pt(i,j,k)* &
491 (1.+dp1(i,j,k)) /delz(i,j,k)) )
493 (1.+dp1(i,j,k))*(1.-q_con(i,j,k))/delz(i,j,k)) )
496 pkz(i,j,k) = exp( kappa*log(rdg*delp(i,j,k)*pt(i,j,k)* &
497 (1.+dp1(i,j,k))/delz(i,j,k)) )
510 kapad(i,j,k)= kappa * (
virqd(q(i,j,k,:))/
vicpqd(q(i,j,k,:)))
511 pkz(i,j,k) = exp(kapad(i,j,k)*log(rdg*
virqd(q(i,j,k,:))*delp(i,j,k)*pt(i,j,k)/delz(i,j,k)))
513 pkz(i,j,k) = exp(kappa*log(rdg*delp(i,j,k)*pt(i,j,k)/delz(i,j,k)))
520 if ( flagstruct%fv_debug )
then 522 call prt_mxm(
'cappa', cappa, is, ie, js, je, ng, npz, 1., gridstruct%area_64, domain)
524 call prt_mxm(
'PS', ps, is, ie, js, je, ng, 1, 0.01, gridstruct%area_64, domain)
525 call prt_mxm(
'T_dyn_b', pt, is, ie, js, je, ng, npz, 1., gridstruct%area_64, domain)
526 if ( .not. hydrostatic)
call prt_mxm(
'delz', delz, is, ie, js, je, ng, npz, 1., gridstruct%area_64, domain)
527 call prt_mxm(
'delp_b ', delp, is, ie, js, je, ng, npz, 0.01, gridstruct%area_64, domain)
528 call prt_mxm(
'pk_b', pk, is, ie, js, je, 0, npz+1, 1.,gridstruct%area_64, domain)
529 call prt_mxm(
'pkz_b', pkz,is, ie, js, je, 0, npz, 1.,gridstruct%area_64, domain)
538 u, v, w, delz, pt, delp, q, dp1, pe, peln, phis, &
539 gridstruct%rsin2, gridstruct%cosa_s, &
540 zvir, cp_air, rdgas, hlv, te_2d, ua, va, teq, &
541 flagstruct%moist_phys, nwat, sphum, liq_wat, rainwat, &
542 ice_wat, snowwat, graupel, hydrostatic, idiag%id_te)
543 if( idiag%id_te>0 )
then 544 used = send_data(idiag%id_te, teq,
fv_time)
550 if( (flagstruct%consv_am.or.idiag%id_amdt>0) .and. (.not.
do_adiabatic_init) )
then 551 call compute_aam(npz, is, ie, js, je, isd, ied, jsd, jed, gridstruct, bd, &
552 ptop, ua, va, u, v, delp, teq, ps2, m_fac)
555 if( .not.flagstruct%RF_fast .and. flagstruct%tau .ne. 0. )
then 556 if ( gridstruct%grid_type<4 )
then 561 call rayleigh_super(abs(bdt), npx, npy, npz, ks, pfull, phis, flagstruct%tau, u, v, w, pt, &
562 ua, va, delz, gridstruct%agrid, cp_air, rdgas, ptop, hydrostatic, &
563 (.not. (neststruct%nested .or. flagstruct%regional)), flagstruct%rf_cutoff, gridstruct, domain, bd)
566 call rayleigh_friction(abs(bdt), npx, npy, npz, ks, pfull, flagstruct%tau, u, v, w, pt, &
567 ua, va, delz, cp_air, rdgas, ptop, hydrostatic, .true., flagstruct%rf_cutoff, gridstruct, domain, bd)
575 if ( flagstruct%adiabatic .and. flagstruct%kord_tm>0 )
then 581 pt(i,j,k) = pt(i,j,k)/pkz(i,j,k)
584 if ( theta_d>0 )
then 587 q(i,j,k,theta_d) = pt(i,j,k)
595 #if defined(CCPP) && defined(__GFORTRAN__) 605 pt(i,j,k) = pt(i,j,k)*(1.+dp1(i,j,k))/pkz(i,j,k)
608 pt(i,j,k) = pt(i,j,k)*(1.+dp1(i,j,k))*(1.-q_con(i,j,k))/pkz(i,j,k)
610 pt(i,j,k) = pt(i,j,k)*(1.+dp1(i,j,k))/pkz(i,j,k)
621 mdt = bdt /
real(k_split)
625 allocate ( dtdt_m(is:ie,js:je,npz) )
642 call start_group_halo_update(i_pack(11), q_con, domain)
644 call start_group_halo_update(i_pack(12), cappa, domain)
648 call start_group_halo_update(i_pack(13), kapad, domain)
650 call start_group_halo_update(i_pack(1), delp, domain, complete=.false.)
651 call start_group_halo_update(i_pack(1), pt, domain, complete=.true.)
653 call start_group_halo_update(i_pack(8), u, v, domain, gridtype=dgrid_ne)
656 #if defined(CCPP) && defined(__GFORTRAN__) 664 dp1(i,j,k) = delp(i,j,k)
669 if ( n_map==k_split ) last_step = .true.
673 call complete_group_halo_update(i_pack(11), domain)
675 call complete_group_halo_update(i_pack(12), domain)
680 call complete_group_halo_update(i_pack(13), domain)
684 call dyn_core(npx, npy, npz, ng, sphum, nq, mdt, n_split, zvir, cp_air, akap, cappa, &
689 u, v, w, delz, pt, q, delp, pe, pk, phis, ws, omga, ptop, pfull, ua, va, &
690 uc, vc, mfx, mfy, cx, cy, pkz, peln, q_con, ak, bk, ks, &
691 gridstruct, flagstruct, neststruct, idiag, bd, &
692 domain, n_map==1, i_pack, last_step, diss_est,time_total)
699 ps(i,j) = delp(i,j,1) *
agrav 703 if( .not. flagstruct%inline_q .and. nq /= 0 )
then 709 if (gridstruct%nested .or. flagstruct%regional)
then 710 call tracer_2d_nested(q, dp1, mfx, mfy, cx, cy, gridstruct, bd, domain, npx, npy, npz, nq, &
711 flagstruct%hord_tr, q_split, mdt, idiag%id_divg, i_pack(10), &
712 flagstruct%nord_tr, flagstruct%trdm2, &
713 k_split, neststruct, parent_grid, flagstruct%lim_fac,flagstruct%regional)
715 if ( flagstruct%z_tracer )
then 716 call tracer_2d_1l(q, dp1, mfx, mfy, cx, cy, gridstruct, bd, domain, npx, npy, npz, nq, &
717 flagstruct%hord_tr, q_split, mdt, idiag%id_divg, i_pack(10), &
718 flagstruct%nord_tr, flagstruct%trdm2, flagstruct%lim_fac,flagstruct%regional)
720 call tracer_2d(q, dp1, mfx, mfy, cx, cy, gridstruct, bd, domain, npx, npy, npz, nq, &
721 flagstruct%hord_tr, q_split, mdt, idiag%id_divg, i_pack(10), &
722 flagstruct%nord_tr, flagstruct%trdm2, flagstruct%lim_fac,flagstruct%regional)
728 if ( flagstruct%hord_tr<8 .and. flagstruct%moist_phys )
then 731 call fill2d(is, ie, js, je, ng, npz, q(isd,jsd,1,liq_wat), delp, gridstruct%area, domain, neststruct%nested, gridstruct%regional, npx, npy)
733 call fill2d(is, ie, js, je, ng, npz, q(isd,jsd,1,rainwat), delp, gridstruct%area, domain, neststruct%nested, gridstruct%regional, npx, npy)
735 call fill2d(is, ie, js, je, ng, npz, q(isd,jsd,1,ice_wat), delp, gridstruct%area, domain, neststruct%nested, gridstruct%regional, npx, npy)
737 call fill2d(is, ie, js, je, ng, npz, q(isd,jsd,1,snowwat), delp, gridstruct%area, domain, neststruct%nested, gridstruct%regional, npx, npy)
739 call fill2d(is, ie, js, je, ng, npz, q(isd,jsd,1,graupel), delp, gridstruct%area, domain, neststruct%nested, gridstruct%regional, npx, npy)
744 if( last_step .and. idiag%id_divg>0 )
then 745 used = send_data(idiag%id_divg, dp1,
fv_time)
746 if(flagstruct%fv_debug)
call prt_mxm(
'divg', dp1, is, ie, js, je, 0, npz, 1.,gridstruct%area_64, domain)
759 kord_tracer(iq) = flagstruct%kord_tr
760 if ( iq==cld_amt ) kord_tracer(iq) = 9
763 do_omega = hydrostatic .and. last_step
766 call avec_timer_start(6)
770 pkz, pk, mdt, bdt, npz, is,ie,js,je, isd,ied,jsd,jed, &
771 nq, nwat, sphum, q_con, u, v, w, delz, pt, q, phis, &
772 zvir, cp_air, akap, cappa, flagstruct%kord_mt, flagstruct%kord_wz, &
773 kord_tracer, flagstruct%kord_tm, peln, te_2d, &
774 ng, ua, va, omga, dp1, ws, fill, reproduce_sum, &
775 idiag%id_mdt>0, dtdt_m, ptop, ak, bk, pfull, gridstruct, domain, &
776 flagstruct%do_sat_adj, hydrostatic, hybrid_z, do_omega, &
780 call avec_timer_stop(6)
784 if ( neststruct%nested .and. .not. last_step)
then 786 0, 0, npx, npy, npz, bd,
real(n_map+1),
real(k_split), &
787 neststruct%cappa_BC, bctype=neststruct%nestbctype )
790 if ( flagstruct%regional .and. .not. last_step)
then 792 call regional_boundary_update(cappa,
'cappa', &
793 isd, ied, jsd, jed, npz, &
795 isd, ied, jsd, jed, &
801 if( .not. hydrostatic )
then 806 omga(i,j,k) = delp(i,j,k)/delz(i,j,k)*w(i,j,k)
814 if(flagstruct%nf_omega>0) &
815 call del2_cubed(omga, 0.18*gridstruct%da_min, gridstruct, domain, npx, npy, npz, flagstruct%nf_omega, bd)
823 #if defined(CCPP) && defined(__GFORTRAN__) 831 dtdt_m(i,j,k) = dtdt_m(i,j,k) *( 86400.0 / bdt)
836 used = send_data(idiag%id_mdt, dtdt_m,
fv_time)
838 deallocate ( dtdt_m )
843 if (cld_amt > 0)
then 844 call neg_adj3(is, ie, js, je, ng, npz, &
845 flagstruct%hydrostatic, &
847 pt, delp, q(isd,jsd,1,sphum), &
848 q(isd,jsd,1,liq_wat), &
849 q(isd,jsd,1,rainwat), &
850 q(isd,jsd,1,ice_wat), &
851 q(isd,jsd,1,snowwat), &
852 q(isd,jsd,1,graupel), &
853 q(isd,jsd,1,cld_amt), flagstruct%check_negative)
855 call neg_adj3(is, ie, js, je, ng, npz, &
856 flagstruct%hydrostatic, &
858 pt, delp, q(isd,jsd,1,sphum), &
859 q(isd,jsd,1,liq_wat), &
860 q(isd,jsd,1,rainwat), &
861 q(isd,jsd,1,ice_wat), &
862 q(isd,jsd,1,snowwat), &
863 q(isd,jsd,1,graupel), check_negative=flagstruct%check_negative)
865 if ( flagstruct%fv_debug )
then 866 call prt_mxm(
'T_dyn_a3', pt, is, ie, js, je, ng, npz, 1., gridstruct%area_64, domain)
867 call prt_mxm(
'SPHUM_dyn', q(isd,jsd,1,sphum ), is, ie, js, je, ng, npz, 1.,gridstruct%area_64, domain)
868 call prt_mxm(
'liq_wat_dyn', q(isd,jsd,1,liq_wat), is, ie, js, je, ng, npz, 1.,gridstruct%area_64, domain)
869 call prt_mxm(
'rainwat_dyn', q(isd,jsd,1,rainwat), is, ie, js, je, ng, npz, 1.,gridstruct%area_64, domain)
870 call prt_mxm(
'ice_wat_dyn', q(isd,jsd,1,ice_wat), is, ie, js, je, ng, npz, 1.,gridstruct%area_64, domain)
871 call prt_mxm(
'snowwat_dyn', q(isd,jsd,1,snowwat), is, ie, js, je, ng, npz, 1.,gridstruct%area_64, domain)
872 call prt_mxm(
'graupel_dyn', q(isd,jsd,1,graupel), is, ie, js, je, ng, npz, 1.,gridstruct%area_64, domain)
877 if (cld_amt > 0)
then 878 call neg_adj2(is, ie, js, je, ng, npz, &
879 flagstruct%hydrostatic, &
881 pt, delp, q(isd,jsd,1,sphum), &
882 q(isd,jsd,1,liq_wat), &
883 q(isd,jsd,1,rainwat), &
884 q(isd,jsd,1,ice_wat), &
885 q(isd,jsd,1,snowwat), &
886 q(isd,jsd,1,cld_amt), flagstruct%check_negative)
888 call neg_adj2(is, ie, js, je, ng, npz, &
889 flagstruct%hydrostatic, &
891 pt, delp, q(isd,jsd,1,sphum), &
892 q(isd,jsd,1,liq_wat), &
893 q(isd,jsd,1,rainwat), &
894 q(isd,jsd,1,ice_wat), &
895 q(isd,jsd,1,snowwat), &
896 check_negative=flagstruct%check_negative)
898 if ( flagstruct%fv_debug )
then 899 call prt_mxm(
'T_dyn_a3', pt, is, ie, js, je, ng, npz, 1., gridstruct%area_64, domain)
900 call prt_mxm(
'SPHUM_dyn', q(isd,jsd,1,sphum ), is, ie, js, je, ng, npz, 1.,gridstruct%area_64, domain)
901 call prt_mxm(
'liq_wat_dyn', q(isd,jsd,1,liq_wat), is, ie, js, je, ng, npz, 1.,gridstruct%area_64, domain)
902 call prt_mxm(
'rainwat_dyn', q(isd,jsd,1,rainwat), is, ie, js, je, ng, npz, 1.,gridstruct%area_64, domain)
903 call prt_mxm(
'ice_wat_dyn', q(isd,jsd,1,ice_wat), is, ie, js, je, ng, npz, 1.,gridstruct%area_64, domain)
904 call prt_mxm(
'snowwat_dyn', q(isd,jsd,1,snowwat), is, ie, js, je, ng, npz, 1.,gridstruct%area_64, domain)
908 if( (flagstruct%consv_am.or.idiag%id_amdt>0.or.idiag%id_aam>0) .and. (.not.
do_adiabatic_init) )
then 909 call compute_aam(npz, is, ie, js, je, isd, ied, jsd, jed, gridstruct, bd, &
910 ptop, ua, va, u, v, delp, te_2d, ps, m_fac)
911 if( idiag%id_aam>0 )
then 912 used = send_data(idiag%id_aam, te_2d,
fv_time)
914 gam =
g_sum( domain, te_2d, is, ie, js, je, ng, gridstruct%area_64, 0)
915 if( is_master() )
write(6,*)
'Total AAM =', gam
920 if( (flagstruct%consv_am.or.idiag%id_amdt>0) .and. (.not.
do_adiabatic_init) )
then 921 #if defined(CCPP) && defined(__GFORTRAN__) 930 te_2d(i,j) = te_2d(i,j)-teq(i,j) + dt2*(ps2(i,j)+ps(i,j))*idiag%zxg(i,j)
933 if( idiag%id_amdt>0 ) used = send_data(idiag%id_amdt, te_2d/bdt,
fv_time)
935 if ( flagstruct%consv_am .or.
prt_minmax )
then 936 amdt =
g_sum( domain, te_2d, is, ie, js, je, ng, gridstruct%area_64, 0, reproduce=.true.)
937 u0 = -radius*amdt/
g_sum( domain, m_fac, is, ie, js, je, ng, gridstruct%area_64, 0,reproduce=.true.)
939 write(6,*)
'Dynamic AM tendency (Hadleys)=', amdt/(bdt*1.e18),
'del-u (per day)=', u0*86400./bdt
942 if( flagstruct%consv_am )
then 946 m_fac(i,j) = u0*cos(gridstruct%agrid(i,j,2))
954 u(i,j,k) = u(i,j,k) + u0*gridstruct%l2c_u(i,j)
959 v(i,j,k) = v(i,j,k) + u0*gridstruct%l2c_v(i,j)
967 npx, npy, npz, 1, gridstruct%grid_type, domain, gridstruct%nested, flagstruct%c2l_ord, bd)
977 if ( flagstruct%fv_debug )
then 978 call prt_mxm(
'UA', ua, is, ie, js, je, ng, npz, 1., gridstruct%area_64, domain)
979 call prt_mxm(
'VA', va, is, ie, js, je, ng, npz, 1., gridstruct%area_64, domain)
980 call prt_mxm(
'TA', pt, is, ie, js, je, ng, npz, 1., gridstruct%area_64, domain)
981 if (.not. hydrostatic)
call prt_mxm(
'W ', w, is, ie, js, je, ng, npz, 1., gridstruct%area_64, domain)
984 if ( flagstruct%range_warn )
then 985 call range_check(
'UA_dyn', ua, is, ie, js, je, ng, npz, gridstruct%agrid, &
987 call range_check(
'VA_dyn', ua, is, ie, js, je, ng, npz, gridstruct%agrid, &
989 call range_check(
'TA_dyn', pt, is, ie, js, je, ng, npz, gridstruct%agrid, &
991 if ( .not. hydrostatic ) &
992 call range_check(
'W_dyn', w, is, ie, js, je, ng, npz, gridstruct%agrid, &
997 end associate ccpp_associate
1003 subroutine rayleigh_fast(dt, npx, npy, npz, pfull, tau, u, v, w, &
1004 ks, dp, ptop, hydrostatic, rf_cutoff, bd)
1006 real,
intent(in):: dt
1007 real,
intent(in):: tau
1008 real,
intent(in):: ptop, rf_cutoff
1009 real,
intent(in),
dimension(npz):: pfull
1010 integer,
intent(in):: npx, npy, npz, ks
1011 logical,
intent(in):: hydrostatic
1012 type(fv_grid_bounds_type),
intent(IN) :: bd
1013 real,
intent(inout):: u(bd%isd:bd%ied ,bd%jsd:bd%jed+1,npz)
1014 real,
intent(inout):: v(bd%isd:bd%ied+1,bd%jsd:bd%jed,npz)
1015 real,
intent(inout):: w(bd%isd: ,bd%jsd: ,1: )
1016 real,
intent(in):: dp(npz)
1018 real(kind=R_GRID):: rff(npz)
1019 real,
parameter:: sday = 86400.
1020 real,
dimension(bd%is:bd%ie+1):: dmv
1021 real,
dimension(bd%is:bd%ie):: dmu
1025 integer :: is, ie, js, je
1026 integer :: isd, ied, jsd, jed
1042 if( is_master() )
write(6,*)
'Fast Rayleigh friction E-folding time (days):' 1044 if ( pfull(k) < rf_cutoff )
then 1045 rff(k) = dt/tau0*sin(0.5*pi*log(rf_cutoff/pfull(k))/log(rf_cutoff/ptop))**2
1047 if( is_master() )
write(6,*) k, 0.01*pfull(k), dt/(rff(k)*sday)
1049 rff(k) = 1.d0 / (1.0d0+rff(k))
1057 if ( pfull(k) < 100.e2 )
then 1064 if( is_master() )
write(6,*)
'k_rf=',
k_rf, 0.01*pfull(
k_rf),
'dm=', dm
1081 dmu(i) = dmu(i) + (1.-
rf(k))*dp(k)*u(i,j,k)
1082 u(i,j,k) =
rf(k)*u(i,j,k)
1086 dmv(i) = dmv(i) + (1.-
rf(k))*dp(k)*v(i,j,k)
1087 v(i,j,k) =
rf(k)*v(i,j,k)
1089 if ( .not. hydrostatic )
then 1091 w(i,j,k) =
rf(k)*w(i,j,k)
1098 dmu(i) = dmu(i) / dm
1102 dmv(i) = dmv(i) / dm
1108 u(i,j,k) = u(i,j,k) + dmu(i)
1112 v(i,j,k) = v(i,j,k) + dmv(i)
1119 end subroutine rayleigh_fast
1124 subroutine rayleigh_super(dt, npx, npy, npz, ks, pm, phis, tau, u, v, w, pt, &
1125 ua, va, delz, agrid, cp, rg, ptop, hydrostatic, &
1126 conserve, rf_cutoff, gridstruct, domain, bd)
1127 real,
intent(in):: dt
1128 real,
intent(in):: tau
1129 real,
intent(in):: cp, rg, ptop, rf_cutoff
1130 real,
intent(in),
dimension(npz):: pm
1131 integer,
intent(in):: npx, npy, npz, ks
1132 logical,
intent(in):: hydrostatic
1133 logical,
intent(in):: conserve
1134 type(fv_grid_bounds_type),
intent(IN) :: bd
1135 real,
intent(inout):: u(bd%isd:bd%ied ,bd%jsd:bd%jed+1,npz)
1136 real,
intent(inout):: v(bd%isd:bd%ied+1,bd%jsd:bd%jed,npz)
1137 real,
intent(inout):: w(bd%isd: ,bd%jsd: ,1: )
1138 real,
intent(inout):: pt(bd%isd:bd%ied,bd%jsd:bd%jed,npz)
1139 real,
intent(inout):: ua(bd%isd:bd%ied,bd%jsd:bd%jed,npz)
1140 real,
intent(inout):: va(bd%isd:bd%ied,bd%jsd:bd%jed,npz)
1141 real,
intent(inout):: delz(bd%isd: ,bd%jsd: ,1: )
1142 real,
intent(in) :: agrid(bd%isd:bd%ied, bd%jsd:bd%jed,2)
1143 real,
intent(in) :: phis(bd%isd:bd%ied,bd%jsd:bd%jed)
1144 type(fv_grid_type),
intent(IN) :: gridstruct
1145 type(domain2d),
intent(INOUT) :: domain
1147 real,
allocatable :: u2f(:,:,:)
1148 real,
parameter:: u0 = 60.
1149 real,
parameter:: sday = 86400.
1153 integer :: is, ie, js, je
1154 integer :: isd, ied, jsd, jed
1165 rcv = 1. / (cp - rg)
1169 allocate ( u00(is:ie, js:je+1,npz) )
1170 allocate ( v00(is:ie+1,js:je ,npz) )
1175 u00(i,j,k) = u(i,j,k)
1180 v00(i,j,k) = v(i,j,k)
1188 tau0 = abs( tau * sday )
1190 if( is_master() )
write(6,*)
'Rayleigh_Super tau=',tau
1195 if( is_master() )
write(6,*) k, 0.01*pm(k)
1197 if( is_master() )
write(6,*)
'Rayleigh friction E-folding time (days):' 1199 if ( pm(k) < rf_cutoff )
then 1201 rf(k) = dt/tau0*( log(rf_cutoff/pm(k)) )**2
1203 rf(k) = dt/tau0*sin(0.5*pi*log(rf_cutoff/pm(k))/log(rf_cutoff/ptop))**2
1205 if( is_master() )
write(6,*) k, 0.01*pm(k), dt/(
rf(k)*sday)
1214 call c2l_ord2(u, v, ua, va, gridstruct, npz, gridstruct%grid_type, bd, gridstruct%nested)
1216 allocate( u2f(isd:ied,jsd:jed,
kmax) )
1221 if ( pm(k) < rf_cutoff )
then 1222 u2f(:,:,k) = 1. / (1.+
rf(k))
1227 call timing_on(
'COMM_TOTAL')
1228 call mpp_update_domains(u2f, domain)
1229 call timing_off(
'COMM_TOTAL')
1237 if ( pm(k) < rf_cutoff )
then 1239 if (.not. hydrostatic)
then 1242 w(i,j,k) = w(i,j,k)/(1.+
rf(k))
1248 u(i,j,k) = (u(i,j,k)+
rf(k)*u00(i,j,k))/(1.+
rf(k))
1253 v(i,j,k) = (v(i,j,k)+
rf(k)*v00(i,j,k))/(1.+
rf(k))
1258 if ( conserve )
then 1259 if ( hydrostatic )
then 1262 pt(i,j,k) = pt(i,j,k) + 0.5*(ua(i,j,k)**2+va(i,j,k)**2)*(1.-u2f(i,j,k)**2)/(cp-rg*ptop/pm(k))
1268 pt(i,j,k) = pt(i,j,k) + 0.5*(ua(i,j,k)**2+va(i,j,k)**2+w(i,j,k)**2)*(1.-u2f(i,j,k)**2)*rcv
1276 u(i,j,k) = 0.5*(u2f(i,j-1,k)+u2f(i,j,k))*u(i,j,k)
1281 v(i,j,k) = 0.5*(u2f(i-1,j,k)+u2f(i,j,k))*v(i,j,k)
1284 if ( .not. hydrostatic )
then 1287 w(i,j,k) = u2f(i,j,k)*w(i,j,k)
1300 subroutine rayleigh_friction(dt, npx, npy, npz, ks, pm, tau, u, v, w, pt, &
1301 ua, va, delz, cp, rg, ptop, hydrostatic, conserve, &
1302 rf_cutoff, gridstruct, domain, bd)
1303 real,
intent(in):: dt
1304 real,
intent(in):: tau
1305 real,
intent(in):: cp, rg, ptop, rf_cutoff
1306 real,
intent(in),
dimension(npz):: pm
1307 integer,
intent(in):: npx, npy, npz, ks
1308 logical,
intent(in):: hydrostatic
1309 logical,
intent(in):: conserve
1310 type(fv_grid_bounds_type),
intent(IN) :: bd
1311 real,
intent(inout):: u(bd%isd:bd%ied ,bd%jsd:bd%jed+1,npz)
1312 real,
intent(inout):: v(bd%isd:bd%ied+1,bd%jsd:bd%jed,npz)
1313 real,
intent(inout):: w(bd%isd: ,bd%jsd: ,1: )
1314 real,
intent(inout):: pt(bd%isd:bd%ied,bd%jsd:bd%jed,npz)
1315 real,
intent(inout):: ua(bd%isd:bd%ied,bd%jsd:bd%jed,npz)
1316 real,
intent(inout):: va(bd%isd:bd%ied,bd%jsd:bd%jed,npz)
1317 real,
intent(inout):: delz(bd%isd: ,bd%jsd: ,1: )
1318 type(fv_grid_type),
intent(IN) :: gridstruct
1319 type(domain2d),
intent(INOUT) :: domain
1321 real,
allocatable :: u2f(:,:,:)
1322 real,
parameter:: sday = 86400.
1323 real,
parameter:: u000 = 4900.
1327 integer :: is, ie, js, je
1328 integer :: isd, ied, jsd, jed
1341 rcv = 1. / (cp - rg)
1345 if( is_master() )
write(6,*)
'Rayleigh friction E-folding time (days):' 1347 if ( pm(k) < rf_cutoff )
then 1348 rf(k) = dt/(tau*sday)*sin(0.5*pi*log(rf_cutoff/pm(k))/log(rf_cutoff/ptop))**2
1349 if( is_master() )
write(6,*) k, 0.01*pm(k), dt/(
rf(k)*sday)
1358 allocate( u2f(isd:ied,jsd:jed,
kmax) )
1360 call c2l_ord2(u, v, ua, va, gridstruct, npz, gridstruct%grid_type, bd, gridstruct%nested)
1364 if ( hydrostatic )
then 1367 u2f(i,j,k) = ua(i,j,k)**2 + va(i,j,k)**2
1373 u2f(i,j,k) = ua(i,j,k)**2 + va(i,j,k)**2 + w(i,j,k)**2
1379 call timing_on(
'COMM_TOTAL')
1380 call mpp_update_domains(u2f, domain)
1381 call timing_off(
'COMM_TOTAL')
1387 if ( conserve )
then 1388 if ( hydrostatic )
then 1391 pt(i,j,k) = pt(i,j,k) + 0.5*u2f(i,j,k)/(cp-rg*ptop/pm(k)) &
1392 * ( 1. - 1./(1.+
rf(k)*sqrt(u2f(i,j,k)/u000))**2 )
1398 delz(i,j,k) = delz(i,j,k) / pt(i,j,k)
1399 pt(i,j,k) = pt(i,j,k) + 0.5*u2f(i,j,k) * rcv &
1400 * ( 1. - 1./(1.+
rf(k)*sqrt(u2f(i,j,k)/u000))**2 )
1401 delz(i,j,k) = delz(i,j,k) * pt(i,j,k)
1409 u2f(i,j,k) =
rf(k)*sqrt(u2f(i,j,k)/u000)
1415 u(i,j,k) = u(i,j,k) / (1.+0.5*(u2f(i,j-1,k)+u2f(i,j,k)))
1420 v(i,j,k) = v(i,j,k) / (1.+0.5*(u2f(i-1,j,k)+u2f(i,j,k)))
1424 if ( .not. hydrostatic )
then 1427 w(i,j,k) = w(i,j,k) / (1.+u2f(i,j,k))
1439 subroutine compute_aam(npz, is, ie, js, je, isd, ied, jsd, jed, gridstruct, bd, &
1440 ptop, ua, va, u, v, delp, aam, ps, m_fac)
1442 integer,
intent(in):: npz
1443 integer,
intent(in):: is, ie, js, je
1444 integer,
intent(in):: isd, ied, jsd, jed
1445 real,
intent(in):: ptop
1446 real,
intent(inout):: u(isd:ied ,jsd:jed+1,npz)
1447 real,
intent(inout):: v(isd:ied+1,jsd:jed,npz)
1448 real,
intent(inout):: delp(isd:ied,jsd:jed,npz)
1449 real,
intent(inout),
dimension(isd:ied,jsd:jed, npz):: ua, va
1450 real,
intent(out):: aam(is:ie,js:je)
1451 real,
intent(out):: m_fac(is:ie,js:je)
1452 real,
intent(out):: ps(isd:ied,jsd:jed)
1453 type(fv_grid_bounds_type),
intent(IN) :: bd
1454 type(fv_grid_type),
intent(IN) :: gridstruct
1456 real,
dimension(is:ie):: r1, r2, dm
1459 call c2l_ord2(u, v, ua, va, gridstruct, npz, gridstruct%grid_type, bd, gridstruct%nested)
1465 r1(i) = radius*cos(gridstruct%agrid(i,j,2))
1474 ps(i,j) = ps(i,j) + dm(i)
1476 aam(i,j) = aam(i,j) + (r2(i)*omega + r1(i)*ua(i,j,k)) * dm(i)
1477 m_fac(i,j) = m_fac(i,j) + dm(i)*r2(i)
subroutine, public init_ijk_mem(i1, i2, j1, j2, km, array, var)
subroutine, public fill2d(is, ie, js, je, ng, km, q, delp, area, domain, nested, regional, npx, npy)
The subroutine 'fill2D' fills in nonphysical negative values in a single scalar field using a two-dim...
logical, public do_adiabatic_init
pure real function, public vicpqd(q)
subroutine, public setup_nested_grid_bcs(npx, npy, npz, zvir, ncnst, u, v, w, pt, delp, delz, q, uc, vc, pkz, nested, inline_q, make_nh, ng, gridstruct, flagstruct, neststruct, nest_timestep, tracer_nest_timestep, domain, bd, nwat)
The subroutine 'setup_nested_grid_BCs' fetches data from the coarse grid to set up the nested-grid bo...
subroutine, public regional_boundary_update(array, bc_vbl_name, lbnd_x, ubnd_x, lbnd_y, ubnd_y, ubnd_z, is, ie, js, je, isd, ied, jsd, jed, fcst_time, index4)
real, public current_time_in_seconds
The module 'fv_mp_mod' is a single program multiple data (SPMD) parallel decompostion/communication m...
subroutine timing_off(blk_name)
The subroutine 'timing_off' stops a timer.
The type 'fv_grid_type' is made up of grid-dependent information from fv_grid_tools and fv_grid_utils...
subroutine, public neg_adj3(is, ie, js, je, ng, kbot, hydrostatic, peln, delz, pt, dp, qv, ql, qr, qi, qs, qg, qa, check_negative)
integer, parameter, public h_stagger
subroutine, public moist_cv(is, ie, isd, ied, jsd, jed, km, j, k, nwat, sphum, liq_wat, rainwat, ice_wat, snowwat, graupel, q, qd, cvm, t1)
The subroutine 'moist_cv' computes the FV3-consistent moist heat capacity under constant volume...
The module 'multi_gases' peforms multi constitutents computations.
subroutine compute_aam(npz, is, ie, js, je, isd, ied, jsd, jed, gridstruct, bd, ptop, ua, va, u, v, delp, aam, ps, m_fac)
The subroutine 'compute_aam' computes vertically (mass) integrated Atmospheric Angular Momentum...
real function, public g_sum(domain, p, ifirst, ilast, jfirst, jlast, ngc, area, mode, reproduce)
The function 'g_sum' is the fast version of 'globalsum'.
subroutine, public compute_total_energy(is, ie, js, je, isd, ied, jsd, jed, km, u, v, w, delz, pt, delp, q, qc, pe, peln, hs, rsin2_l, cosa_s_l, r_vir, cp, rg, hlv, te_2d, ua, va, teq, moist_phys, nwat, sphum, liq_wat, rainwat, ice_wat, snowwat, graupel, hydrostatic, id_te)
The subroutine 'compute_total_energy' performs the FV3-consistent computation of the global total ene...
The module 'fv_dynamics' is the top-level routine for the dynamical core.
The module 'fv_sg' performs FV sub-grid mixing.
pure real function, public virqd(q)
The module fv_nwp_nudge contains routines for nudging to input analyses. note This module is currentl...
subroutine, public prt_mxm(qname, q, is, ie, js, je, n_g, km, fac, area, domain)
The module 'fv_timing' contains FV3 timers.
pure real function, public virq(q)
The module 'boundary' contains utility routines for grid nesting and boundary conditions.
subroutine, public lagrangian_to_eulerian(last_step, consv, ps, pe, delp, pkz, pk, mdt, pdt, km, is, ie, js, je, isd, ied, jsd, jed, nq, nwat, sphum, q_con, u, v, w, delz, pt, q, hs, r_vir, cp, akap, cappa, kord_mt, kord_wz, kord_tr, kord_tm, peln, te0_2d, ng, ua, va, omga, te, ws, fill, reproduce_sum, out_dt, dtdt, ptop, ak, bk, pfull, gridstruct, domain, do_sat_adj, hydrostatic, hybrid_z, do_omega, adiabatic, do_adiabatic_init)
The subroutine 'Lagrangian_to_Eulerian' remaps deformed Lagrangian layers back to the reference Euler...
subroutine, public fv_dynamics(npx, npy, npz, nq_tot, ng, bdt, consv_te, fill, reproduce_sum, kappa, cp_air, zvir, ptop, ks, ncnst, n_split, q_split, u, v, w, delz, hydrostatic, pt, delp, q, ps, pe, pk, peln, pkz, phis, q_con, omga, ua, va, uc, vc, ak, bk, mfx, mfy, cx, cy, ze0, hybrid_z, gridstruct, flagstruct, neststruct, idiag, bd, parent_grid, domain, diss_est, time_total)
The module 'fv_tracer2d.F90' performs sub-cycled tracer advection.
subroutine, public set_regional_bcs(delp, delz, w, pt ifdef USE_COND
The module 'fv_mapz' contains the vertical mapping routines .
subroutine, public dyn_core(npx, npy, npz, ng, sphum, nq, bdt, n_split, zvir, cp, akap, cappa, grav, hydrostatic, u, v, w, delz, pt, q, delp, pe, pk, phis, ws, omga, ptop, pfull, ua, va, uc, vc, mfx, mfy, cx, cy, pkz, peln, q_con, ak, bk, ks, gridstruct, flagstruct, neststruct, idiag, bd, domain, init_step, i_pack, end_step, diss_est, time_total)
logical, public prt_minmax
The module 'fv_arrays' contains the 'fv_atmos_type' and associated datatypes.
subroutine, public tracer_2d(q, dp1, mfx, mfy, cx, cy, gridstruct, bd, domain, npx, npy, npz, nq, hord, q_split, dt, id_divg, q_pack, nord_tr, trdm, lim_fac, regional)
The subroutine 'tracer_2d' is the standard routine for sub-cycled tracer advection.
subroutine, public moist_cp(is, ie, isd, ied, jsd, jed, km, j, k, nwat, sphum, liq_wat, rainwat, ice_wat, snowwat, graupel, q, qd, cpm, t1)
The subroutine 'moist_cp' computes the FV3-consistent moist heat capacity under constant pressure...
subroutine, public del2_cubed(q, cd, gridstruct, domain, npx, npy, km, nmax, bd)
The subroutine 'del2-cubed' filters the omega field for the physics.
type(time_type), public fv_time
subroutine rayleigh_friction(dt, npx, npy, npz, ks, pm, tau, u, v, w, pt, ua, va, delz, cp, rg, ptop, hydrostatic, conserve, rf_cutoff, gridstruct, domain, bd)
subroutine timing_on(blk_name)
The subroutine 'timing_on' starts a timer.
subroutine, public range_check(qname, q, is, ie, js, je, n_g, km, pos, q_low, q_hi, bad_range)
integer, parameter, public v_stagger
subroutine, public c2l_ord2(u, v, ua, va, gridstruct, km, grid_type, bd, do_halo)
@ The module 'fv_diagnostics' contains routines to compute diagnosic fields.
The module 'fv_grid_utils' contains routines for setting up and computing grid-related quantities...
The module 'dyn_core' peforms the Lagrangian acoustic dynamics described by .
subroutine rayleigh_super(dt, npx, npy, npz, ks, pm, phis, tau, u, v, w, pt, ua, va, delz, agrid, cp, rg, ptop, hydrostatic, conserve, rf_cutoff, gridstruct, domain, bd)
subroutine, public tracer_2d_1l(q, dp1, mfx, mfy, cx, cy, gridstruct, bd, domain, npx, npy, npz, nq, hord, q_split, dt, id_divg, q_pack, nord_tr, trdm, lim_fac, regional)
The subroutine 'tracer_2d_1L' performs 2-D horizontal-to-lagrangian transport.
subroutine, public nested_grid_bc_apply_intt(var_nest, istag, jstag, npx, npy, npz, bd, step, split, BC, bctype)
The subroutine 'nested_grid_BC_apply_intT' performs linear interpolation or extrapolation in time for...
subroutine, public neg_adj2(is, ie, js, je, ng, kbot, hydrostatic, peln, delz, pt, dp, qv, ql, qr, qi, qs, qa, check_negative)
subroutine, public cubed_to_latlon(u, v, ua, va, gridstruct, npx, npy, km, mode, grid_type, domain, nested, c2l_ord, bd)
real, dimension(:), allocatable rf
subroutine, public tracer_2d_nested(q, dp1, mfx, mfy, cx, cy, gridstruct, bd, domain, npx, npy, npz, nq, hord, q_split, dt, id_divg, q_pack, nord_tr, trdm, k_split, neststruct, parent_grid, lim_fac, regional)
integer, parameter, public u_stagger
The module 'fv_nesting' is a collection of routines pertaining to grid nesting .