INTERFACE:
subroutine do_advection_3d(dt,f,uu,vv,ww,hu,hv,ho,hn, & split,hscheme,vscheme,AH,tag, & hires,advres)DESCRIPTION:
Here, advection terms for all three-dimensional state variables are calculated by means of a finite-volume approach (an exception is the possibility to directly calculate the momentum advection by a one-step three-dimensional upstream scheme, see uv_advection_3d) and the advection step is carried out as a fractional advection time step. Those 3D variables may be defined on T-, U-, V- and W-points. The latter option is interesting for turbulent quantities. Inside this advection routine, it does not matter, wehre the advected variable is located on the grid. All finite volume fluxes and geometric coefficients need to be calculated before do_advection_3d is called.
Originally, this 3D advection routine has been written for tracer equations. There, after multiplying the layer-integrated and transformed to curvilinear coordinates tracer equation (40) with , the advective terms in this equation are discretised as follows.
First advection term in (40):
Second advection term in (40):
Vertical advective fluxes in (40):
The interfacial concentrations are calculated according to upwind or higher order directional split schemes, which are discussed in detail below and in sections 7.4.2 and 8.6.4.
However, as said above, in the same way these routines may be applied to quantities on U-, V-, and W-points, if the transports are properly calculated.
There are various combinations of advection schemes possible.
The options for split are:
split = NOSPLIT: | no split (one 3D uvw step) |
split = FULLSPLIT: | full step splitting (u + v + w) |
split = HALFSPLIT: | half step splitting (u/2 + v/2 + w + v/2 + u/2) |
split = HVSPLIT: | hor./ver. splitting (uv + w) |
The options for the horizontal scheme hscheme 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 |
The options for the vertical scheme vscheme are:
scheme = NOADV: | advection disabled |
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) |
With the compiler option SLICE_MODEL, the advection in meridional direction is not executed.
USES:
use halo_zones, only: update_3d_halo,wait_halo,D_TAG,H_TAG,U_TAG,V_TAG use getm_timers, only: tic,toc,TIM_ADV3D,TIM_ADV3DH IMPLICIT NONEINPUT PARAMETERS:
REALTYPE,intent(in) :: dt,AH REALTYPE,dimension(I3DFIELD),intent(in) :: uu,vv,ww,ho,hn,hu,hv integer,intent(in) :: split,hscheme,vscheme,tagINPUT/OUTPUT PARAMETERS:
REALTYPE,dimension(I3DFIELD),intent(inout) :: fOUTPUT PARAMETERS:
REALTYPE,dimension(I3DFIELD),target,intent(out),optional :: hires,advresLOCAL VARIABLES:
type(t_adv_grid),pointer :: adv_grid REALTYPE,dimension(I3DFIELD),target :: fi,hi,adv REALTYPE,dimension(:,:,:),pointer,contiguous :: p_hi,p_adv integer :: tag2d,i,j,k