INTERFACE:
subroutine do_advection(dt,f,U,V,DU,DV,Do,Dn,split,scheme,AH,tag, & Dires,advres)DESCRIPTION:
Laterally advects a 2D quantity. The location of the quantity on the grid (either T-, U- or V-points) must be specified by the argument tag. The transports through the interfaces of the corresponding Finite-Volumes and their different height information (all relative to the given quantity) must be provided as well. Depending on split and scheme several fractional steps (Strang splitting) with different options for the calculation of the interfacial fluxes are carried out.
The options for split are:
split = NOSPLIT: | no split (one 2D uv step) |
split = FULLSPLIT: | full step splitting (u + v) |
split = HALFSPLIT: | half step splitting (u/2 + v + u/2) |
The options for scheme are:
scheme = NOADV: | advection disabled |
scheme = UPSTREAM: | first-order upstream (monotone) |
scheme = UPSTREAM_2DH: | 2DH upstream with forced monotonicity |
scheme = P2: | third-order polynomial (non-monotone) |
scheme = SUPERBEE: | second-order TVD (monotone) |
scheme = MUSCL: | second-order TVD (monotone) |
scheme = P2_PDM: | third-order ULTIMATE-QUICKEST (monotone) |
scheme = J7: | 2DH Arakawa J7 |
scheme = FCT: | 2DH FCT with forced monotonicity |
scheme = P2_2DH: | 2DH P2 with forced monotonicity |
With the compiler option SLICE_MODEL, the advection in meridional direction is not executed.
USES:
use halo_zones, only: update_2d_halo,wait_halo,D_TAG,H_TAG,U_TAG,V_TAG use getm_timers, only: tic,toc,TIM_ADV,TIM_ADVH IMPLICIT NONEINPUT PARAMETERS:
REALTYPE,intent(in) :: dt,AH REALTYPE,dimension(E2DFIELD),intent(in) :: U,V,Do,Dn,DU,DV integer,intent(in) :: split,scheme,tagINPUT/OUTPUT PARAMETERS:
REALTYPE,dimension(E2DFIELD),intent(inout) :: fOUTPUT PARAMETERS:
REALTYPE,dimension(E2DFIELD),target,intent(out),optional :: Dires,advresLOCAL VARIABLES:
type(t_adv_grid),pointer :: adv_grid REALTYPE,dimension(E2DFIELD),target :: fi,Di,adv REALTYPE,dimension(:,:),pointer,contiguous :: p_Di,p_adv integer :: i,j