105 use mpp_domains_mod
, only: mpp_update_domains
106 use mpp_domains_mod
, only: mpp_global_field
107 use field_manager_mod
, only: model_atmos
108 use tracer_manager_mod
, only: get_tracer_index
110 use mpp_domains_mod
, only: mpp_get_data_domain, mpp_get_compute_domain, mpp_get_global_domain
111 use mpp_domains_mod
, only: dgrid_ne, mpp_update_domains, domain2d
113 use mpp_mod
, only: mpp_sync_self, mpp_sync, mpp_send, mpp_recv, mpp_error, fatal
114 use mpp_domains_mod
, only: mpp_global_sum, bitwise_efp_sum, bitwise_exact_sum
117 use fv_mp_mod, only: is, ie, js, je, isd, ied, jsd, jed, isc, iec, jsc, jec
122 use constants_mod
, only: grav, pi=>pi_8, radius, hlv, rdgas, cp_air, rvgas, cp_vapor, kappa
133 real,
allocatable ::
rf(:),
rw(:)
158 u, v, w, pt, delp, delz,q, uc, vc, pkz, &
159 nested, inline_q, make_nh, ng, &
160 gridstruct, flagstruct, neststruct, &
161 nest_timestep, tracer_nest_timestep, &
166 real,
intent(IN) :: zvir
168 integer,
intent(IN) :: npx, npy, npz
169 integer,
intent(IN) :: ncnst, ng, nwat
170 logical,
intent(IN) :: inline_q, make_nh,nested
172 real,
intent(inout),
dimension(bd%isd:bd%ied ,bd%jsd:bd%jed+1,npz) :: u
173 real,
intent(inout),
dimension(bd%isd:bd%ied+1,bd%jsd:bd%jed ,npz) :: v
174 real,
intent(inout) :: w( bd%isd: ,bd%jsd: ,1:)
175 real,
intent(inout) :: pt( bd%isd:bd%ied ,bd%jsd:bd%jed ,npz)
176 real,
intent(inout) :: delp(bd%isd:bd%ied ,bd%jsd:bd%jed ,npz)
177 real,
intent(inout) :: delz(bd%isd: ,bd%jsd: ,1:)
178 real,
intent(inout) :: q( bd%isd:bd%ied ,bd%jsd:bd%jed ,npz, ncnst)
179 real,
intent(inout) :: uc(bd%isd:bd%ied+1,bd%jsd:bd%jed ,npz)
180 real,
intent(inout) :: vc(bd%isd:bd%ied ,bd%jsd:bd%jed+1,npz)
181 real,
intent(inout) :: pkz (bd%is:bd%ie,bd%js:bd%je,npz)
182 integer,
intent(INOUT) :: nest_timestep, tracer_nest_timestep
187 type(domain2d),
intent(INOUT) :: domain
188 real :: divg(bd%isd:bd%ied+1,bd%jsd:bd%jed+1, npz)
189 real :: ua(bd%isd:bd%ied,bd%jsd:bd%jed)
190 real :: va(bd%isd:bd%ied,bd%jsd:bd%jed)
192 real :: pkz_coarse( bd%isd:bd%ied ,bd%jsd:bd%jed ,npz)
193 integer :: i,j,k,n,p, sphum
199 logical,
pointer :: child_grids(:)
201 integer :: is, ie, js, je
202 integer :: isd, ied, jsd, jed
213 child_grids => neststruct%child_grids
220 if (.not. inline_q) tracer_nest_timestep = 0
223 if (neststruct%nested .and. (.not. (neststruct%first_step) .or. make_nh) )
then 225 call set_bcs_t0(ncnst, flagstruct%hydrostatic, neststruct)
233 if (any(neststruct%child_grids))
then 237 call mpp_update_domains(u, v, &
238 domain, gridtype=dgrid_ne, complete=.true.)
244 call d2c_setup(u(isd,jsd,k), v(isd,jsd,k), &
246 uc(isd,jsd,k), vc(isd,jsd,k), flagstruct%nord>0, &
247 isd,ied,jsd,jed, is,ie,js,je, npx,npy, &
248 gridstruct%grid_type, gridstruct%nested, &
249 gridstruct%se_corner, gridstruct%sw_corner, &
250 gridstruct%ne_corner, gridstruct%nw_corner, &
251 gridstruct%rsin_u, gridstruct%rsin_v, &
252 gridstruct%cosa_s, gridstruct%rsin2, flagstruct%regional )
254 call divergence_corner_nest(u(isd,jsd,k), v(isd,jsd,k), ua, va, divg(isd,jsd,k), gridstruct, flagstruct, bd)
256 call divergence_corner(u(isd,jsd,k), v(isd,jsd,k), ua, va, divg(isd,jsd,k), gridstruct, flagstruct, bd)
262 if (flagstruct%hydrostatic)
then 267 pkz_coarse(i,j,k) = pkz(i,j,k)
274 if (neststruct%nested)
then 275 if (.not.
allocated(
q_buf))
then 276 allocate(
q_buf(ncnst))
289 if (flagstruct%hydrostatic)
then 290 call allocate_fv_nest_bc_type(pkz_bc,is,ie,js,je,isd,ied,jsd,jed,npx,npy,npz,ng,0,0,0,.false.)
315 do p=1,
size(child_grids)
316 if (child_grids(p))
then 324 if (flagstruct%hydrostatic)
then 343 neststruct%ind_h, neststruct%wt_h, 0, 0, npx, npy, npz, bd, &
344 neststruct%delp_BC,
delp_buf, pd_in=do_pd)
347 neststruct%ind_h, neststruct%wt_h, 0, 0, npx, npy, npz, bd, &
348 neststruct%q_BC(n),
q_buf(n), pd_in=do_pd)
352 neststruct%ind_h, neststruct%wt_h, 0, 0, npx, npy, npz, bd, &
355 sphum = get_tracer_index(model_atmos,
'sphum')
356 if (flagstruct%hydrostatic)
then 358 neststruct%ind_h, neststruct%wt_h, 0, 0, npx, npy, npz, bd, &
360 call setup_pt_bc(neststruct%pt_BC, pkz_bc, neststruct%q_BC(sphum), npx, npy, npz, zvir, bd)
363 neststruct%ind_h, neststruct%wt_h, 0, 0, npx, npy, npz, bd, &
364 neststruct%w_BC,
w_buf)
366 neststruct%ind_h, neststruct%wt_h, 0, 0, npx, npy, npz, bd, &
369 call setup_pt_nh_bc(neststruct%pt_BC, neststruct%delp_BC, neststruct%delz_BC, &
370 neststruct%q_BC(sphum), neststruct%q_BC, ncnst, &
372 neststruct%q_con_BC, &
374 neststruct%cappa_BC, &
377 npx, npy, npz, zvir, bd)
381 neststruct%ind_u, neststruct%wt_u, 0, 1, npx, npy, npz, bd, &
382 neststruct%u_BC,
u_buf)
384 neststruct%ind_u, neststruct%wt_u, 0, 1, npx, npy, npz, bd, &
387 neststruct%ind_v, neststruct%wt_v, 1, 0, npx, npy, npz, bd, &
388 neststruct%v_BC,
v_buf)
390 neststruct%ind_v, neststruct%wt_v, 1, 0, npx, npy, npz, bd, &
394 neststruct%ind_b, neststruct%wt_b, 1, 1, npx, npy, npz, bd, &
398 if (neststruct%first_step)
then 399 if (neststruct%nested)
call set_bcs_t0(ncnst, flagstruct%hydrostatic, neststruct)
400 neststruct%first_step = .false.
401 if (.not. flagstruct%hydrostatic) flagstruct%make_nh= .false.
402 else if (flagstruct%make_nh)
then 404 flagstruct%make_nh= .false.
421 subroutine setup_pt_bc(pt_BC, pkz_BC, sphum_BC, npx, npy, npz, zvir, bd)
426 integer,
intent(IN) :: npx, npy, npz
427 real,
intent(IN) :: zvir
429 real,
dimension(:,:,:),
pointer :: ptBC, pkzBC, sphumBC
431 integer :: i,j,k, istart, iend
433 integer :: is, ie, js, je
434 integer :: isd, ied, jsd, jed
446 ptbc => pt_bc%west_t1
447 pkzbc => pkz_bc%west_t1
448 sphumbc => sphum_bc%west_t1
453 ptbc(i,j,k) = ptbc(i,j,k)/pkzbc(i,j,k)*(1.+zvir*sphumbc(i,j,k))
460 ptbc => pt_bc%south_t1
461 pkzbc => pkz_bc%south_t1
462 sphumbc => sphum_bc%south_t1
468 if (ie == npx-1)
then 478 ptbc(i,j,k) = ptbc(i,j,k)/pkzbc(i,j,k) * &
479 (1.+zvir*sphumbc(i,j,k))
486 if (ie == npx-1)
then 487 ptbc => pt_bc%east_t1
488 pkzbc => pkz_bc%east_t1
489 sphumbc => sphum_bc%east_t1
494 ptbc(i,j,k) = ptbc(i,j,k)/pkzbc(i,j,k) * &
495 (1.+zvir*sphumbc(i,j,k))
501 if (je == npy-1)
then 502 ptbc => pt_bc%north_t1
503 pkzbc => pkz_bc%north_t1
504 sphumbc => sphum_bc%north_t1
510 if (ie == npx-1)
then 520 ptbc(i,j,k) = ptbc(i,j,k)/pkzbc(i,j,k) * &
521 (1.+zvir*sphumbc(i,j,k))
529 subroutine setup_pt_nh_bc(pt_BC, delp_BC, delz_BC, sphum_BC, q_BC, nq, &
536 npx, npy, npz, zvir, bd)
541 integer,
intent(IN) :: nq
549 integer,
intent(IN) :: npx, npy, npz
550 real,
intent(IN) :: zvir
552 real,
parameter:: c_liq = 4185.5
553 real,
parameter:: c_ice = 1972.
554 real,
parameter:: cv_vap = cp_vapor - rvgas
556 real,
dimension(:,:,:),
pointer :: ptBC, sphumBC, qconBC, delpBC, delzBC, cappaBC
557 real,
dimension(:,:,:),
pointer :: liq_watBC_west, ice_watBC_west, rainwatBC_west, snowwatBC_west, graupelBC_west
558 real,
dimension(:,:,:),
pointer :: liq_watBC_east, ice_watBC_east, rainwatBC_east, snowwatBC_east, graupelBC_east
559 real,
dimension(:,:,:),
pointer :: liq_watBC_north, ice_watBC_north, rainwatBC_north, snowwatBC_north, graupelBC_north
560 real,
dimension(:,:,:),
pointer :: liq_watBC_south, ice_watBC_south, rainwatBC_south, snowwatBC_south, graupelBC_south
562 real :: dp1, q_liq, q_sol, q_con = 0., cvm, pkz, rdg, cv_air
564 integer :: i,j,k, istart, iend
565 integer :: liq_wat, ice_wat, rainwat, snowwat, graupel
566 real,
parameter:: tice = 273.16
567 real,
parameter:: t_i0 = 15.
569 integer :: is, ie, js, je
570 integer :: isd, ied, jsd, jed
582 cv_air = cp_air - rdgas
584 liq_wat = get_tracer_index(model_atmos,
'liq_wat')
585 ice_wat = get_tracer_index(model_atmos,
'ice_wat')
586 rainwat = get_tracer_index(model_atmos,
'rainwat')
587 snowwat = get_tracer_index(model_atmos,
'snowwat')
588 graupel = get_tracer_index(model_atmos,
'graupel')
592 allocate(
dum_west(isd:0,jsd:jed,npz))
616 if (ie == npx-1)
then 618 allocate(
dum_east(npx:ied,jsd:jed,npz))
629 if (je == npy-1)
then 643 if (liq_wat > 0)
then 644 liq_watbc_west => q_bc(liq_wat)%west_t1
645 liq_watbc_east => q_bc(liq_wat)%east_t1
646 liq_watbc_north => q_bc(liq_wat)%north_t1
647 liq_watbc_south => q_bc(liq_wat)%south_t1
654 if (ice_wat > 0)
then 655 ice_watbc_west => q_bc(ice_wat)%west_t1
656 ice_watbc_east => q_bc(ice_wat)%east_t1
657 ice_watbc_north => q_bc(ice_wat)%north_t1
658 ice_watbc_south => q_bc(ice_wat)%south_t1
665 if (rainwat > 0)
then 666 rainwatbc_west => q_bc(rainwat)%west_t1
667 rainwatbc_east => q_bc(rainwat)%east_t1
668 rainwatbc_north => q_bc(rainwat)%north_t1
669 rainwatbc_south => q_bc(rainwat)%south_t1
676 if (snowwat > 0)
then 677 snowwatbc_west => q_bc(snowwat)%west_t1
678 snowwatbc_east => q_bc(snowwat)%east_t1
679 snowwatbc_north => q_bc(snowwat)%north_t1
680 snowwatbc_south => q_bc(snowwat)%south_t1
687 if (graupel > 0)
then 688 graupelbc_west => q_bc(graupel)%west_t1
689 graupelbc_east => q_bc(graupel)%east_t1
690 graupelbc_north => q_bc(graupel)%north_t1
691 graupelbc_south => q_bc(graupel)%south_t1
700 ptbc => pt_bc%west_t1
701 sphumbc => sphum_bc%west_t1
703 qconbc => q_con_bc%west_t1
705 cappabc => cappa_bc%west_t1
708 delpbc => delp_bc%west_t1
709 delzbc => delz_bc%west_t1
717 dp1 = zvir*sphumbc(i,j,k)
720 q_con = liq_watbc_west(i,j,k)
721 q_sol = q_con*max(min((tice-ptbc(i,j,k))/t_i0,1.),0.)
722 q_liq = q_con - q_sol
724 q_liq = liq_watbc_west(i,j,k) + rainwatbc_west(i,j,k)
725 q_sol = ice_watbc_west(i,j,k) + snowwatbc_west(i,j,k) + graupelbc_west(i,j,k)
726 q_con = q_liq + q_sol
728 qconbc(i,j,k) = q_con
730 cvm = (1.-(sphumbc(i,j,k)+q_con))*cv_air+sphumbc(i,j,k)*cv_vap+q_liq*c_liq+q_sol*c_ice
731 cappabc(i,j,k) = rdgas/(rdgas + cvm/(1.+dp1))
732 pkz = exp( cappabc(i,j,k)*log(rdg*delpbc(i,j,k)*ptbc(i,j,k) * &
733 (1.+dp1)*(1.-q_con)/delzbc(i,j,k)))
735 pkz = exp( kappa*log(rdg*delpbc(i,j,k)*ptbc(i,j,k) * &
736 (1.+dp1)*(1.-q_con)/delzbc(i,j,k)))
738 ptbc(i,j,k) = ptbc(i,j,k)*(1.+dp1)*(1.-q_con)/pkz
740 pkz = exp( kappa*log(rdg*delpbc(i,j,k)*ptbc(i,j,k) * &
741 (1.+dp1)/delzbc(i,j,k)))
742 ptbc(i,j,k) = ptbc(i,j,k)*(1.+dp1)/pkz
751 ptbc => pt_bc%south_t1
752 sphumbc => sphum_bc%south_t1
754 qconbc => q_con_bc%south_t1
756 cappabc => cappa_bc%south_t1
759 delpbc => delp_bc%south_t1
760 delzbc => delz_bc%south_t1
766 if (ie == npx-1)
then 780 dp1 = zvir*sphumbc(i,j,k)
783 q_con = liq_watbc_south(i,j,k)
784 q_sol = q_con*max(min((tice-ptbc(i,j,k))/t_i0,1.),0.)
785 q_liq = q_con - q_sol
787 q_liq = liq_watbc_south(i,j,k) + rainwatbc_south(i,j,k)
788 q_sol = ice_watbc_south(i,j,k) + snowwatbc_south(i,j,k) + graupelbc_south(i,j,k)
789 q_con = q_liq + q_sol
791 qconbc(i,j,k) = q_con
793 cvm = (1.-(sphumbc(i,j,k)+q_con))*cv_air+sphumbc(i,j,k)*cv_vap+q_liq*c_liq+q_sol*c_ice
794 cappabc(i,j,k) = rdgas/(rdgas + cvm/(1.+dp1))
795 pkz = exp( cappabc(i,j,k)*log(rdg*delpbc(i,j,k)*ptbc(i,j,k) * &
796 (1.+dp1)*(1.-q_con)/delzbc(i,j,k)))
798 pkz = exp( kappa*log(rdg*delpbc(i,j,k)*ptbc(i,j,k) * &
799 (1.+dp1)*(1.-q_con)/delzbc(i,j,k)))
801 ptbc(i,j,k) = ptbc(i,j,k)*(1.+dp1)*(1.-q_con)/pkz
803 pkz = exp( kappa*log(rdg*delpbc(i,j,k)*ptbc(i,j,k) * &
804 (1.+dp1)/delzbc(i,j,k)))
805 ptbc(i,j,k) = ptbc(i,j,k)*(1.+dp1)/pkz
813 if (ie == npx-1)
then 814 ptbc => pt_bc%east_t1
815 sphumbc => sphum_bc%east_t1
817 qconbc => q_con_bc%east_t1
819 cappabc => cappa_bc%east_t1
822 delpbc => delp_bc%east_t1
823 delzbc => delz_bc%east_t1
831 dp1 = zvir*sphumbc(i,j,k)
834 q_con = liq_watbc_east(i,j,k)
835 q_sol = q_con*max(min((tice-ptbc(i,j,k))/t_i0,1.),0.)
836 q_liq = q_con - q_sol
838 q_liq = liq_watbc_east(i,j,k) + rainwatbc_east(i,j,k)
839 q_sol = ice_watbc_east(i,j,k) + snowwatbc_east(i,j,k) + graupelbc_east(i,j,k)
840 q_con = q_liq + q_sol
842 qconbc(i,j,k) = q_con
844 cvm = (1.-(sphumbc(i,j,k)+q_con))*cv_air+sphumbc(i,j,k)*cv_vap+q_liq*c_liq+q_sol*c_ice
845 cappabc(i,j,k) = rdgas/(rdgas + cvm/(1.+dp1))
846 pkz = exp( cappabc(i,j,k)*log(rdg*delpbc(i,j,k)*ptbc(i,j,k) * &
847 (1.+dp1)*(1.-q_con)/delzbc(i,j,k)))
849 pkz = exp( kappa*log(rdg*delpbc(i,j,k)*ptbc(i,j,k) * &
850 (1.+dp1)*(1.-q_con)/delzbc(i,j,k)))
852 ptbc(i,j,k) = ptbc(i,j,k)*(1.+dp1)*(1.-q_con)/pkz
854 pkz = exp( kappa*log(rdg*delpbc(i,j,k)*ptbc(i,j,k) * &
855 (1.+dp1)/delzbc(i,j,k)))
856 ptbc(i,j,k) = ptbc(i,j,k)*(1.+dp1)/pkz
863 if (je == npy-1)
then 864 ptbc => pt_bc%north_t1
865 sphumbc => sphum_bc%north_t1
867 qconbc => q_con_bc%north_t1
869 cappabc => cappa_bc%north_t1
872 delpbc => delp_bc%north_t1
873 delzbc => delz_bc%north_t1
879 if (ie == npx-1)
then 892 dp1 = zvir*sphumbc(i,j,k)
895 q_con = liq_watbc_north(i,j,k)
896 q_sol = q_con*max(min((tice-ptbc(i,j,k))/t_i0,1.),0.)
897 q_liq = q_con - q_sol
899 q_liq = liq_watbc_north(i,j,k) + rainwatbc_north(i,j,k)
900 q_sol = ice_watbc_north(i,j,k) + snowwatbc_north(i,j,k) + graupelbc_north(i,j,k)
901 q_con = q_liq + q_sol
903 qconbc(i,j,k) = q_con
905 cvm = (1.-(sphumbc(i,j,k)+q_con))*cv_air+sphumbc(i,j,k)*cv_vap+q_liq*c_liq+q_sol*c_ice
906 cappabc(i,j,k) = rdgas/(rdgas + cvm/(1.+dp1))
907 pkz = exp( cappabc(i,j,k)*log(rdg*delpbc(i,j,k)*ptbc(i,j,k) * &
908 (1.+dp1)*(1.-q_con)/delzbc(i,j,k)))
910 pkz = exp( kappa*log(rdg*delpbc(i,j,k)*ptbc(i,j,k) * &
911 (1.+dp1)*(1.-q_con)/delzbc(i,j,k)))
913 ptbc(i,j,k) = ptbc(i,j,k)*(1.+dp1)*(1.-q_con)/pkz
915 pkz = exp( kappa*log(rdg*delpbc(i,j,k)*ptbc(i,j,k) * &
916 (1.+dp1)/delzbc(i,j,k)))
917 ptbc(i,j,k) = ptbc(i,j,k)*(1.+dp1)/pkz
934 neststruct%delz_BC%east_t0 = neststruct%delz_BC%east_t1
935 neststruct%delz_BC%west_t0 = neststruct%delz_BC%west_t1
936 neststruct%delz_BC%north_t0 = neststruct%delz_BC%north_t1
937 neststruct%delz_BC%south_t0 = neststruct%delz_BC%south_t1
939 neststruct%w_BC%east_t0 = neststruct%w_BC%east_t1
940 neststruct%w_BC%west_t0 = neststruct%w_BC%west_t1
941 neststruct%w_BC%north_t0 = neststruct%w_BC%north_t1
942 neststruct%w_BC%south_t0 = neststruct%w_BC%south_t1
947 subroutine set_bcs_t0(ncnst, hydrostatic, neststruct)
949 integer,
intent(IN) :: ncnst
950 logical,
intent(IN) :: hydrostatic
955 neststruct%delp_BC%east_t0 = neststruct%delp_BC%east_t1
956 neststruct%delp_BC%west_t0 = neststruct%delp_BC%west_t1
957 neststruct%delp_BC%north_t0 = neststruct%delp_BC%north_t1
958 neststruct%delp_BC%south_t0 = neststruct%delp_BC%south_t1
960 neststruct%q_BC(n)%east_t0 = neststruct%q_BC(n)%east_t1
961 neststruct%q_BC(n)%west_t0 = neststruct%q_BC(n)%west_t1
962 neststruct%q_BC(n)%north_t0 = neststruct%q_BC(n)%north_t1
963 neststruct%q_BC(n)%south_t0 = neststruct%q_BC(n)%south_t1
966 neststruct%pt_BC%east_t0 = neststruct%pt_BC%east_t1
967 neststruct%pt_BC%west_t0 = neststruct%pt_BC%west_t1
968 neststruct%pt_BC%north_t0 = neststruct%pt_BC%north_t1
969 neststruct%pt_BC%south_t0 = neststruct%pt_BC%south_t1
970 neststruct%pt_BC%east_t0 = neststruct%pt_BC%east_t1
971 neststruct%pt_BC%west_t0 = neststruct%pt_BC%west_t1
972 neststruct%pt_BC%north_t0 = neststruct%pt_BC%north_t1
973 neststruct%pt_BC%south_t0 = neststruct%pt_BC%south_t1
976 neststruct%q_con_BC%east_t0 = neststruct%q_con_BC%east_t1
977 neststruct%q_con_BC%west_t0 = neststruct%q_con_BC%west_t1
978 neststruct%q_con_BC%north_t0 = neststruct%q_con_BC%north_t1
979 neststruct%q_con_BC%south_t0 = neststruct%q_con_BC%south_t1
981 neststruct%cappa_BC%east_t0 = neststruct%cappa_BC%east_t1
982 neststruct%cappa_BC%west_t0 = neststruct%cappa_BC%west_t1
983 neststruct%cappa_BC%north_t0 = neststruct%cappa_BC%north_t1
984 neststruct%cappa_BC%south_t0 = neststruct%cappa_BC%south_t1
988 if (.not. hydrostatic)
then 992 neststruct%u_BC%east_t0 = neststruct%u_BC%east_t1
993 neststruct%u_BC%west_t0 = neststruct%u_BC%west_t1
994 neststruct%u_BC%north_t0 = neststruct%u_BC%north_t1
995 neststruct%u_BC%south_t0 = neststruct%u_BC%south_t1
996 neststruct%v_BC%east_t0 = neststruct%v_BC%east_t1
997 neststruct%v_BC%west_t0 = neststruct%v_BC%west_t1
998 neststruct%v_BC%north_t0 = neststruct%v_BC%north_t1
999 neststruct%v_BC%south_t0 = neststruct%v_BC%south_t1
1002 neststruct%vc_BC%east_t0 = neststruct%vc_BC%east_t1
1003 neststruct%vc_BC%west_t0 = neststruct%vc_BC%west_t1
1004 neststruct%vc_BC%north_t0 = neststruct%vc_BC%north_t1
1005 neststruct%vc_BC%south_t0 = neststruct%vc_BC%south_t1
1006 neststruct%uc_BC%east_t0 = neststruct%uc_BC%east_t1
1007 neststruct%uc_BC%west_t0 = neststruct%uc_BC%west_t1
1008 neststruct%uc_BC%north_t0 = neststruct%uc_BC%north_t1
1009 neststruct%uc_BC%south_t0 = neststruct%uc_BC%south_t1
1011 neststruct%divg_BC%east_t0 = neststruct%divg_BC%east_t1
1012 neststruct%divg_BC%west_t0 = neststruct%divg_BC%west_t1
1013 neststruct%divg_BC%north_t0 = neststruct%divg_BC%north_t1
1014 neststruct%divg_BC%south_t0 = neststruct%divg_BC%south_t1
1046 integer,
intent(IN) :: ngrids
1047 logical,
intent(IN) :: grids_on_this_pe(ngrids)
1048 real,
intent(IN) :: zvir
1050 integer :: n, p, sphum
1053 if (ngrids > 1)
then 1058 if (atm(n)%neststruct%twowaynest )
then 1059 if (grids_on_this_pe(n) .or. grids_on_this_pe(atm(n)%parent_grid%grid_number))
then 1060 sphum = get_tracer_index(model_atmos,
'sphum')
1062 atm(n)%ncnst, sphum, atm(n)%u, atm(n)%v, atm(n)%w, atm(n)%omga, &
1063 atm(n)%pt, atm(n)%delp, atm(n)%q, atm(n)%uc, atm(n)%vc, &
1064 atm(n)%pkz, atm(n)%delz, atm(n)%ps, atm(n)%ptop, &
1065 atm(n)%gridstruct, atm(n)%flagstruct, atm(n)%neststruct, atm(n)%parent_grid, atm(n)%bd, .false.)
1073 if (atm(n)%neststruct%parent_of_twoway .and. grids_on_this_pe(n))
then 1075 atm(n)%u, atm(n)%v, atm(n)%w, atm(n)%delz, &
1076 atm(n)%pt, atm(n)%delp, atm(n)%q, &
1077 atm(n)%ps, atm(n)%pe, atm(n)%pk, atm(n)%peln, atm(n)%pkz, &
1078 atm(n)%phis, atm(n)%ua, atm(n)%va, &
1079 atm(n)%ptop, atm(n)%gridstruct, atm(n)%flagstruct, &
1080 atm(n)%domain, atm(n)%bd)
1091 u, v, w, omga, pt, delp, q, &
1092 uc, vc, pkz, delz, ps, ptop, &
1093 gridstruct, flagstruct, neststruct, &
1094 parent_grid, bd, conv_theta_in)
1096 real,
intent(IN) :: zvir, ptop
1098 integer,
intent(IN) :: npx, npy, npz
1099 integer,
intent(IN) :: ncnst, sphum
1100 logical,
intent(IN),
OPTIONAL :: conv_theta_in
1103 real,
intent(inout),
dimension(bd%isd:bd%ied ,bd%jsd:bd%jed+1,npz) :: u
1104 real,
intent(inout),
dimension(bd%isd:bd%ied+1,bd%jsd:bd%jed ,npz) :: v
1105 real,
intent(inout) :: w( bd%isd: ,bd%jsd: ,1: )
1106 real,
intent(inout) :: omga(bd%isd:bd%ied,bd%jsd:bd%jed,npz)
1107 real,
intent(inout) :: pt( bd%isd:bd%ied ,bd%jsd:bd%jed ,npz)
1108 real,
intent(inout) :: delp(bd%isd:bd%ied ,bd%jsd:bd%jed ,npz)
1109 real,
intent(inout) :: q( bd%isd:bd%ied ,bd%jsd:bd%jed ,npz, ncnst)
1110 real,
intent(inout) :: uc(bd%isd:bd%ied+1,bd%jsd:bd%jed ,npz)
1111 real,
intent(inout) :: vc(bd%isd:bd%ied ,bd%jsd:bd%jed+1,npz)
1113 real,
intent(inout) :: pkz (bd%is:bd%ie,bd%js:bd%je,npz)
1114 real,
intent(inout) :: delz(bd%isd: ,bd%jsd: ,1: )
1115 real,
intent(inout) :: ps (bd%isd:bd%ied ,bd%jsd:bd%jed)
1123 real,
allocatable :: t_nest(:,:,:), ps0(:,:)
1125 integer :: isd_p, ied_p, jsd_p, jed_p, isc_p, iec_p, jsc_p, jec_p
1126 integer :: isg, ieg, jsg,jeg, npx_p, npy_p
1127 integer :: istart, iend
1128 real :: qmass_b, qmass_a, fix = 1.
1129 logical :: used, conv_theta=.true.
1131 real :: qdp( bd%isd:bd%ied ,bd%jsd:bd%jed ,npz)
1132 real,
allocatable :: qdp_coarse(:,:,:)
1133 real(kind=f_p),
allocatable :: q_diff(:,:,:)
1134 real :: L_sum_b(npz), L_sum_a(npz)
1137 integer :: is, ie, js, je
1138 integer :: isd, ied, jsd, jed
1139 integer :: isu, ieu, jsu, jeu
1149 isu = neststruct%isu
1150 ieu = neststruct%ieu
1151 jsu = neststruct%jsu
1152 jeu = neststruct%jeu
1154 upoff = neststruct%upoff
1158 if (
present(conv_theta_in)) conv_theta = conv_theta_in
1160 if ((.not. neststruct%parent_proc) .and. (.not. neststruct%child_proc))
return 1162 call mpp_get_data_domain( parent_grid%domain, &
1163 isd_p, ied_p, jsd_p, jed_p )
1164 call mpp_get_compute_domain( parent_grid%domain, &
1165 isc_p, iec_p, jsc_p, jec_p )
1170 if (neststruct%nestupdate < 3)
then 1173 neststruct%ind_update_h, gridstruct%dx, gridstruct%dy, gridstruct%area, &
1174 isd_p, ied_p, jsd_p, jed_p, isd, ied, jsd, jed, &
1175 neststruct%isu, neststruct%ieu, neststruct%jsu, neststruct%jeu, &
1176 npx, npy, npz, 0, 0, &
1177 neststruct%refinement, neststruct%nestupdate, upoff, 0, &
1178 neststruct%parent_proc, neststruct%child_proc, parent_grid)
1183 if (neststruct%parent_proc)
then 1187 parent_grid%ps(i,j) = &
1188 parent_grid%delp(i,j,1)/grav
1198 if (neststruct%nestupdate /= 3 .and. neststruct%nestupdate /= 7 .and. neststruct%nestupdate /= 8)
then 1200 allocate(qdp_coarse(isd_p:ied_p,jsd_p:jed_p,npz))
1201 if (parent_grid%flagstruct%nwat > 0)
then 1202 allocate(q_diff(isd_p:ied_p,jsd_p:jed_p,npz))
1206 do n=1,parent_grid%flagstruct%nwat
1209 if (neststruct%child_proc)
then 1213 qdp(i,j,k) = q(i,j,k,n)*delp(i,j,k)
1221 if (neststruct%parent_proc)
then 1226 qdp_coarse(i,j,k) = parent_grid%q(i,j,k,n)*parent_grid%delp(i,j,k)
1230 call level_sum(qdp_coarse, parent_grid%gridstruct%area, parent_grid%domain, &
1231 parent_grid%bd, npz, l_sum_b)
1235 if (neststruct%parent_proc)
then 1236 if (n <= parent_grid%flagstruct%nwat)
then 1240 q_diff(i,j,k) = q_diff(i,j,k) - qdp_coarse(i,j,k)
1248 neststruct%ind_update_h, gridstruct%dx, gridstruct%dy, gridstruct%area, &
1249 isd_p, ied_p, jsd_p, jed_p, isd, ied, jsd, jed, &
1250 neststruct%isu, neststruct%ieu, neststruct%jsu, neststruct%jeu, &
1251 npx, npy, npz, 0, 0, &
1252 neststruct%refinement, neststruct%nestupdate, upoff, 0, neststruct%parent_proc, neststruct%child_proc, parent_grid)
1256 if (neststruct%parent_proc)
then 1257 call level_sum(qdp_coarse, parent_grid%gridstruct%area, parent_grid%domain, &
1258 parent_grid%bd, npz, l_sum_a)
1260 if (l_sum_a(k) > 0.)
then 1261 fix = l_sum_b(k)/l_sum_a(k)
1265 parent_grid%q(i,j,k,n) = qdp_coarse(i,j,k)*fix
1272 if (neststruct%parent_proc)
then 1273 if (n <= parent_grid%flagstruct%nwat)
then 1277 q_diff(i,j,k) = q_diff(i,j,k) + parent_grid%q(i,j,k,n)
1286 if (neststruct%parent_proc)
then 1287 if (parent_grid%flagstruct%nwat > 0)
then 1291 parent_grid%delp(i,j,k) = parent_grid%delp(i,j,k) + q_diff(i,j,k)
1297 do n=1,parent_grid%flagstruct%nwat
1301 parent_grid%q(i,j,k,n) = parent_grid%q(i,j,k,n)/parent_grid%delp(i,j,k)
1308 deallocate(qdp_coarse)
1309 if (
allocated(q_diff))
deallocate(q_diff)
1314 if (neststruct%nestupdate /= 3 .and. neststruct%nestupdate /= 8)
then 1316 if (conv_theta)
then 1318 if (neststruct%child_proc)
then 1322 allocate(t_nest(isd:ied,jsd:jed,1:npz))
1327 t_nest(i,j,k) = pt(i,j,k)*pkz(i,j,k)/(1.+zvir*q(i,j,k,sphum))
1335 t_nest, neststruct%nest_domain, &
1336 neststruct%ind_update_h, gridstruct%dx, gridstruct%dy, gridstruct%area, &
1337 isd_p, ied_p, jsd_p, jed_p, isd, ied, jsd, jed, &
1338 neststruct%isu, neststruct%ieu, neststruct%jsu, neststruct%jeu, &
1339 npx, npy, npz, 0, 0, &
1340 neststruct%refinement, neststruct%nestupdate, upoff, 0, neststruct%parent_proc, neststruct%child_proc, parent_grid)
1344 pt, neststruct%nest_domain, &
1345 neststruct%ind_update_h, gridstruct%dx, gridstruct%dy, gridstruct%area, &
1346 isd_p, ied_p, jsd_p, jed_p, isd, ied, jsd, jed, &
1347 neststruct%isu, neststruct%ieu, neststruct%jsu, neststruct%jeu, &
1348 npx, npy, npz, 0, 0, &
1349 neststruct%refinement, neststruct%nestupdate, upoff, 0, neststruct%parent_proc, neststruct%child_proc, parent_grid)
1356 if (.not. flagstruct%hydrostatic)
then 1359 neststruct%ind_update_h, gridstruct%dx, gridstruct%dy, gridstruct%area, &
1360 isd_p, ied_p, jsd_p, jed_p, isd, ied, jsd, jed, &
1361 neststruct%isu, neststruct%ieu, neststruct%jsu, neststruct%jeu, &
1362 npx, npy, npz, 0, 0, &
1363 neststruct%refinement, neststruct%nestupdate, upoff, 0, neststruct%parent_proc, neststruct%child_proc, parent_grid)
1379 neststruct%ind_update_h, gridstruct%dx, gridstruct%dy, gridstruct%area, &
1380 isd_p, ied_p, jsd_p, jed_p, isd, ied, jsd, jed, &
1381 neststruct%isu, neststruct%ieu, neststruct%jsu, neststruct%jeu, &
1382 npx, npy, npz, 0, 1, &
1383 neststruct%refinement, neststruct%nestupdate, upoff, 0, neststruct%parent_proc, neststruct%child_proc, parent_grid)
1386 neststruct%ind_update_h, gridstruct%dx, gridstruct%dy, gridstruct%area, &
1387 isd_p, ied_p, jsd_p, jed_p, isd, ied, jsd, jed, &
1388 neststruct%isu, neststruct%ieu, neststruct%jsu, neststruct%jeu, &
1389 npx, npy, npz, 1, 0, &
1390 neststruct%refinement, neststruct%nestupdate, upoff, 0, neststruct%parent_proc, neststruct%child_proc, parent_grid)
1396 if (neststruct%nestupdate >= 5 .and. npz > 4)
then 1403 allocate(ps0(isd_p:ied_p,jsd_p:jed_p))
1404 if (neststruct%parent_proc)
then 1406 parent_grid%ps = parent_grid%ptop
1412 parent_grid%ps(i,j) = parent_grid%ps(i,j) + &
1413 parent_grid%delp(i,j,k)
1418 ps0 = parent_grid%ps
1421 if (neststruct%child_proc)
then 1428 ps(i,j) = ps(i,j) + delp(i,j,k)
1435 neststruct%ind_update_h, gridstruct%dx, gridstruct%dy, gridstruct%area, &
1436 isd_p, ied_p, jsd_p, jed_p, isd, ied, jsd, jed, &
1437 neststruct%isu, neststruct%ieu, neststruct%jsu, neststruct%jeu, &
1439 neststruct%refinement, neststruct%nestupdate, upoff, 0, neststruct%parent_proc, neststruct%child_proc, parent_grid)
1445 if (neststruct%parent_proc)
then 1446 call mpp_update_domains(parent_grid%ps, parent_grid%domain, complete=.false.)
1447 call mpp_update_domains(ps0, parent_grid%domain, complete=.true.)
1453 if (parent_grid%tile == neststruct%parent_tile)
then 1455 if (neststruct%parent_proc)
then 1461 if (.not. parent_grid%flagstruct%remap_t)
then 1466 parent_grid%pt(i,j,k) = &
1467 parent_grid%pt(i,j,k)/parent_grid%pkz(i,j,k)*&
1468 (1.+zvir*parent_grid%q(i,j,k,sphum))
1474 parent_grid%ps, parent_grid%delp, &
1475 parent_grid%pt, parent_grid%q, parent_grid%w, &
1476 parent_grid%flagstruct%hydrostatic, &
1477 npz, ps0, zvir, parent_grid%ptop, ncnst, &
1478 parent_grid%flagstruct%kord_tm, parent_grid%flagstruct%kord_tr, &
1479 parent_grid%flagstruct%kord_wz, &
1480 isc_p, iec_p, jsc_p, jec_p, isd_p, ied_p, jsd_p, jed_p, .false. )
1481 if (.not. parent_grid%flagstruct%remap_t)
then 1486 parent_grid%pt(i,j,k) = &
1487 parent_grid%pt(i,j,k)*parent_grid%pkz(i,j,k) / &
1488 (1.+zvir*parent_grid%q(i,j,k,sphum))
1497 parent_grid%v, npz, ps0, parent_grid%flagstruct%kord_mt, &
1498 isc_p, iec_p, jsc_p, jec_p, isd_p, ied_p, jsd_p, jed_p, parent_grid%ptop)
1504 if (
allocated(ps0))
deallocate(ps0)
1512 subroutine level_sum(q, area, domain, bd, npz, L_sum)
1514 integer,
intent(IN) :: npz
1516 real,
intent(in) :: area( bd%isd:bd%ied ,bd%jsd:bd%jed)
1517 real,
intent(in) :: q( bd%isd:bd%ied ,bd%jsd:bd%jed ,npz)
1518 real,
intent(OUT) :: L_sum( npz )
1519 type(domain2d),
intent(IN) :: domain
1521 integer :: i, j, k, n
1529 qa = qa + q(i,j,k)*area(i,j)
1532 call mp_reduce_sum(qa)
1542 u, v, w, delz, pt, delp, q, &
1543 ps, pe, pk, peln, pkz, phis, ua, va, &
1544 ptop, gridstruct, flagstruct, &
1548 real,
intent(IN) :: ptop
1550 integer,
intent(IN) :: ng, npx, npy, npz
1551 integer,
intent(IN) :: ncnst
1553 real,
intent(inout),
dimension(bd%isd:bd%ied ,bd%jsd:bd%jed+1,npz) :: u
1554 real,
intent(inout),
dimension(bd%isd:bd%ied+1,bd%jsd:bd%jed ,npz) :: v
1555 real,
intent(inout) :: w( bd%isd: ,bd%jsd: ,1: )
1556 real,
intent(inout) :: pt( bd%isd:bd%ied ,bd%jsd:bd%jed ,npz)
1557 real,
intent(inout) :: delp(bd%isd:bd%ied ,bd%jsd:bd%jed ,npz)
1558 real,
intent(inout) :: q( bd%isd:bd%ied ,bd%jsd:bd%jed ,npz, ncnst)
1559 real,
intent(inout) :: delz(bd%isd: ,bd%jsd: ,1: )
1566 real,
intent(inout) :: ps (bd%isd:bd%ied ,bd%jsd:bd%jed)
1567 real,
intent(inout) :: pe (bd%is-1:bd%ie+1, npz+1,bd%js-1:bd%je+1)
1568 real,
intent(inout) :: pk (bd%is:bd%ie,bd%js:bd%je, npz+1)
1569 real,
intent(inout) :: peln(bd%is:bd%ie,npz+1,bd%js:bd%je)
1570 real,
intent(inout) :: pkz (bd%is:bd%ie,bd%js:bd%je,npz)
1575 real,
intent(inout) :: phis(bd%isd:bd%ied,bd%jsd:bd%jed)
1577 real,
intent(inout),
dimension(bd%isd:bd%ied ,bd%jsd:bd%jed ,npz):: ua, va
1580 type(domain2d),
intent(INOUT) :: domain
1582 logical :: bad_range
1584 integer :: is, ie, js, je
1585 integer :: isd, ied, jsd, jed
1597 gridstruct, npx, npy, npz, &
1598 1, gridstruct%grid_type, domain, &
1599 gridstruct%nested, flagstruct%c2l_ord, bd)
1613 q, ng, flagstruct%ncnst, gridstruct%area_64, 0., &
1615 flagstruct%moist_phys, flagstruct%hydrostatic, &
1616 flagstruct%nwat, domain, .false.)
1620 if (flagstruct%range_warn)
then 1621 call range_check(
'TA update', pt, is, ie, js, je, ng, npz, gridstruct%agrid, 130., 350., bad_range)
1622 call range_check(
'UA update', ua, is, ie, js, je, ng, npz, gridstruct%agrid, -220., 250., bad_range)
1623 call range_check(
'VA update', va, is, ie, js, je, ng, npz, gridstruct%agrid, -220., 220., bad_range)
1624 if (.not. flagstruct%hydrostatic)
then 1625 call range_check(
'W update', w, is, ie, js, je, ng, npz, gridstruct%agrid, -50., 100., bad_range)
1636 subroutine update_remap_tqw( npz, ak, bk, ps, delp, t, q, w, hydrostatic, &
1637 kmd, ps0, zvir, ptop, nq, kord_tm, kord_tr, kord_wz, &
1638 is, ie, js, je, isd, ied, jsd, jed, do_q)
1639 integer,
intent(in):: npz, kmd, nq, kord_tm, kord_tr, kord_wz
1640 real,
intent(in):: zvir, ptop
1641 real,
intent(in):: ak(npz+1), bk(npz+1)
1642 real,
intent(in),
dimension(isd:ied,jsd:jed):: ps0
1643 real,
intent(in),
dimension(isd:ied,jsd:jed):: ps
1644 real,
intent(in),
dimension(isd:ied,jsd:jed,npz):: delp
1645 real,
intent(inout),
dimension(isd:ied,jsd:jed,npz):: t, w
1646 real,
intent(inout),
dimension(isd:ied,jsd:jed,npz,nq):: q
1647 integer,
intent(in) :: is, ie, js, je, isd, ied, jsd, jed
1648 logical,
intent(in) :: hydrostatic, do_q
1650 real,
dimension(is:ie,kmd):: tp, qp
1651 real,
dimension(is:ie,kmd+1):: pe0, pn0
1652 real,
dimension(is:ie,npz):: qn1
1653 real,
dimension(is:ie,npz+1):: pe1, pn1
1663 pe0(i,k) = ak(k) + bk(k)*ps0(i,j)
1664 pn0(i,k) = log(pe0(i,k))
1669 pe1(i,k) = ak(k) + bk(k)*ps(i,j)
1670 pn1(i,k) = log(pe1(i,k))
1677 qp(i,k) = q(i,j,k,iq)
1680 call mappm(kmd, pe0, qp, npz, pe1, qn1, is,ie, 0, kord_tr, ptop)
1683 q(i,j,k,iq) = qn1(i,k)
1695 call mappm(kmd, pn0, tp, npz, pn1, qn1, is,ie, 1, abs(kord_tm), ptop)
1703 if (.not. hydrostatic)
then 1711 call mappm(kmd, pe0, tp, npz, pe1, qn1, is,ie, -1, kord_wz, ptop)
1725 subroutine update_remap_uv(npz, ak, bk, ps, u, v, kmd, ps0, kord_mt, &
1726 is, ie, js, je, isd, ied, jsd, jed, ptop)
1727 integer,
intent(in):: npz
1728 real,
intent(in):: ak(npz+1), bk(npz+1)
1729 real,
intent(in):: ps(isd:ied,jsd:jed)
1730 real,
intent(inout),
dimension(isd:ied,jsd:jed+1,npz):: u
1731 real,
intent(inout),
dimension(isd:ied+1,jsd:jed,npz):: v
1733 integer,
intent(in):: kmd, kord_mt
1734 real,
intent(IN) :: ptop
1735 real,
intent(in):: ps0(isd:ied,jsd:jed)
1736 integer,
intent(in) :: is, ie, js, je, isd, ied, jsd, jed
1739 real,
dimension(is:ie+1,kmd+1):: pe0
1740 real,
dimension(is:ie+1,npz+1):: pe1
1741 real,
dimension(is:ie+1,kmd):: qt
1742 real,
dimension(is:ie+1,npz):: qn1
1756 pe0(i,k) = ak(k) + bk(k)*0.5*(ps0(i,j)+ps0(i,j-1))
1764 pe1(i,k) = ak(k) + bk(k)*0.5*(ps(i,j)+ps(i,j-1))
1777 call mappm(kmd, pe0(is:ie,:), qt(is:ie,:), npz, pe1(is:ie,:), qn1(is:ie,:), is,ie, -1, kord_mt, ptop)
1797 pe0(i,k) = ak(k) + bk(k)*0.5*(ps0(i,j)+ps0(i-1,j))
1805 pe1(i,k) = ak(k) + bk(k)*0.5*(ps(i,j)+ps(i-1,j))
1818 call mappm(kmd, pe0(is:ie+1,:), qt(is:ie+1,:), npz, pe1(is:ie+1,:), qn1(is:ie+1,:), is,ie+1, -1, 8, ptop)
real, dimension(:,:,:), allocatable, target dum_east
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 divergence_corner_nest(u, v, ua, va, divg_d, gridstruct, flagstruct, bd)
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.
type(fv_nest_bc_type_3d) divg_buf
The type 'fv_grid_type' is made up of grid-dependent information from fv_grid_tools and fv_grid_utils...
real, dimension(:,:,:), allocatable, target dum_south
subroutine, public neg_adj3(is, ie, js, je, ng, kbot, hydrostatic, peln, delz, pt, dp, qv, ql, qr, qi, qs, qg, qa, check_negative)
type(fv_nest_bc_type_3d) delz_buf
subroutine, public nested_grid_bc_send(var_coarse, nest_domain, istag, jstag)
The subroutine 'nested_grid_BC_send' sends coarse-grid data to create boundary conditions.
The interface'update_coarse_grid_mpp'contains subroutines that fetch data from the nested grid and in...
subroutine set_nh_bcs_t0(neststruct)
type(fv_nest_bc_type_3d) vc_buf
subroutine set_bcs_t0(ncnst, hydrostatic, neststruct)
subroutine, public mappm(km, pe1, q1, kn, pe2, q2, i1, i2, iv, kord, ptop)
The subroutine 'mappm' is a general-purpose routine for remapping one set of vertical levels to anoth...
real, public sphum_ll_fix
The module 'sw_core' advances the forward step of the Lagrangian dynamics as described by ...
type(fv_nest_bc_type_3d) uc_buf
subroutine, public divergence_corner(u, v, ua, va, divg_d, gridstruct, flagstruct, bd)
The subroutine 'divergence_corner' computes the cell-mean divergence on the "dual grid"...
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'.
The module 'fv_sg' performs FV sub-grid mixing.
'allocate_fv_nest_BC_type' is an interface to subroutines that allocate the 'fv_nest_BC_type' structu...
type(fv_nest_bc_type_3d) v_buf
subroutine twoway_nest_update(npx, npy, npz, zvir, ncnst, sphum, u, v, w, omga, pt, delp, q, uc, vc, pkz, delz, ps, ptop, gridstruct, flagstruct, neststruct, parent_grid, bd, conv_theta_in)
type(fv_nest_bc_type_3d) delp_buf
The module 'fv_timing' contains FV3 timers.
The module 'boundary' contains utility routines for grid nesting and boundary conditions.
subroutine update_remap_tqw(npz, ak, bk, ps, delp, t, q, w, hydrostatic, kmd, ps0, zvir, ptop, nq, kord_tm, kord_tr, kord_wz, is, ie, js, je, isd, ied, jsd, jed, do_q)
The subroutine 'update_remap_tqw' remaps (interpolated) nested-grid data to the coarse-grid's vertica...
subroutine, public nested_grid_bc_recv(nest_domain, istag, jstag, npz, bd, nest_BC_buffers)
subroutine 'nested_grid_BC_recv' receives coarse-grid data to create boundary conditions.
real, parameter, public ptop_min
integer, parameter, public f_p
subroutine setup_pt_bc(pt_BC, pkz_BC, sphum_BC, npx, npy, npz, zvir, bd)
subroutine, public d2a_setup(u, v, ua, va, dord4, isd, ied, jsd, jed, is, ie, js, je, npx, npy, grid_type, nested, cosa_s, rsin2, regional)
The module 'fv_mapz' contains the vertical mapping routines .
real, dimension(:,:,:), allocatable, target dum_west
type(fv_nest_bc_type_3d) w_buf
subroutine after_twoway_nest_update(npx, npy, npz, ng, ncnst, u, v, w, delz, pt, delp, q, ps, pe, pk, peln, pkz, phis, ua, va, ptop, gridstruct, flagstruct, domain, bd)
The module 'fv_arrays' contains the 'fv_atmos_type' and associated datatypes.
type(fv_nest_bc_type_3d), dimension(:), allocatable q_buf
type(fv_nest_bc_type_3d) pkz_buf
subroutine update_remap_uv(npz, ak, bk, ps, u, v, kmd, ps0, kord_mt, is, ie, js, je, isd, ied, jsd, jed, ptop)
real, dimension(:), allocatable rw
subroutine, public twoway_nesting(Atm, ngrids, grids_on_this_pe, zvir)
The subroutine'twoway_nesting' performs a two-way update of nested-grid data onto the parent grid...
type(fv_nest_bc_type_3d) pt_buf
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)
@ 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...
real, dimension(:,:,:), allocatable, target dum_north
subroutine level_sum(q, area, domain, bd, npz, L_sum)
subroutine, public cubed_to_latlon(u, v, ua, va, gridstruct, npx, npy, km, mode, grid_type, domain, nested, c2l_ord, bd)
subroutine setup_pt_nh_bc(pt_BC, delp_BC, delz_BC, sphum_BC, q_BC, nq, ifdef USE_COND
subroutine, public nested_grid_bc_save_proc(nest_domain, ind, wt, istag, jstag, npx, npy, npz, bd, nest_BC, nest_BC_buffers, pd_in)
The subroutine 'nested_grid_BC_save_proc' saves data received by 'nested_grid_BC_recv' into the datat...
real, dimension(:,:), allocatable te_2d_coarse
real, dimension(:), allocatable rf
subroutine, public d2c_setup(u, v, ua, va, uc, vc, dord4, isd, ied, jsd, jed, is, ie, js, je, npx, npy, grid_type, nested, se_corner, sw_corner, ne_corner, nw_corner, rsin_u, rsin_v, cosa_s, rsin2, regional)
subroutine, public p_var(km, ifirst, ilast, jfirst, jlast, ptop, ptop_min, delp, delz, pt, ps, pe, peln, pk, pkz, cappa, q, ng, nq, area, dry_mass, adjust_dry_mass, mountain, moist_phys, hydrostatic, nwat, domain, make_nh)
the subroutine 'p_var' computes auxiliary pressure variables for a hydrostatic state.
real, dimension(:,:,:), allocatable dp1_coarse
The module 'fv_nesting' is a collection of routines pertaining to grid nesting .
type(fv_nest_bc_type_3d) u_buf