waLBerla 7.2
Loading...
Searching...
No Matches
CellwiseSweep.h File Reference

Detailed Description

Classes

class  walberla::lbm::CellwiseSweep< LatticeModel_T, Filter_T, DensityVelocityIn_T, DensityVelocityOut_T, Enable >
 Class for performing the streaming and collision step of the LBM. More...
 

Namespaces

namespace  walberla
 Storage for detected contacts which can be used to perform actions for all contacts, e.g.
 
namespace  walberla::lbm
 

Macros

#define WALBERLA_LBM_CELLWISE_SWEEP_CLASS_HEAD_AND_STREAM(specialization)
 
#define WALBERLA_LBM_CELLWISE_SWEEP_STREAM_COLLIDE_HEAD(specialization)
 
#define WALBERLA_LBM_CELLWISE_SWEEP_STREAM_COLLIDE_FOOT()   src->swapDataPointers( dst ); }
 
#define WALBERLA_LBM_CELLWISE_SWEEP_COLLIDE_HEAD(specialization)
 
#define WALBERLA_LBM_CELLWISE_SWEEP_COLLIDE_FOOT()   }
 
#define WALBERLA_LBM_CELLWISE_SWEEP_D2Q9_STREAM_COLLIDE_PULL()
 
#define WALBERLA_LBM_CELLWISE_SWEEP_D2Q9_COLLIDE_GET()
 
#define WALBERLA_LBM_CELLWISE_SWEEP_D2Q9_DENSITY_VELOCITY_INCOMP()
 
#define WALBERLA_LBM_CELLWISE_SWEEP_D2Q9_DENSITY_VELOCITY_COMP()
 
#define WALBERLA_LBM_CELLWISE_SWEEP_D3Q19_STREAM_COLLIDE_PULL()
 
#define WALBERLA_LBM_CELLWISE_SWEEP_D3Q19_COLLIDE_GET()
 
#define WALBERLA_LBM_CELLWISE_SWEEP_D3Q19_DENSITY_VELOCITY_INCOMP()
 
#define WALBERLA_LBM_CELLWISE_SWEEP_D3Q19_DENSITY_VELOCITY_COMP()
 
#define WALBERLA_LBM_CELLWISE_SWEEP_D3Q27_STREAM_COLLIDE_PULL()
 
#define WALBERLA_LBM_CELLWISE_SWEEP_D3Q27_COLLIDE_GET()
 
#define WALBERLA_LBM_CELLWISE_SWEEP_D3Q27_DENSITY_VELOCITY_INCOMP()
 
#define WALBERLA_LBM_CELLWISE_SWEEP_D3Q27_DENSITY_VELOCITY_COMP()
 

Functions

template<typename LatticeModel_T , typename Filter_T , typename DensityVelocityIn_T , typename DensityVelocityOut_T >
shared_ptr< CellwiseSweep< LatticeModel_T, Filter_T, DensityVelocityIn_T, DensityVelocityOut_T > > walberla::lbm::makeCellwiseSweep (const BlockDataID &pdfFieldId, const Filter_T &filter, const DensityVelocityIn_T &densityVelocityIn, const DensityVelocityOut_T &densityVelocityOut)
 
template<typename LatticeModel_T , typename Filter_T , typename DensityVelocityIn_T , typename DensityVelocityOut_T >
shared_ptr< CellwiseSweep< LatticeModel_T, Filter_T, DensityVelocityIn_T, DensityVelocityOut_T > > walberla::lbm::makeCellwiseSweep (const BlockDataID &src, const BlockDataID &dst, const Filter_T &filter, const DensityVelocityIn_T &densityVelocityIn, const DensityVelocityOut_T &densityVelocityOut)
 
template<typename LatticeModel_T >
shared_ptr< CellwiseSweep< LatticeModel_T, walberla::field::DefaultEvaluationFilter, DefaultDensityEquilibriumVelocityCalculation, DefaultDensityVelocityCallback > > walberla::lbm::makeCellwiseSweep (const BlockDataID &pdfFieldId)
 
template<typename LatticeModel_T >
shared_ptr< CellwiseSweep< LatticeModel_T, walberla::field::DefaultEvaluationFilter, DefaultDensityEquilibriumVelocityCalculation, DefaultDensityVelocityCallback > > walberla::lbm::makeCellwiseSweep (const BlockDataID &src, const BlockDataID &dst)
 
template<typename LatticeModel_T , typename FlagField_T >
shared_ptr< CellwiseSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, DefaultDensityEquilibriumVelocityCalculation, DefaultDensityVelocityCallback > > walberla::lbm::makeCellwiseSweep (const BlockDataID &pdfFieldId, const ConstBlockDataID &flagFieldId, const Set< FlagUID > &cellsToEvaluate)
 
template<typename LatticeModel_T , typename FlagField_T >
shared_ptr< CellwiseSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, DefaultDensityEquilibriumVelocityCalculation, DefaultDensityVelocityCallback > > walberla::lbm::makeCellwiseSweep (const BlockDataID &src, const BlockDataID &dst, const ConstBlockDataID &flagFieldId, const Set< FlagUID > &cellsToEvaluate)
 
template<typename LatticeModel_T , typename FlagField_T , typename VelocityField_T >
shared_ptr< CellwiseSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, DefaultDensityEquilibriumVelocityCalculation, VelocityCallback< VelocityField_T > > > walberla::lbm::makeCellwiseSweep (const BlockDataID &pdfFieldId, const ConstBlockDataID &flagFieldId, const Set< FlagUID > &cellsToEvaluate, const BlockDataID &velocityFieldId)
 
template<typename LatticeModel_T , typename FlagField_T , typename VelocityField_T >
shared_ptr< CellwiseSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, DefaultDensityEquilibriumVelocityCalculation, VelocityCallback< VelocityField_T > > > walberla::lbm::makeCellwiseSweep (const BlockDataID &src, const BlockDataID &dst, const ConstBlockDataID &flagFieldId, const Set< FlagUID > &cellsToEvaluate, const BlockDataID &velocityFieldId)
 
template<typename LatticeModel_T , typename FlagField_T , typename VelocityField_T , typename DensityField_T >
shared_ptr< CellwiseSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, DefaultDensityEquilibriumVelocityCalculation, DensityVelocityCallback< VelocityField_T, DensityField_T > > > walberla::lbm::makeCellwiseSweep (const BlockDataID &pdfFieldId, const ConstBlockDataID &flagFieldId, const Set< FlagUID > &cellsToEvaluate, const BlockDataID &velocityFieldId, const BlockDataID &densityFieldId)
 
template<typename LatticeModel_T , typename FlagField_T , typename VelocityField_T , typename DensityField_T >
shared_ptr< CellwiseSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, DefaultDensityEquilibriumVelocityCalculation, DensityVelocityCallback< VelocityField_T, DensityField_T > > > walberla::lbm::makeCellwiseSweep (const BlockDataID &src, const BlockDataID &dst, const ConstBlockDataID &flagFieldId, const Set< FlagUID > &cellsToEvaluate, const BlockDataID &velocityFieldId, const BlockDataID &densityFieldId)
 
template<typename LatticeModel_T , typename VelocityField_T , typename Filter_T , typename DensityVelocityOut_T >
shared_ptr< CellwiseSweep< LatticeModel_T, Filter_T, AdvectionDiffusionDensityEquilibriumVelocityCalculation< VelocityField_T >, DensityVelocityOut_T > > walberla::lbm::makeCellwiseAdvectionDiffusionSweep (const BlockDataID &pdfFieldId, const ConstBlockDataID &velocityFieldId, const Filter_T &filter, const DensityVelocityOut_T &densityVelocityOut)
 
template<typename LatticeModel_T , typename VelocityField_T , typename Filter_T , typename DensityVelocityOut_T >
shared_ptr< CellwiseSweep< LatticeModel_T, Filter_T, AdvectionDiffusionDensityEquilibriumVelocityCalculation< VelocityField_T >, DensityVelocityOut_T > > walberla::lbm::makeCellwiseAdvectionDiffusionSweep (const BlockDataID &src, const BlockDataID &dst, const ConstBlockDataID &velocityFieldId, const Filter_T &filter, const DensityVelocityOut_T &densityVelocityOut)
 
template<typename LatticeModel_T , typename VelocityField_T >
shared_ptr< CellwiseSweep< LatticeModel_T, walberla::field::DefaultEvaluationFilter, AdvectionDiffusionDensityEquilibriumVelocityCalculation< VelocityField_T >, DefaultDensityVelocityCallback > > walberla::lbm::makeCellwiseAdvectionDiffusionSweep (const BlockDataID &pdfFieldId, const ConstBlockDataID &velocityFieldId)
 
template<typename LatticeModel_T , typename VelocityField_T >
shared_ptr< CellwiseSweep< LatticeModel_T, walberla::field::DefaultEvaluationFilter, AdvectionDiffusionDensityEquilibriumVelocityCalculation< VelocityField_T >, DefaultDensityVelocityCallback > > walberla::lbm::makeCellwiseAdvectionDiffusionSweep (const BlockDataID &src, const BlockDataID &dst, const ConstBlockDataID &velocityFieldId)
 
template<typename LatticeModel_T , typename VelocityField_T , typename FlagField_T >
shared_ptr< CellwiseSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, AdvectionDiffusionDensityEquilibriumVelocityCalculation< VelocityField_T >, DefaultDensityVelocityCallback > > walberla::lbm::makeCellwiseAdvectionDiffusionSweep (const BlockDataID &pdfFieldId, const ConstBlockDataID &velocityFieldId, const ConstBlockDataID &flagFieldId, const Set< FlagUID > &cellsToEvaluate)
 
template<typename LatticeModel_T , typename VelocityField_T , typename FlagField_T >
shared_ptr< CellwiseSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, AdvectionDiffusionDensityEquilibriumVelocityCalculation< VelocityField_T >, DefaultDensityVelocityCallback > > walberla::lbm::makeCellwiseAdvectionDiffusionSweep (const BlockDataID &src, const BlockDataID &dst, const ConstBlockDataID &velocityFieldId, const ConstBlockDataID &flagFieldId, const Set< FlagUID > &cellsToEvaluate)
 
template<typename LatticeModel_T , typename VelocityField_T , typename FlagField_T , typename VelocityFieldOut_T >
shared_ptr< CellwiseSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, AdvectionDiffusionDensityEquilibriumVelocityCalculation< VelocityField_T >, VelocityCallback< VelocityFieldOut_T > > > walberla::lbm::makeCellwiseAdvectionDiffusionSweep (const BlockDataID &pdfFieldId, const ConstBlockDataID &velocityFieldId, const ConstBlockDataID &flagFieldId, const Set< FlagUID > &cellsToEvaluate, const BlockDataID &velocityFieldOutId)
 
template<typename LatticeModel_T , typename VelocityField_T , typename FlagField_T , typename VelocityFieldOut_T >
shared_ptr< CellwiseSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, AdvectionDiffusionDensityEquilibriumVelocityCalculation< VelocityField_T >, VelocityCallback< VelocityFieldOut_T > > > walberla::lbm::makeCellwiseAdvectionDiffusionSweep (const BlockDataID &src, const BlockDataID &dst, const ConstBlockDataID &velocityFieldId, const ConstBlockDataID &flagFieldId, const Set< FlagUID > &cellsToEvaluate, const BlockDataID &velocityFieldOutId)
 
template<typename LatticeModel_T , typename VelocityField_T , typename FlagField_T , typename VelocityFieldOut_T , typename DensityFieldOut_T >
shared_ptr< CellwiseSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, AdvectionDiffusionDensityEquilibriumVelocityCalculation< VelocityField_T >, DensityVelocityCallback< VelocityFieldOut_T, DensityFieldOut_T > > > walberla::lbm::makeCellwiseAdvectionDiffusionSweep (const BlockDataID &pdfFieldId, const ConstBlockDataID &velocityFieldId, const ConstBlockDataID &flagFieldId, const Set< FlagUID > &cellsToEvaluate, const BlockDataID &velocityFieldOutId, const BlockDataID &densityFieldOutId)
 
template<typename LatticeModel_T , typename VelocityField_T , typename FlagField_T , typename VelocityFieldOut_T , typename DensityFieldOut_T >
shared_ptr< CellwiseSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, AdvectionDiffusionDensityEquilibriumVelocityCalculation< VelocityField_T >, DensityVelocityCallback< VelocityFieldOut_T, DensityFieldOut_T > > > walberla::lbm::makeCellwiseAdvectionDiffusionSweep (const BlockDataID &src, const BlockDataID &dst, const ConstBlockDataID &velocityFieldId, const ConstBlockDataID &flagFieldId, const Set< FlagUID > &cellsToEvaluate, const BlockDataID &velocityFieldOutId, const BlockDataID &densityFieldOutId)
 

Macro Definition Documentation

◆ WALBERLA_LBM_CELLWISE_SWEEP_CLASS_HEAD_AND_STREAM

#define WALBERLA_LBM_CELLWISE_SWEEP_CLASS_HEAD_AND_STREAM ( specialization)

◆ WALBERLA_LBM_CELLWISE_SWEEP_COLLIDE_FOOT

#define WALBERLA_LBM_CELLWISE_SWEEP_COLLIDE_FOOT ( )    }

◆ WALBERLA_LBM_CELLWISE_SWEEP_COLLIDE_HEAD

#define WALBERLA_LBM_CELLWISE_SWEEP_COLLIDE_HEAD ( specialization)
Value:
template< typename LatticeModel_T, typename Filter_T, typename DensityVelocityIn_T, typename DensityVelocityOut_T > \
void CellwiseSweep< LatticeModel_T, Filter_T, DensityVelocityIn_T, DensityVelocityOut_T, typename std::enable_if_t< specialization > \
>::collide( IBlock * const block, const uint_t numberOfGhostLayersToInclude ) \
{ \
PdfField_T * src = this->getSrcField( block ); \
WALBERLA_ASSERT_GREATER_EQUAL( src->nrOfGhostLayers(), numberOfGhostLayersToInclude ); \
\
const auto & lm = src->latticeModel(); \
\
this->filter( *block ); \
this->densityVelocityIn( *block ); \
this->densityVelocityOut( *block );

◆ WALBERLA_LBM_CELLWISE_SWEEP_D2Q9_COLLIDE_GET

#define WALBERLA_LBM_CELLWISE_SWEEP_D2Q9_COLLIDE_GET ( )
Value:
const real_t vC = src->get( x, y, z, Stencil_T::idx[C] ); \
const real_t vN = src->get( x, y, z, Stencil_T::idx[N] ); \
const real_t vS = src->get( x, y, z, Stencil_T::idx[S] ); \
const real_t vW = src->get( x, y, z, Stencil_T::idx[W] ); \
const real_t vE = src->get( x, y, z, Stencil_T::idx[E] ); \
const real_t vNW = src->get( x, y, z, Stencil_T::idx[NW] ); \
const real_t vNE = src->get( x, y, z, Stencil_T::idx[NE] ); \
const real_t vSW = src->get( x, y, z, Stencil_T::idx[SW] ); \
const real_t vSE = src->get( x, y, z, Stencil_T::idx[SE] );

◆ WALBERLA_LBM_CELLWISE_SWEEP_D2Q9_DENSITY_VELOCITY_COMP

#define WALBERLA_LBM_CELLWISE_SWEEP_D2Q9_DENSITY_VELOCITY_COMP ( )
Value:
const real_t velXTerm = vE + vNE + vSE; \
const real_t velYTerm = vN + vNW; \
\
const real_t rho = vC + vS + vW + vB + vSW + velXTerm + velYTerm; \
const real_t invRho = real_t(1.0) / rho; \
\
const real_t velX = invRho * ( velXTerm - vW - vNW - vSW ); \
const real_t velY = invRho * ( velYTerm + vNE - vS - vSW - vSE );

◆ WALBERLA_LBM_CELLWISE_SWEEP_D2Q9_DENSITY_VELOCITY_INCOMP

#define WALBERLA_LBM_CELLWISE_SWEEP_D2Q9_DENSITY_VELOCITY_INCOMP ( )
Value:
const real_t velXTerm = vE + vNE + vSE; \
const real_t velYTerm = vN + vNW; \
\
const real_t rho = vC + vS + vW + vSW + velXTerm + velYTerm; \
\
const real_t velX = velXTerm - vW - vNW - vSW; \
const real_t velY = velYTerm + vNE - vS - vSW - vSE;

◆ WALBERLA_LBM_CELLWISE_SWEEP_D2Q9_STREAM_COLLIDE_PULL

#define WALBERLA_LBM_CELLWISE_SWEEP_D2Q9_STREAM_COLLIDE_PULL ( )
Value:
const real_t vC = src->get( x , y , z, Stencil_T::idx[C] ); \
const real_t vN = src->get( x , y-1, z, Stencil_T::idx[N] ); \
const real_t vS = src->get( x , y+1, z, Stencil_T::idx[S] ); \
const real_t vW = src->get( x+1, y , z, Stencil_T::idx[W] ); \
const real_t vE = src->get( x-1, y , z, Stencil_T::idx[E] ); \
const real_t vNW = src->get( x+1, y-1, z, Stencil_T::idx[NW] ); \
const real_t vNE = src->get( x-1, y-1, z, Stencil_T::idx[NE] ); \
const real_t vSW = src->get( x+1, y+1, z, Stencil_T::idx[SW] ); \
const real_t vSE = src->get( x-1, y+1, z, Stencil_T::idx[SE] );

◆ WALBERLA_LBM_CELLWISE_SWEEP_D3Q19_COLLIDE_GET

#define WALBERLA_LBM_CELLWISE_SWEEP_D3Q19_COLLIDE_GET ( )
Value:
const real_t vC = src->get( x, y, z, Stencil_T::idx[C] ); \
const real_t vN = src->get( x, y, z, Stencil_T::idx[N] ); \
const real_t vS = src->get( x, y, z, Stencil_T::idx[S] ); \
const real_t vW = src->get( x, y, z, Stencil_T::idx[W] ); \
const real_t vE = src->get( x, y, z, Stencil_T::idx[E] ); \
const real_t vT = src->get( x, y, z, Stencil_T::idx[T] ); \
const real_t vB = src->get( x, y, z, Stencil_T::idx[B] ); \
const real_t vNW = src->get( x, y, z, Stencil_T::idx[NW] ); \
const real_t vNE = src->get( x, y, z, Stencil_T::idx[NE] ); \
const real_t vSW = src->get( x, y, z, Stencil_T::idx[SW] ); \
const real_t vSE = src->get( x, y, z, Stencil_T::idx[SE] ); \
const real_t vTN = src->get( x, y, z, Stencil_T::idx[TN] ); \
const real_t vTS = src->get( x, y, z, Stencil_T::idx[TS] ); \
const real_t vTW = src->get( x, y, z, Stencil_T::idx[TW] ); \
const real_t vTE = src->get( x, y, z, Stencil_T::idx[TE] ); \
const real_t vBN = src->get( x, y, z, Stencil_T::idx[BN] ); \
const real_t vBS = src->get( x, y, z, Stencil_T::idx[BS] ); \
const real_t vBW = src->get( x, y, z, Stencil_T::idx[BW] ); \
const real_t vBE = src->get( x, y, z, Stencil_T::idx[BE] );

◆ WALBERLA_LBM_CELLWISE_SWEEP_D3Q19_DENSITY_VELOCITY_COMP

#define WALBERLA_LBM_CELLWISE_SWEEP_D3Q19_DENSITY_VELOCITY_COMP ( )
Value:
const real_t velXTerm = vE + vNE + vSE + vTE + vBE; \
const real_t velYTerm = vN + vNW + vTN + vBN; \
const real_t velZTerm = vT + vTS + vTW; \
\
const real_t rho = vC + vS + vW + vB + vSW + vBS + vBW + velXTerm + velYTerm + velZTerm; \
const real_t invRho = real_t(1.0) / rho; \
\
const real_t velX = invRho * ( velXTerm - vW - vNW - vSW - vTW - vBW ); \
const real_t velY = invRho * ( velYTerm + vNE - vS - vSW - vSE - vTS - vBS ); \
const real_t velZ = invRho * ( velZTerm + vTN + vTE - vB - vBN - vBS - vBW - vBE );

◆ WALBERLA_LBM_CELLWISE_SWEEP_D3Q19_DENSITY_VELOCITY_INCOMP

#define WALBERLA_LBM_CELLWISE_SWEEP_D3Q19_DENSITY_VELOCITY_INCOMP ( )
Value:
const real_t velXTerm = vE + vNE + vSE + vTE + vBE; \
const real_t velYTerm = vN + vNW + vTN + vBN; \
const real_t velZTerm = vT + vTS + vTW; \
\
const real_t rho = vC + vS + vW + vB + vSW + vBS + vBW + velXTerm + velYTerm + velZTerm; \
\
const real_t velX = velXTerm - vW - vNW - vSW - vTW - vBW; \
const real_t velY = velYTerm + vNE - vS - vSW - vSE - vTS - vBS; \
const real_t velZ = velZTerm + vTN + vTE - vB - vBN - vBS - vBW - vBE;

◆ WALBERLA_LBM_CELLWISE_SWEEP_D3Q19_STREAM_COLLIDE_PULL

#define WALBERLA_LBM_CELLWISE_SWEEP_D3Q19_STREAM_COLLIDE_PULL ( )
Value:
const real_t vC = src->get( x , y , z , Stencil_T::idx[C] ); \
const real_t vN = src->get( x , y-1, z , Stencil_T::idx[N] ); \
const real_t vS = src->get( x , y+1, z , Stencil_T::idx[S] ); \
const real_t vW = src->get( x+1, y , z , Stencil_T::idx[W] ); \
const real_t vE = src->get( x-1, y , z , Stencil_T::idx[E] ); \
const real_t vT = src->get( x , y , z-1, Stencil_T::idx[T] ); \
const real_t vB = src->get( x , y , z+1, Stencil_T::idx[B] ); \
const real_t vNW = src->get( x+1, y-1, z , Stencil_T::idx[NW] ); \
const real_t vNE = src->get( x-1, y-1, z , Stencil_T::idx[NE] ); \
const real_t vSW = src->get( x+1, y+1, z , Stencil_T::idx[SW] ); \
const real_t vSE = src->get( x-1, y+1, z , Stencil_T::idx[SE] ); \
const real_t vTN = src->get( x , y-1, z-1, Stencil_T::idx[TN] ); \
const real_t vTS = src->get( x , y+1, z-1, Stencil_T::idx[TS] ); \
const real_t vTW = src->get( x+1, y , z-1, Stencil_T::idx[TW] ); \
const real_t vTE = src->get( x-1, y , z-1, Stencil_T::idx[TE] ); \
const real_t vBN = src->get( x , y-1, z+1, Stencil_T::idx[BN] ); \
const real_t vBS = src->get( x , y+1, z+1, Stencil_T::idx[BS] ); \
const real_t vBW = src->get( x+1, y , z+1, Stencil_T::idx[BW] ); \
const real_t vBE = src->get( x-1, y , z+1, Stencil_T::idx[BE] );

◆ WALBERLA_LBM_CELLWISE_SWEEP_D3Q27_COLLIDE_GET

#define WALBERLA_LBM_CELLWISE_SWEEP_D3Q27_COLLIDE_GET ( )
Value:
const real_t vC = src->get( x, y, z, Stencil_T::idx[C] ); \
const real_t vN = src->get( x, y, z, Stencil_T::idx[N] ); \
const real_t vS = src->get( x, y, z, Stencil_T::idx[S] ); \
const real_t vW = src->get( x, y, z, Stencil_T::idx[W] ); \
const real_t vE = src->get( x, y, z, Stencil_T::idx[E] ); \
const real_t vT = src->get( x, y, z, Stencil_T::idx[T] ); \
const real_t vB = src->get( x, y, z, Stencil_T::idx[B] ); \
const real_t vNW = src->get( x, y, z, Stencil_T::idx[NW] ); \
const real_t vNE = src->get( x, y, z, Stencil_T::idx[NE] ); \
const real_t vSW = src->get( x, y, z, Stencil_T::idx[SW] ); \
const real_t vSE = src->get( x, y, z, Stencil_T::idx[SE] ); \
const real_t vTN = src->get( x, y, z, Stencil_T::idx[TN] ); \
const real_t vTS = src->get( x, y, z, Stencil_T::idx[TS] ); \
const real_t vTW = src->get( x, y, z, Stencil_T::idx[TW] ); \
const real_t vTE = src->get( x, y, z, Stencil_T::idx[TE] ); \
const real_t vBN = src->get( x, y, z, Stencil_T::idx[BN] ); \
const real_t vBS = src->get( x, y, z, Stencil_T::idx[BS] ); \
const real_t vBW = src->get( x, y, z, Stencil_T::idx[BW] ); \
const real_t vBE = src->get( x, y, z, Stencil_T::idx[BE] ); \
const real_t vTNE = src->get( x, y, z, Stencil_T::idx[TNE] ); \
const real_t vTNW = src->get( x, y, z, Stencil_T::idx[TNW] ); \
const real_t vTSE = src->get( x, y, z, Stencil_T::idx[TSE] ); \
const real_t vTSW = src->get( x, y, z, Stencil_T::idx[TSW] ); \
const real_t vBNE = src->get( x, y, z, Stencil_T::idx[BNE] ); \
const real_t vBNW = src->get( x, y, z, Stencil_T::idx[BNW] ); \
const real_t vBSE = src->get( x, y, z, Stencil_T::idx[BSE] ); \
const real_t vBSW = src->get( x, y, z, Stencil_T::idx[BSW] );

◆ WALBERLA_LBM_CELLWISE_SWEEP_D3Q27_DENSITY_VELOCITY_COMP

#define WALBERLA_LBM_CELLWISE_SWEEP_D3Q27_DENSITY_VELOCITY_COMP ( )
Value:
const real_t velXTerm = vE + vNE + vSE + vTE + vBE + vTNE + vTSE + vBNE + vBSE; \
const real_t velYTerm = vN + vNW + vTN + vBN + vTNW + vBNW; \
const real_t velZTerm = vT + vTS + vTW + vTSW; \
\
const real_t rho = vC + vS + vW + vB + vSW + vBS + vBW + vBSW + velXTerm + velYTerm + velZTerm; \
const real_t invRho = real_t(1.0) / rho; \
\
const real_t velX = invRho * ( velXTerm - vW - vNW - vSW - vTW - vBW - vTNW - vTSW - vBNW - vBSW ); \
const real_t velY = invRho * ( velYTerm + vNE + vTNE + vBNE - vS - vSW - vSE - vTS - vBS - vTSE - vTSW - vBSE - vBSW ); \
const real_t velZ = invRho * ( velZTerm + vTN + vTE + vTNE + vTNW + vTSE - vB - vBN - vBS - vBW - vBE - vBNE - vBNW - vBSE - vBSW );

◆ WALBERLA_LBM_CELLWISE_SWEEP_D3Q27_DENSITY_VELOCITY_INCOMP

#define WALBERLA_LBM_CELLWISE_SWEEP_D3Q27_DENSITY_VELOCITY_INCOMP ( )
Value:
const real_t velXTerm = vE + vNE + vSE + vTE + vBE + vTNE + vTSE + vBNE + vBSE; \
const real_t velYTerm = vN + vNW + vTN + vBN + vTNW + vBNW; \
const real_t velZTerm = vT + vTS + vTW + vTSW; \
\
const real_t rho = vC + vS + vW + vB + vSW + vBS + vBW + vBSW + velXTerm + velYTerm + velZTerm; \
\
const real_t velX = velXTerm - vW - vNW - vSW - vTW - vBW - vTNW - vTSW - vBNW - vBSW ; \
const real_t velY = velYTerm + vNE + vTNE + vBNE - vS - vSW - vSE - vTS - vBS - vTSE - vTSW - vBSE - vBSW; \
const real_t velZ = velZTerm + vTN + vTE + vTNE + vTNW + vTSE - vB - vBN - vBS - vBW - vBE - vBNE - vBNW - vBSE - vBSW;

◆ WALBERLA_LBM_CELLWISE_SWEEP_D3Q27_STREAM_COLLIDE_PULL

#define WALBERLA_LBM_CELLWISE_SWEEP_D3Q27_STREAM_COLLIDE_PULL ( )
Value:
const real_t vC = src->get( x , y , z , Stencil_T::idx[C] ); \
const real_t vN = src->get( x , y-1, z , Stencil_T::idx[N] ); \
const real_t vS = src->get( x , y+1, z , Stencil_T::idx[S] ); \
const real_t vW = src->get( x+1, y , z , Stencil_T::idx[W] ); \
const real_t vE = src->get( x-1, y , z , Stencil_T::idx[E] ); \
const real_t vT = src->get( x , y , z-1, Stencil_T::idx[T] ); \
const real_t vB = src->get( x , y , z+1, Stencil_T::idx[B] ); \
const real_t vNW = src->get( x+1, y-1, z , Stencil_T::idx[NW] ); \
const real_t vNE = src->get( x-1, y-1, z , Stencil_T::idx[NE] ); \
const real_t vSW = src->get( x+1, y+1, z , Stencil_T::idx[SW] ); \
const real_t vSE = src->get( x-1, y+1, z , Stencil_T::idx[SE] ); \
const real_t vTN = src->get( x , y-1, z-1, Stencil_T::idx[TN] ); \
const real_t vTS = src->get( x , y+1, z-1, Stencil_T::idx[TS] ); \
const real_t vTW = src->get( x+1, y , z-1, Stencil_T::idx[TW] ); \
const real_t vTE = src->get( x-1, y , z-1, Stencil_T::idx[TE] ); \
const real_t vBN = src->get( x , y-1, z+1, Stencil_T::idx[BN] ); \
const real_t vBS = src->get( x , y+1, z+1, Stencil_T::idx[BS] ); \
const real_t vBW = src->get( x+1, y , z+1, Stencil_T::idx[BW] ); \
const real_t vBE = src->get( x-1, y , z+1, Stencil_T::idx[BE] ); \
const real_t vTNE = src->get( x-1, y-1, z-1, Stencil_T::idx[TNE] ); \
const real_t vTNW = src->get( x+1, y-1, z-1, Stencil_T::idx[TNW] ); \
const real_t vTSE = src->get( x-1, y+1, z-1, Stencil_T::idx[TSE] ); \
const real_t vTSW = src->get( x+1, y+1, z-1, Stencil_T::idx[TSW] ); \
const real_t vBNE = src->get( x-1, y-1, z+1, Stencil_T::idx[BNE] ); \
const real_t vBNW = src->get( x+1, y-1, z+1, Stencil_T::idx[BNW] ); \
const real_t vBSE = src->get( x-1, y+1, z+1, Stencil_T::idx[BSE] ); \
const real_t vBSW = src->get( x+1, y+1, z+1, Stencil_T::idx[BSW] );

◆ WALBERLA_LBM_CELLWISE_SWEEP_STREAM_COLLIDE_FOOT

#define WALBERLA_LBM_CELLWISE_SWEEP_STREAM_COLLIDE_FOOT ( )    src->swapDataPointers( dst ); }

◆ WALBERLA_LBM_CELLWISE_SWEEP_STREAM_COLLIDE_HEAD

#define WALBERLA_LBM_CELLWISE_SWEEP_STREAM_COLLIDE_HEAD ( specialization)
Value:
template< typename LatticeModel_T, typename Filter_T, typename DensityVelocityIn_T, typename DensityVelocityOut_T > \
void CellwiseSweep< LatticeModel_T, Filter_T, DensityVelocityIn_T, DensityVelocityOut_T, typename std::enable_if_t< specialization > \
>::streamCollide( IBlock * const block, const uint_t numberOfGhostLayersToInclude ) \
{ \
PdfField_T * src( nullptr ); \
PdfField_T * dst( nullptr ); \
\
this->getFields( block, src, dst ); \
\
WALBERLA_ASSERT_GREATER( src->nrOfGhostLayers(), numberOfGhostLayersToInclude ); \
WALBERLA_ASSERT_GREATER_EQUAL( dst->nrOfGhostLayers(), numberOfGhostLayersToInclude ); \
\
const auto & lm = src->latticeModel(); \
dst->resetLatticeModel( lm ); /* required so that member functions for getting density and equilibrium velocity can be called for dst! */ \
\
this->filter( *block ); \
this->densityVelocityIn( *block ); \
this->densityVelocityOut( *block );