INTERFACE:
module advectionDESCRIPTION:
This module does lateral advection of scalars. It follows the same convention as the other modules in 'getm'. The module is initialised by calling 'init_advection()'. In the time-loop 'do_advection()' is called. 'do_advection' is a wrapper routine which - dependent on the actual advection scheme chosen - makes calls to the appropriate subroutines, which may be done as one-step or multiple-step schemes. The actual subroutines are coded in external FORTRAN files. New advection schemes are easily implemented - at least from a program point of view - since only this module needs to be changed. Additional work arrays can easily be added following the stencil given below. To add a new advection scheme three things must be done:
use domain, only: imin,imax,jmin,jmax
IMPLICIT NONE
private
PUBLIC DATA MEMBERS:
public init_advection,do_advection,print_adv_settings
public adv_split_u,adv_split_v,adv_upstream_2dh,adv_arakawa_j7_2dh,adv_fct_2dh
public adv_interfacial_reconstruction
type, public :: t_adv_grid
logical,dimension(:,:),pointer,contiguous :: mask_uflux,mask_vflux,mask_xflux
logical,dimension(:,:),pointer,contiguous :: mask_uupdate,mask_vupdate
logical,dimension(:,:),pointer,contiguous :: mask_finalise
integer,dimension(:,:),pointer,contiguous :: az
#if defined(SPHERICAL) || defined(CURVILINEAR)
REALTYPE,dimension(:,:),pointer,contiguous :: dxu,dyu,dxv,dyv,arcd1
#endif
end type t_adv_grid
type(t_adv_grid),public,target :: adv_gridH,adv_gridU,adv_gridV
integer,public,parameter :: NOSPLIT=0,FULLSPLIT=1,HALFSPLIT=2
character(len=64),public,parameter :: adv_splits(0:2) = &
(/"no split: one 2D uv step ", &
"full step splitting: u + v ", &
"half step splitting: u/2 + v + u/2"/)
integer,public,parameter :: NOADV=0,UPSTREAM=1,UPSTREAM_2DH=2
integer,public,parameter :: P2=3,SUPERBEE=4,MUSCL=5,P2_PDM=6
integer,public,parameter :: J7=7,FCT=8,P2_2DH=9
character(len=64),public,parameter :: adv_schemes(0:9) = &
(/"advection disabled ", &
"upstream advection (first-order, monotone) ", &
"2DH-upstream advection with forced monotonicity", &
"P2 advection (third-order, non-monotone) ", &
"TVD-Superbee advection (second-order, monotone)", &
"TVD-MUSCL advection (second-order, monotone) ", &
"TVD-P2-PDM advection (third-order, monotone) ", &
"2DH-J7 advection (Arakawa and Lamb, 1977) ", &
"2DH-FCT advection ", &
"2DH-P2 advection "/)
LOCAL VARIABLES:
#ifdef STATIC
logical,dimension(E2DFIELD),target :: mask_updateH
logical,dimension(E2DFIELD),target :: mask_uflux,mask_vflux,mask_xflux
logical,dimension(E2DFIELD),target :: mask_uupdateU,mask_vupdateV
#else
logical,dimension(:,:),allocatable,target :: mask_updateH
logical,dimension(:,:),allocatable,target :: mask_uflux,mask_vflux,mask_xflux
logical,dimension(:,:),allocatable,target :: mask_uupdateU,mask_vupdateV
#endif
REVISION HISTORY:
Original author(s): Knut Klingbeil