applySrc_absorbLayer Subroutine

public subroutine applySrc_absorbLayer(fun, inState, outState, neigh, auxField, nPdfSize, iLevel, varSys, time, phyConvFac, derVarPos)

Update state with source variable "absorb_layer". absorb_layer is used to absorb the flow and gradually reduce the flow quantities like pressure and velocity to a fixed value. It is based on: Xu, H., & Sagaut, P. (2013). Analysis of the absorbing layers for the weakly-compressible lattice Boltzmann methods. Journal of Computational Physics, 245(x), 14–42. 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.

This subroutine's interface must match the abstract interface definition proc_apply_source in derived/mus_source_type_module.f90 in order to be callable via applySrc function pointer.

Arguments

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

Description of method to apply source terms

real(kind=rk), intent(in) :: inState(:)

input pdf vector

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

output pdf vector

integer, intent(in) :: neigh(:)

connectivity Array corresponding to state vector

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

auxField array

integer, intent(in) :: nPdfSize

number of elements in state Array

integer, intent(in) :: iLevel

current level

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

variable system

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

Point in time at which to evaluate the variable.

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_varSys_data_type), private, pointer:: fPtr
type(mus_scheme_type), private, pointer:: scheme
real(kind=rk), private :: spongeField(fun%elemLvl(iLevel)%nElems*5)
real(kind=rk), private :: dens
real(kind=rk), private :: vel(3)
real(kind=rk), private :: ucx
real(kind=rk), private :: uMinusCX(3)
real(kind=rk), private :: sponge_velTerm
integer, private :: nElems
integer, private :: iElem
integer, private :: iDir
integer, private :: QQ
integer, private :: nScalars
integer, private :: posInTotal
integer, private :: statePos
integer, private :: elemOff
integer, private :: vel_pos(3)
integer, private :: dens_pos
real(kind=rk), private :: sigma
real(kind=rk), private :: dens_ref
real(kind=rk), private :: vel_ref(3)
real(kind=rk), private :: sponge_vel(3)
real(kind=rk), private :: sponge_dens
real(kind=rk), private :: inv_rho_phy
real(kind=rk), private :: inv_vel_phy
real(kind=rk), private :: omega
real(kind=rk), private :: omega_fac