68 use fv_mp_mod, only: group_halo_update_type
69 use fv_mp_mod, only: start_group_halo_update, complete_group_halo_update
70 use mpp_domains_mod
, only: mpp_update_domains, cgrid_ne, domain2d
76 use mpp_mod
, only: mpp_error, fatal, mpp_broadcast, mpp_send, mpp_recv, mpp_sum, mpp_max
92 subroutine tracer_2d_1l(q, dp1, mfx, mfy, cx, cy, gridstruct, bd, domain, npx, npy, npz, &
93 nq, hord, q_split, dt, id_divg, q_pack, nord_tr, trdm, lim_fac, regional)
96 integer,
intent(IN) :: npx
97 integer,
intent(IN) :: npy
98 integer,
intent(IN) :: npz
99 integer,
intent(IN) :: nq
100 integer,
intent(IN) :: hord, nord_tr
101 integer,
intent(IN) :: q_split
102 integer,
intent(IN) :: id_divg
103 real ,
intent(IN) :: dt, trdm
104 real ,
intent(IN) :: lim_fac
105 logical,
intent(IN) :: regional
106 type(group_halo_update_type),
intent(inout) :: q_pack
107 real ,
intent(INOUT) :: q(bd%isd:bd%ied,bd%jsd:bd%jed,npz,nq)
108 real ,
intent(INOUT) :: dp1(bd%isd:bd%ied,bd%jsd:bd%jed,npz)
109 real ,
intent(INOUT) :: mfx(bd%is:bd%ie+1,bd%js:bd%je, npz)
110 real ,
intent(INOUT) :: mfy(bd%is:bd%ie ,bd%js:bd%je+1,npz)
111 real ,
intent(INOUT) :: cx(bd%is:bd%ie+1,bd%jsd:bd%jed ,npz)
112 real ,
intent(INOUT) :: cy(bd%isd:bd%ied,bd%js :bd%je +1,npz)
114 type(domain2d),
intent(INOUT) :: domain
117 real :: qn2(bd%isd:bd%ied,bd%jsd:bd%jed,nq)
118 real :: dp2(bd%is:bd%ie,bd%js:bd%je)
119 real :: fx(bd%is:bd%ie+1,bd%js:bd%je )
120 real :: fy(bd%is:bd%ie , bd%js:bd%je+1)
121 real :: ra_x(bd%is:bd%ie,bd%jsd:bd%jed)
122 real :: ra_y(bd%isd:bd%ied,bd%js:bd%je)
123 real :: xfx(bd%is:bd%ie+1,bd%jsd:bd%jed ,npz)
124 real :: yfx(bd%isd:bd%ied,bd%js: bd%je+1, npz)
128 integer :: i,j,k,it,iq
130 real,
pointer,
dimension(:,:) :: area, rarea
131 real,
pointer,
dimension(:,:,:) :: sin_sg
132 real,
pointer,
dimension(:,:) :: dxa, dya, dx, dy
134 integer :: is, ie, js, je
135 integer :: isd, ied, jsd, jed
146 area => gridstruct%area
147 rarea => gridstruct%rarea
149 sin_sg => gridstruct%sin_sg
150 dxa => gridstruct%dxa
151 dya => gridstruct%dya
160 if (cx(i,j,k) > 0.)
then 161 xfx(i,j,k) = cx(i,j,k)*dxa(i-1,j)*dy(i,j)*sin_sg(i-1,j,3)
163 xfx(i,j,k) = cx(i,j,k)*dxa(i, j)*dy(i,j)*sin_sg(i, j,1)
169 if (cy(i,j,k) > 0.)
then 170 yfx(i,j,k) = cy(i,j,k)*dya(i,j-1)*dx(i,j)*sin_sg(i,j-1,4)
172 yfx(i,j,k) = cy(i,j,k)*dya(i,j )*dx(i,j)*sin_sg(i,j, 2)
178 if ( k < npz/6 )
then 181 cmax(k) = max( cmax(k), abs(cx(i,j,k)), abs(cy(i,j,k)) )
187 cmax(k) = max( cmax(k), max(abs(cx(i,j,k)),abs(cy(i,j,k)))+1.-sin_sg(i,j,5) )
193 call mp_reduce_max(cmax,npz)
200 nsplt = int(1. + cmax(k))
201 if ( nsplt > 1 )
then 202 frac = 1. /
real(nsplt)
205 cx(i,j,k) = cx(i,j,k) * frac
206 xfx(i,j,k) = xfx(i,j,k) * frac
211 mfx(i,j,k) = mfx(i,j,k) * frac
216 cy(i,j,k) = cy(i,j,k) * frac
217 yfx(i,j,k) = yfx(i,j,k) * frac
222 mfy(i,j,k) = mfy(i,j,k) * frac
230 call complete_group_halo_update(q_pack, domain)
240 ra_x(i,j) = area(i,j) + xfx(i,j,k) - xfx(i+1,j,k)
242 if ( j>=js .and. j<=je )
then 244 ra_y(i,j) = area(i,j) + yfx(i,j,k) - yfx(i,j+1,k)
249 nsplt = int(1. + cmax(k))
255 dp2(i,j) = dp1(i,j,k) + (mfx(i,j,k)-mfx(i+1,j,k)+mfy(i,j,k)-mfy(i,j+1,k))*rarea(i,j)
263 if ( nsplt /= 1 )
then 267 qn2(i,j,iq) = q(i,j,k,iq)
271 call fv_tp_2d(qn2(isd,jsd,iq), cx(is,jsd,k), cy(isd,js,k), &
272 npx, npy, hord, fx, fy, xfx(is,jsd,k), yfx(isd,js,k), &
273 gridstruct, bd, ra_x, ra_y, lim_fac, regional, mfx=mfx(is,js,k), mfy=mfy(is,js,k))
274 if ( it < nsplt )
then 277 qn2(i,j,iq) = (qn2(i,j,iq)*dp1(i,j,k)+(fx(i,j)-fx(i+1,j)+fy(i,j)-fy(i,j+1))*rarea(i,j))/dp2(i,j)
283 q(i,j,k,iq) = (qn2(i,j,iq)*dp1(i,j,k)+(fx(i,j)-fx(i+1,j)+fy(i,j)-fy(i,j+1))*rarea(i,j))/dp2(i,j)
288 call fv_tp_2d(q(isd,jsd,k,iq), cx(is,jsd,k), cy(isd,js,k), &
289 npx, npy, hord, fx, fy, xfx(is,jsd,k), yfx(isd,js,k), &
290 gridstruct, bd, ra_x, ra_y, lim_fac, regional, mfx=mfx(is,js,k), mfy=mfy(is,js,k))
293 q(i,j,k,iq) = (q(i,j,k,iq)*dp1(i,j,k)+(fx(i,j)-fx(i+1,j)+fy(i,j)-fy(i,j+1))*rarea(i,j))/dp2(i,j)
299 if ( it < nsplt )
then 302 dp1(i,j,k) = dp2(i,j)
307 call mpp_update_domains(qn2, domain)
317 subroutine tracer_2d(q, dp1, mfx, mfy, cx, cy, gridstruct, bd, domain, npx, npy, npz, &
318 nq, hord, q_split, dt, id_divg, q_pack, nord_tr, trdm, lim_fac, regional)
321 integer,
intent(IN) :: npx
322 integer,
intent(IN) :: npy
323 integer,
intent(IN) :: npz
324 integer,
intent(IN) :: nq
325 integer,
intent(IN) :: hord, nord_tr
326 integer,
intent(IN) :: q_split
327 integer,
intent(IN) :: id_divg
328 real ,
intent(IN) :: dt, trdm
329 real ,
intent(IN) :: lim_fac
330 logical,
intent(IN) :: regional
331 type(group_halo_update_type),
intent(inout) :: q_pack
332 real ,
intent(INOUT) :: q(bd%isd:bd%ied,bd%jsd:bd%jed,npz,nq)
333 real ,
intent(INOUT) :: dp1(bd%isd:bd%ied,bd%jsd:bd%jed,npz)
334 real ,
intent(INOUT) :: mfx(bd%is:bd%ie+1,bd%js:bd%je, npz)
335 real ,
intent(INOUT) :: mfy(bd%is:bd%ie ,bd%js:bd%je+1,npz)
336 real ,
intent(INOUT) :: cx(bd%is:bd%ie+1,bd%jsd:bd%jed ,npz)
337 real ,
intent(INOUT) :: cy(bd%isd:bd%ied,bd%js :bd%je +1,npz)
339 type(domain2d),
intent(INOUT) :: domain
342 real :: dp2(bd%is:bd%ie,bd%js:bd%je)
343 real :: fx(bd%is:bd%ie+1,bd%js:bd%je )
344 real :: fy(bd%is:bd%ie , bd%js:bd%je+1)
345 real :: ra_x(bd%is:bd%ie,bd%jsd:bd%jed)
346 real :: ra_y(bd%isd:bd%ied,bd%js:bd%je)
347 real :: xfx(bd%is:bd%ie+1,bd%jsd:bd%jed ,npz)
348 real :: yfx(bd%isd:bd%ied,bd%js: bd%je+1, npz)
352 integer :: ksplt(npz)
354 integer :: i,j,k,it,iq
356 real,
pointer,
dimension(:,:) :: area, rarea
357 real,
pointer,
dimension(:,:,:) :: sin_sg
358 real,
pointer,
dimension(:,:) :: dxa, dya, dx, dy
360 integer :: is, ie, js, je
361 integer :: isd, ied, jsd, jed
372 area => gridstruct%area
373 rarea => gridstruct%rarea
375 sin_sg => gridstruct%sin_sg
376 dxa => gridstruct%dxa
377 dya => gridstruct%dya
386 if (cx(i,j,k) > 0.)
then 387 xfx(i,j,k) = cx(i,j,k)*dxa(i-1,j)*dy(i,j)*sin_sg(i-1,j,3)
389 xfx(i,j,k) = cx(i,j,k)*dxa(i,j)*dy(i,j)*sin_sg(i,j,1)
395 if (cy(i,j,k) > 0.)
then 396 yfx(i,j,k) = cy(i,j,k)*dya(i,j-1)*dx(i,j)*sin_sg(i,j-1,4)
398 yfx(i,j,k) = cy(i,j,k)*dya(i,j)*dx(i,j)*sin_sg(i,j,2)
403 if ( q_split == 0 )
then 405 if ( k < npz/6 )
then 408 cmax(k) = max( cmax(k), abs(cx(i,j,k)), abs(cy(i,j,k)) )
414 cmax(k) = max( cmax(k), max(abs(cx(i,j,k)),abs(cy(i,j,k)))+1.-sin_sg(i,j,5) )
426 if ( q_split == 0 )
then 427 call mp_reduce_max(cmax,npz)
432 c_global = max(cmax(k), c_global)
435 nsplt = int(1. + c_global)
436 if ( is_master() .and. nsplt > 4 )
write(*,*)
'Tracer_2d_split=', nsplt, c_global
443 if( nsplt /= 1 )
then 451 ksplt(k) = int(1. + cmax(k))
453 frac = 1. /
real(ksplt(k))
457 cx(i,j,k) = cx(i,j,k) * frac
458 xfx(i,j,k) = xfx(i,j,k) * frac
463 mfx(i,j,k) = mfx(i,j,k) * frac
469 cy(i,j,k) = cy(i,j,k) * frac
470 yfx(i,j,k) = yfx(i,j,k) * frac
475 mfy(i,j,k) = mfy(i,j,k) * frac
485 call complete_group_halo_update(q_pack, domain)
494 if ( it .le. ksplt(k) )
then 498 dp2(i,j) = dp1(i,j,k) + (mfx(i,j,k)-mfx(i+1,j,k)+mfy(i,j,k)-mfy(i,j+1,k))*rarea(i,j)
504 ra_x(i,j) = area(i,j) + xfx(i,j,k) - xfx(i+1,j,k)
509 ra_y(i,j) = area(i,j) + yfx(i,j,k) - yfx(i,j+1,k)
514 if ( it==1 .and. trdm>1.e-4 )
then 515 call fv_tp_2d(q(isd,jsd,k,iq), cx(is,jsd,k), cy(isd,js,k), &
516 npx, npy, hord, fx, fy, xfx(is,jsd,k), yfx(isd,js,k), &
517 gridstruct, bd, ra_x, ra_y, lim_fac, regional, mfx=mfx(is,js,k), mfy=mfy(is,js,k), &
518 mass=dp1(isd,jsd,k), nord=nord_tr, damp_c=trdm)
520 call fv_tp_2d(q(isd,jsd,k,iq), cx(is,jsd,k), cy(isd,js,k), &
521 npx, npy, hord, fx, fy, xfx(is,jsd,k), yfx(isd,js,k), &
522 gridstruct, bd, ra_x, ra_y, lim_fac, regional, mfx=mfx(is,js,k), mfy=mfy(is,js,k))
526 q(i,j,k,iq) = ( q(i,j,k,iq)*dp1(i,j,k) + &
527 (fx(i,j)-fx(i+1,j)+fy(i,j)-fy(i,j+1))*rarea(i,j) )/dp2(i,j)
532 if ( it /= nsplt )
then 535 dp1(i,j,k) = dp2(i,j)
544 if ( it /= nsplt )
then 547 call start_group_halo_update(q_pack, q, domain)
558 subroutine tracer_2d_nested(q, dp1, mfx, mfy, cx, cy, gridstruct, bd, domain, npx, npy, npz, &
559 nq, hord, q_split, dt, id_divg, q_pack, nord_tr, trdm, &
560 k_split, neststruct, parent_grid, lim_fac, regional)
563 integer,
intent(IN) :: npx
564 integer,
intent(IN) :: npy
565 integer,
intent(IN) :: npz
566 integer,
intent(IN) :: nq
567 integer,
intent(IN) :: hord, nord_tr
568 integer,
intent(IN) :: q_split, k_split
569 integer,
intent(IN) :: id_divg
570 real ,
intent(IN) :: dt, trdm
571 real ,
intent(IN) :: lim_fac
572 logical,
intent(IN) :: regional
573 type(group_halo_update_type),
intent(inout) :: q_pack
574 real ,
intent(INOUT) :: q(bd%isd:bd%ied,bd%jsd:bd%jed,npz,nq)
575 real ,
intent(INOUT) :: dp1(bd%isd:bd%ied,bd%jsd:bd%jed,npz)
576 real ,
intent(INOUT) :: mfx(bd%is:bd%ie+1,bd%js:bd%je, npz)
577 real ,
intent(INOUT) :: mfy(bd%is:bd%ie ,bd%js:bd%je+1,npz)
578 real ,
intent(INOUT) :: cx(bd%is:bd%ie+1,bd%jsd:bd%jed ,npz)
579 real ,
intent(INOUT) :: cy(bd%isd:bd%ied,bd%js :bd%je +1,npz)
583 type(domain2d),
intent(INOUT) :: domain
586 real :: dp2(bd%is:bd%ie,bd%js:bd%je)
587 real :: fx(bd%is:bd%ie+1,bd%js:bd%je )
588 real :: fy(bd%is:bd%ie , bd%js:bd%je+1)
589 real :: ra_x(bd%is:bd%ie,bd%jsd:bd%jed)
590 real :: ra_y(bd%isd:bd%ied,bd%js:bd%je)
591 real :: xfx(bd%is:bd%ie+1,bd%jsd:bd%jed ,npz)
592 real :: yfx(bd%isd:bd%ied,bd%js: bd%je+1, npz)
597 real :: recip_nsplt,reg_bc_update_time
598 integer :: nsplt, nsplt_parent, msg_split_steps = 1
599 integer :: i,j,k,it,iq
601 real,
pointer,
dimension(:,:) :: area, rarea
602 real,
pointer,
dimension(:,:,:) :: sin_sg
603 real,
pointer,
dimension(:,:) :: dxa, dya, dx, dy
605 integer :: is, ie, js, je
606 integer :: isd, ied, jsd, jed
617 area => gridstruct%area
618 rarea => gridstruct%rarea
620 sin_sg => gridstruct%sin_sg
621 dxa => gridstruct%dxa
622 dya => gridstruct%dya
631 if (cx(i,j,k) > 0.)
then 632 xfx(i,j,k) = cx(i,j,k)*dxa(i-1,j)*dy(i,j)*sin_sg(i-1,j,3)
634 xfx(i,j,k) = cx(i,j,k)*dxa(i,j)*dy(i,j)*sin_sg(i,j,1)
640 if (cy(i,j,k) > 0.)
then 641 yfx(i,j,k) = cy(i,j,k)*dya(i,j-1)*dx(i,j)*sin_sg(i,j-1,4)
643 yfx(i,j,k) = cy(i,j,k)*dya(i,j)*dx(i,j)*sin_sg(i,j,2)
650 if ( q_split == 0 )
then 661 cmax_t = max( abs(cx(i,j,k)), abs(cy(i,j,k)) )
662 cmax(k) = max( cmax_t, cmax(k) )
668 cmax_t = max(abs(cx(i,j,k)), abs(cy(i,j,k))) + 1.-sin_sg(i,j,5)
669 cmax(k) = max( cmax_t, cmax(k) )
674 call mp_reduce_max(cmax,npz)
680 c_global = max(cmax(k), c_global)
683 nsplt = int(1. + c_global)
684 if ( is_master() .and. nsplt > 3 )
write(*,*)
'Tracer_2d_split=', nsplt, c_global
687 if (gridstruct%nested .and. neststruct%nestbctype > 1) msg_split_steps = max(q_split/parent_grid%flagstruct%q_split,1)
692 frac = 1. /
real(nsplt)
693 recip_nsplt = 1. /
real(nsplt)
695 if( nsplt /= 1 )
then 700 cx(i,j,k) = cx(i,j,k) * frac
701 xfx(i,j,k) = xfx(i,j,k) * frac
706 mfx(i,j,k) = mfx(i,j,k) * frac
712 cy(i,j,k) = cy(i,j,k) * frac
713 yfx(i,j,k) = yfx(i,j,k) * frac
719 mfy(i,j,k) = mfy(i,j,k) * frac
727 if ( gridstruct%nested )
then 728 neststruct%tracer_nest_timestep = neststruct%tracer_nest_timestep + 1
732 call complete_group_halo_update(q_pack, domain)
736 if (gridstruct%nested)
then 739 0, 0, npx, npy, npz, bd, &
740 real(neststruct%tracer_nest_timestep)+
real(nsplt*k_split),
real(nsplt*k_split), &
741 neststruct%q_BC(iq), bctype=neststruct%nestbctype )
748 call regional_boundary_update(q(:,:,:,iq),
'q', &
749 isd, ied, jsd, jed, npz, &
751 isd, ied, jsd, jed, &
752 reg_bc_update_time, &
765 dp2(i,j) = dp1(i,j,k) + (mfx(i,j,k)-mfx(i+1,j,k)+mfy(i,j,k)-mfy(i,j+1,k))*rarea(i,j)
771 ra_x(i,j) = area(i,j) + xfx(i,j,k) - xfx(i+1,j,k)
776 ra_y(i,j) = area(i,j) + yfx(i,j,k) - yfx(i,j+1,k)
781 if ( it==1 .and. trdm>1.e-4 )
then 782 call fv_tp_2d(q(isd,jsd,k,iq), cx(is,jsd,k), cy(isd,js,k), &
783 npx, npy, hord, fx, fy, xfx(is,jsd,k), yfx(isd,js,k), &
784 gridstruct, bd, ra_x, ra_y, lim_fac, regional, mfx=mfx(is,js,k), mfy=mfy(is,js,k), &
785 mass=dp1(isd,jsd,k), nord=nord_tr, damp_c=trdm)
787 call fv_tp_2d(q(isd,jsd,k,iq), cx(is,jsd,k), cy(isd,js,k), &
788 npx, npy, hord, fx, fy, xfx(is,jsd,k), yfx(isd,js,k), &
789 gridstruct, bd, ra_x, ra_y, lim_fac, regional, mfx=mfx(is,js,k), mfy=mfy(is,js,k))
793 q(i,j,k,iq) = ( q(i,j,k,iq)*dp1(i,j,k) + &
794 (fx(i,j)-fx(i+1,j)+fy(i,j)-fy(i,j+1))*rarea(i,j) )/dp2(i,j)
800 if ( it /= nsplt )
then 803 call start_group_halo_update(q_pack, q, domain)
808 if ( gridstruct%nested )
then 813 0, 0, npx, npy, npz, bd, &
814 real(neststruct%tracer_nest_timestep),
real(nsplt*k_split), &
815 neststruct%q_BC(iq), bctype=neststruct%nestbctype )
824 if ( id_divg > 0 )
then 831 dp1(i,j,k) = (xfx(i+1,j,k)-xfx(i,j,k) + yfx(i,j+1,k)-yfx(i,j,k))*rarea(i,j)*rdt
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...
real, dimension(:,:,:), allocatable nest_fx_west_accum
real, dimension(:,:,:), allocatable nest_fx_south_accum
The module 'fv_timing' contains FV3 timers.
The module 'boundary' contains utility routines for grid nesting and boundary conditions.
The module 'tp_core' is a collection of routines to support FV transport.
The module 'fv_tracer2d.F90' performs sub-cycled tracer advection.
subroutine, public copy_corners(q, npx, npy, dir, nested, bd, sw_corner, se_corner, nw_corner, ne_corner)
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 timing_on(blk_name)
The subroutine 'timing_on' starts a timer.
integer, parameter, public ng
real, dimension(:,:,:), allocatable nest_fx_east_accum
subroutine, public fv_tp_2d(q, crx, cry, npx, npy, hord, fx, fy, xfx, yfx, gridstruct, bd, ra_x, ra_y, lim_fac, regional, mfx, mfy, mass, nord, damp_c)
The subroutine 'fv_tp_2d' contains the FV advection scheme .
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...
real, dimension(:,:,:), allocatable nest_fx_north_accum
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)