mus_addHRRCorrToAuxField_fluid_D3Q19 Subroutine

public subroutine mus_addHRRCorrToAuxField_fluid_D3Q19(fun, auxField, iLevel, time, varSys, phyConvFac, derVarPos)

This routine add sponge density and velocity field to density and velocity in auxField. Density and velocity in far field are computed by time average.

Reference: Jacob, J.; Sagaut, P. (2019): Solid wall and open boundary conditions in hybrid recursive regularized lattice Boltzmann method for compressible flows. In Physics of Fluids 31 (12), p. 126103.

Arguments

TypeIntentOptionalAttributesName
class(mus_source_op_type), intent(inout) :: fun

Description of method to update source

real(kind=rk), intent(inout) :: auxField(:)

output auxField array

integer, intent(in) :: iLevel

current level

type(tem_time_type), intent(in) :: time

current timing information

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

variable system definition

type(mus_convertFac_type), intent(in) :: phyConvFac

Physics conversion factor for current level

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

position of derived quantities in varsys


Contents


Variables

TypeVisibilityAttributesNameInitial
type(mus_gradData_type), private, pointer:: gradData

gradient data

type(mus_varSys_data_type), private, pointer:: fPtr
type(mus_scheme_type), private, pointer:: scheme
type(mus_scheme_layout_type), private, pointer:: layout
integer, private :: dens_pos
integer, private :: vel_pos(3)
integer, private :: QQ
integer, private :: iDir
integer, private :: nSolve
integer, private :: nDims
integer, private :: iElem
integer, private :: elemOff
integer, private :: nChunks
integer, private :: iChunks
integer, private :: nChunkElems
integer, private :: low_bound
integer, private :: elemPos
real(kind=rk), private :: gradRHOU3(3,vlen)
real(kind=rk), private :: SCorr
real(kind=rk), private :: gradRHOUVZ(3,vlen)
real(kind=rk), private :: H2xx
real(kind=rk), private :: H2xy
real(kind=rk), private :: H2yy
real(kind=rk), private :: H2zz
real(kind=rk), private :: H2yz
real(kind=rk), private :: H2xz
real(kind=rk), private :: c_x
real(kind=rk), private :: c_y
real(kind=rk), private :: c_z
real(kind=rk), private :: dens
real(kind=rk), private :: vel(3)