mus_Poisson_advRel_d2q9 Subroutine

public subroutine mus_Poisson_advRel_d2q9(fieldProp, inState, outState, auxField, neigh, nElems, nSolve, level, layout, params, varSys, derVarPos)

Advection relaxation routine for the poisson equation with an explicit calculation of all equilibrium quantities. Slow and simple. The right hand side of equation is added as a source term in mus_apply_sourceTerms routine

This subroutine interface must match the abstract interface definition kernel in scheme/mus_scheme_type_module.f90 in order to be callable via compute function pointer.

\todo KM 20170821: change omega for multilevel Generic fetching step: Streaming for pull Local copy for push

Arguments

TypeIntentOptionalAttributesName
type(mus_field_prop_type), intent(in) :: fieldProp(:)

Array of field properties (fluid or species)

real(kind=rk), intent(in) :: inState(nElems*varSys%nScalars)

input pdf vector

real(kind=rk), intent(out) :: outState(nElems*varSys%nScalars)

output pdf vector

real(kind=rk), intent(inout) :: auxField(nElems*varSys%nAuxScalars)

Auxiliary field computed from pre-collision state Is updated with correct velocity field for multicomponent models

integer, intent(in) :: neigh(nElems*layout%fStencil%QQ)

connectivity vector

integer, intent(in) :: nElems

number of elements in state Array

integer, intent(in) :: nSolve

number of elements solved in kernel

integer, intent(in) :: level

current level

type(mus_scheme_layout_type), intent(in) :: layout

current layout

type(mus_param_type), intent(in) :: params

global parameters

type(tem_varSys_type), intent(in) :: varSys

variable system definition

type(mus_derVarPos_type), intent(in) :: derVarPos(:)

position of derived quantities in varsys for all fields


Contents


Variables

TypeVisibilityAttributesNameInitial
integer, private :: iElem
real(kind=rk), private :: pdfTmp(QQ_9)
real(kind=rk), private :: omega
real(kind=rk), private :: omega_fac
real(kind=rk), private :: om_pot
real(kind=rk), private :: fac14
real(kind=rk), private :: fac58
real(kind=rk), private :: fac9