INTERFACE:
subroutine adv_split_u(dt,f,fi,Di,adv,U,DU, & #if defined(SPHERICAL) || defined(CURVILINEAR) dxu,dyu,arcd1, & #endif splitfac,scheme,AH, & mask_flux,mask_update) Note (KK): Keep in sync with interface in advection.F90DESCRIPTION:
Executes an advection step in zonal direction for a 2D quantity. The 1D advection equation
is accompanied by an fractional step for the 1D continuity equation
Here, and denote values before and after this operation, respectively, denote intermediate values when other 1D advection steps come after this and denotes intermediate values when other 1D advection steps came before this. Furthermore, when this -directional split step is repeated during the total time step (Strang splitting), the time step denotes a fraction of the full time step.
The interfacial fluxes are calculated according to the third-order polynomial scheme (so-called P scheme), cast in Lax-Wendroff form by:
with the Courant number and
where
(76) |
and
(77) |
It should be noted, that for the original Lax-Wendroff scheme and for the first-oder upstream scheme can be recovered.
In order to obtain a monotonic and positive scheme, the factors are limited in the following way:
and, equivalently, for . This so-called PDM-limiter has been described in detail by Leonard (1991), who named the PDM-limited P scheme also ULTIMATE QUICKEST (quadratic upstream interpolation for convective kinematics with estimated stream terms).
Some simpler limiters which do not exploit the third-order polynomial properties of the discretisation (74) have been listed by Zalezak (1987). Among those are the MUSCL scheme by van Leer (1979),
and the Superbee scheme by Roe (1985),
The selector for the schemes is scheme:
scheme = UPSTREAM: | first-order upstream (monotone) |
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) |
Furthermore, the horizontal diffusion in zonal direction with the constant diffusion coefficient AH is carried out here by means of a central difference second-order scheme. USES:
use domain, only: imin,imax,jmin,jmax #if !( defined(SPHERICAL) || defined(CURVILINEAR) ) use domain, only: dx,dy,ard1 #endif use advection, only: adv_interfacial_reconstruction use advection, only: UPSTREAM $ use omp_lib IMPLICIT NONEINPUT PARAMETERS:
Note (KK): in general dxu, dyu and mask_flux do only have valid data within (_IRANGE_HALO_-1,_JRANGE_HALO_). In some cases the original field extension may even be _IRANGE_HALO_. Then explicit declared array bounds _IRANGE_HALO_-1 require a provision of the corresponding subarray and will cause copying of the non-contiguously data into a temporarily array. Therefore they are declared as pointers here. This however requires, that the provided pointers already carry the correct bounds. REALTYPE,intent(in) :: dt,splitfac,AH REALTYPE,dimension(E2DFIELD),intent(in) :: f,U,DU #if defined(SPHERICAL) || defined(CURVILINEAR) REALTYPE,dimension(:,:),pointer,contiguous,intent(in) :: dxu,dyu REALTYPE,dimension(E2DFIELD),intent(in) :: arcd1 #endif integer,intent(in) :: scheme logical,dimension(:,:),pointer,contiguous,intent(in) :: mask_flux logical,dimension(E2DFIELD),intent(in) :: mask_updateINPUT/OUTPUT PARAMETERS:
REALTYPE,dimension(E2DFIELD),intent(inout) :: fi,Di,advLOCAL VARIABLES:
REALTYPE,dimension(E2DFIELD) :: uflux logical :: use_limiter,use_AH integer :: i,j,isub REALTYPE :: dti,Dio,advn,cfl,fuu,fu,fdREVISION HISTORY:
Original author(s): Hans Burchard & Karsten Bolding