walberla::lbm_mesapd_coupling::psm Namespace Reference

Classes

class  BoxWithOverlap
 
class  CylindricalBoundaryWithOverlap
 
class  EllipsoidWithOverlap
 
class  HalfSpaceWithOverlap
 
struct  OverlapFractionFunctor
 
class  ParticleAndVolumeFractionMapping
 
class  PSMSweep
 LBM sweep for the partially saturated cells method. More...
 
class  SphereWithOverlap
 

Typedefs

using ParticleAndVolumeFraction_T = std::pair< id_t, real_t >
 
using ParticleAndVolumeFractionField_T = GhostLayerField< std::vector< ParticleAndVolumeFraction_T >, 1 >
 

Functions

template<int Weighting_T>
real_t calculateWeighting (const real_t &, const real_t &)
 
template<>
real_t calculateWeighting< 1 > (const real_t &epsilon, const real_t &)
 
template<>
real_t calculateWeighting< 2 > (const real_t &epsilon, const real_t &tau)
 
template<typename LatticeModel_T , typename Filter_T , typename DensityVelocityIn_T , typename DensityVelocityOut_T , int SolidCollision_T, int Weighting_T, typename ParticleAccessor_T >
shared_ptr< PSMSweep< LatticeModel_T, Filter_T, DensityVelocityIn_T, DensityVelocityOut_T, SolidCollision_T, Weighting_T, ParticleAccessor_T > > makePSMSweep (const BlockDataID &pdfFieldID, const BlockDataID &particleAndVolumeFractionFieldID, const shared_ptr< StructuredBlockStorage > &blockStorage, const shared_ptr< ParticleAccessor_T > &ac, 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 , int SolidCollision_T, int Weighting_T, typename ParticleAccessor_T >
shared_ptr< PSMSweep< LatticeModel_T, Filter_T, DensityVelocityIn_T, DensityVelocityOut_T, SolidCollision_T, Weighting_T, ParticleAccessor_T > > makePSMSweep (const BlockDataID &srcID, const BlockDataID &dstID, const BlockDataID &particleAndVolumeFractionFieldID, const shared_ptr< StructuredBlockStorage > &blockStorage, const shared_ptr< ParticleAccessor_T > &ac, const Filter_T &filter, const DensityVelocityIn_T &densityVelocityIn, const DensityVelocityOut_T &densityVelocityOut)
 
template<typename LatticeModel_T , int SolidCollision_T, int Weighting_T, typename ParticleAccessor_T >
shared_ptr< PSMSweep< LatticeModel_T, walberla::field::DefaultEvaluationFilter, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::DefaultDensityVelocityCallback, SolidCollision_T, Weighting_T, ParticleAccessor_T > > makePSMSweep (const BlockDataID &pdfFieldID, const BlockDataID &particleAndVolumeFractionFieldID, const shared_ptr< StructuredBlockStorage > &blockStorage, const shared_ptr< ParticleAccessor_T > &ac)
 
template<typename LatticeModel_T , int SolidCollision_T, int Weighting_T, typename ParticleAccessor_T >
shared_ptr< PSMSweep< LatticeModel_T, walberla::field::DefaultEvaluationFilter, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::DefaultDensityVelocityCallback, SolidCollision_T, Weighting_T, ParticleAccessor_T > > makePSMSweep (const BlockDataID &srcID, const BlockDataID &dstID, const BlockDataID &particleAndVolumeFractionFieldID, const shared_ptr< StructuredBlockStorage > &blockStorage, const shared_ptr< ParticleAccessor_T > &ac)
 
template<typename LatticeModel_T , typename FlagField_T , int SolidCollision_T, int Weighting_T, typename ParticleAccessor_T >
shared_ptr< PSMSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::DefaultDensityVelocityCallback, SolidCollision_T, Weighting_T, ParticleAccessor_T > > makePSMSweep (const BlockDataID &pdfFieldID, const BlockDataID &particleAndVolumeFractionFieldID, const shared_ptr< StructuredBlockStorage > &blockStorage, const shared_ptr< ParticleAccessor_T > &ac, const ConstBlockDataID &flagFieldID, const Set< FlagUID > &cellsToEvaluate)
 
template<typename LatticeModel_T , typename FlagField_T , int SolidCollision_T, int Weighting_T, typename ParticleAccessor_T >
shared_ptr< PSMSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::DefaultDensityVelocityCallback, SolidCollision_T, Weighting_T, ParticleAccessor_T > > makePSMSweep (const BlockDataID &srcID, const BlockDataID &dstID, const BlockDataID &particleAndVolumeFractionFieldID, const shared_ptr< StructuredBlockStorage > &blockStorage, const shared_ptr< ParticleAccessor_T > &ac, const ConstBlockDataID &flagFieldID, const Set< FlagUID > &cellsToEvaluate)
 
template<typename LatticeModel_T , typename FlagField_T , typename VelocityField_T , int SolidCollision_T, int Weighting_T, typename ParticleAccessor_T >
shared_ptr< PSMSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::VelocityCallback< VelocityField_T >, SolidCollision_T, Weighting_T, ParticleAccessor_T > > makePSMSweep (const BlockDataID &pdfFieldID, const BlockDataID &particleAndVolumeFractionFieldID, const shared_ptr< StructuredBlockStorage > &blockStorage, const shared_ptr< ParticleAccessor_T > &ac, const ConstBlockDataID &flagFieldID, const Set< FlagUID > &cellsToEvaluate, const BlockDataID &velocityFieldID)
 
template<typename LatticeModel_T , typename FlagField_T , typename VelocityField_T , int SolidCollision_T, int Weighting_T, typename ParticleAccessor_T >
shared_ptr< PSMSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::VelocityCallback< VelocityField_T >, SolidCollision_T, Weighting_T, ParticleAccessor_T > > makePSMSweep (const BlockDataID &srcID, const BlockDataID &dstID, const BlockDataID &particleAndVolumeFractionFieldID, const shared_ptr< StructuredBlockStorage > &blockStorage, const shared_ptr< ParticleAccessor_T > &ac, const ConstBlockDataID &flagFieldID, const Set< FlagUID > &cellsToEvaluate, const BlockDataID &velocityFieldID)
 
template<typename LatticeModel_T , typename FlagField_T , typename VelocityField_T , typename DensityField_T , int SolidCollision_T, int Weighting_T, typename ParticleAccessor_T >
shared_ptr< PSMSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::DensityVelocityCallback< VelocityField_T, DensityField_T >, SolidCollision_T, Weighting_T, ParticleAccessor_T > > makePSMSweep (const BlockDataID &pdfFieldID, const BlockDataID &particleAndVolumeFractionFieldID, const shared_ptr< StructuredBlockStorage > &blockStorage, const shared_ptr< ParticleAccessor_T > &ac, 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 , int SolidCollision_T, int Weighting_T, typename ParticleAccessor_T >
shared_ptr< PSMSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::DensityVelocityCallback< VelocityField_T, DensityField_T >, SolidCollision_T, Weighting_T, ParticleAccessor_T > > makePSMSweep (const BlockDataID &srcID, const BlockDataID &dstID, const BlockDataID &particleAndVolumeFractionFieldID, const shared_ptr< StructuredBlockStorage > &blockStorage, const shared_ptr< ParticleAccessor_T > &ac, const ConstBlockDataID &flagFieldID, const Set< FlagUID > &cellsToEvaluate, const BlockDataID &velocityFieldID, const BlockDataID &densityFieldID)
 
template<typename LatticeModel_T , int Weighting_T, typename ParticleAccessor_T >
Vector3< real_tgetPSMMacroscopicVelocity (const IBlock &block, lbm::PdfField< LatticeModel_T > *pdfField, GhostLayerField< std::vector< std::pair< walberla::size_t, real_t > >, 1 > *particleAndVolumeFractionField, StructuredBlockStorage &blockStorage, const Cell &cell, const ParticleAccessor_T &ac)
 Evaluates the macroscopic velocity for a given cell when using the PSM method. More...
 
template<typename LatticeModel_T , int Weighting_T, typename ParticleAccessor_T >
void initializeDomainForPSM (StructuredBlockStorage &blockStorage, const BlockDataID &pdfFieldID, const BlockDataID &particleAndVolumeFractionFieldID, const ParticleAccessor_T &ac)
 Initializes the PDF field inside the bodies according to the velocities of the bodies. More...
 

Typedef Documentation

◆ ParticleAndVolumeFraction_T

◆ ParticleAndVolumeFractionField_T

Function Documentation

◆ calculateWeighting()

template<int Weighting_T>
real_t walberla::lbm_mesapd_coupling::psm::calculateWeighting ( const real_t ,
const real_t  
)

◆ calculateWeighting< 1 >()

template<>
real_t walberla::lbm_mesapd_coupling::psm::calculateWeighting< 1 > ( const real_t epsilon,
const real_t  
)
inline

◆ calculateWeighting< 2 >()

template<>
real_t walberla::lbm_mesapd_coupling::psm::calculateWeighting< 2 > ( const real_t epsilon,
const real_t tau 
)
inline

◆ getPSMMacroscopicVelocity()

template<typename LatticeModel_T , int Weighting_T, typename ParticleAccessor_T >
Vector3< real_t > walberla::lbm_mesapd_coupling::psm::getPSMMacroscopicVelocity ( const IBlock block,
lbm::PdfField< LatticeModel_T > *  pdfField,
GhostLayerField< std::vector< std::pair< walberla::size_t, real_t > >, 1 > *  particleAndVolumeFractionField,
StructuredBlockStorage blockStorage,
const Cell cell,
const ParticleAccessor_T &  ac 
)

Evaluates the macroscopic velocity for a given cell when using the PSM method.

returns the cell local macroscopic velocity for the PSM method cell local velocity = \( (1-\sum_s B_s) * u + \sum_s B_s * v_s \) u = fluid velocity v_s = velocity of object s at cell center

lbm::PdfField< LatticeModel_T > is typically typedef'ed as PdfField_T; GhostLayerField< std::vector< std::pair< walberla::id_t, real_t > >, 1 > is typically typedef'ed as ParticleAndVolumeFractionField_T;

Weighting_T is like in the PSMSweep.

◆ initializeDomainForPSM()

template<typename LatticeModel_T , int Weighting_T, typename ParticleAccessor_T >
void walberla::lbm_mesapd_coupling::psm::initializeDomainForPSM ( StructuredBlockStorage blockStorage,
const BlockDataID pdfFieldID,
const BlockDataID particleAndVolumeFractionFieldID,
const ParticleAccessor_T &  ac 
)

Initializes the PDF field inside the bodies according to the velocities of the bodies.

As the Partially Saturated Cells method relies on executing the LBM sweep also inside the bodies, it is good practice (and for some PSM variants also required) to initialize the PDF field ( i.e. the velocity ) in agreement with possible initial velocities of the bodies. This is also the case in the presence of external forces acting on the fluid, as these will often shift the macroscopic velocities during the initialization of the PDF field.

Note, that the ParticleAndVolumeFractionMapping for PSM has to be called first to have a valid field.

Only the velocity of cells intersecting with bodies is set, pure fluid cells remain unchanged.

◆ makePSMSweep() [1/10]

template<typename LatticeModel_T , int SolidCollision_T, int Weighting_T, typename ParticleAccessor_T >
shared_ptr< PSMSweep< LatticeModel_T, walberla::field::DefaultEvaluationFilter, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::DefaultDensityVelocityCallback, SolidCollision_T, Weighting_T, ParticleAccessor_T > > walberla::lbm_mesapd_coupling::psm::makePSMSweep ( const BlockDataID pdfFieldID,
const BlockDataID particleAndVolumeFractionFieldID,
const shared_ptr< StructuredBlockStorage > &  blockStorage,
const shared_ptr< ParticleAccessor_T > &  ac 
)

◆ makePSMSweep() [2/10]

template<typename LatticeModel_T , typename FlagField_T , int SolidCollision_T, int Weighting_T, typename ParticleAccessor_T >
shared_ptr< PSMSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::DefaultDensityVelocityCallback, SolidCollision_T, Weighting_T, ParticleAccessor_T > > walberla::lbm_mesapd_coupling::psm::makePSMSweep ( const BlockDataID pdfFieldID,
const BlockDataID particleAndVolumeFractionFieldID,
const shared_ptr< StructuredBlockStorage > &  blockStorage,
const shared_ptr< ParticleAccessor_T > &  ac,
const ConstBlockDataID flagFieldID,
const Set< FlagUID > &  cellsToEvaluate 
)

◆ makePSMSweep() [3/10]

template<typename LatticeModel_T , typename FlagField_T , typename VelocityField_T , int SolidCollision_T, int Weighting_T, typename ParticleAccessor_T >
shared_ptr< PSMSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::VelocityCallback< VelocityField_T >, SolidCollision_T, Weighting_T, ParticleAccessor_T > > walberla::lbm_mesapd_coupling::psm::makePSMSweep ( const BlockDataID pdfFieldID,
const BlockDataID particleAndVolumeFractionFieldID,
const shared_ptr< StructuredBlockStorage > &  blockStorage,
const shared_ptr< ParticleAccessor_T > &  ac,
const ConstBlockDataID flagFieldID,
const Set< FlagUID > &  cellsToEvaluate,
const BlockDataID velocityFieldID 
)

◆ makePSMSweep() [4/10]

template<typename LatticeModel_T , typename FlagField_T , typename VelocityField_T , typename DensityField_T , int SolidCollision_T, int Weighting_T, typename ParticleAccessor_T >
shared_ptr< PSMSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::DensityVelocityCallback< VelocityField_T, DensityField_T >, SolidCollision_T, Weighting_T, ParticleAccessor_T > > walberla::lbm_mesapd_coupling::psm::makePSMSweep ( const BlockDataID pdfFieldID,
const BlockDataID particleAndVolumeFractionFieldID,
const shared_ptr< StructuredBlockStorage > &  blockStorage,
const shared_ptr< ParticleAccessor_T > &  ac,
const ConstBlockDataID flagFieldID,
const Set< FlagUID > &  cellsToEvaluate,
const BlockDataID velocityFieldID,
const BlockDataID densityFieldID 
)

◆ makePSMSweep() [5/10]

template<typename LatticeModel_T , typename Filter_T , typename DensityVelocityIn_T , typename DensityVelocityOut_T , int SolidCollision_T, int Weighting_T, typename ParticleAccessor_T >
shared_ptr< PSMSweep< LatticeModel_T, Filter_T, DensityVelocityIn_T, DensityVelocityOut_T, SolidCollision_T, Weighting_T, ParticleAccessor_T > > walberla::lbm_mesapd_coupling::psm::makePSMSweep ( const BlockDataID pdfFieldID,
const BlockDataID particleAndVolumeFractionFieldID,
const shared_ptr< StructuredBlockStorage > &  blockStorage,
const shared_ptr< ParticleAccessor_T > &  ac,
const Filter_T &  filter,
const DensityVelocityIn_T &  densityVelocityIn,
const DensityVelocityOut_T &  densityVelocityOut 
)

◆ makePSMSweep() [6/10]

template<typename LatticeModel_T , int SolidCollision_T, int Weighting_T, typename ParticleAccessor_T >
shared_ptr< PSMSweep< LatticeModel_T, walberla::field::DefaultEvaluationFilter, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::DefaultDensityVelocityCallback, SolidCollision_T, Weighting_T, ParticleAccessor_T > > walberla::lbm_mesapd_coupling::psm::makePSMSweep ( const BlockDataID srcID,
const BlockDataID dstID,
const BlockDataID particleAndVolumeFractionFieldID,
const shared_ptr< StructuredBlockStorage > &  blockStorage,
const shared_ptr< ParticleAccessor_T > &  ac 
)

◆ makePSMSweep() [7/10]

template<typename LatticeModel_T , typename FlagField_T , int SolidCollision_T, int Weighting_T, typename ParticleAccessor_T >
shared_ptr< PSMSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::DefaultDensityVelocityCallback, SolidCollision_T, Weighting_T, ParticleAccessor_T > > walberla::lbm_mesapd_coupling::psm::makePSMSweep ( const BlockDataID srcID,
const BlockDataID dstID,
const BlockDataID particleAndVolumeFractionFieldID,
const shared_ptr< StructuredBlockStorage > &  blockStorage,
const shared_ptr< ParticleAccessor_T > &  ac,
const ConstBlockDataID flagFieldID,
const Set< FlagUID > &  cellsToEvaluate 
)

◆ makePSMSweep() [8/10]

template<typename LatticeModel_T , typename FlagField_T , typename VelocityField_T , int SolidCollision_T, int Weighting_T, typename ParticleAccessor_T >
shared_ptr< PSMSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::VelocityCallback< VelocityField_T >, SolidCollision_T, Weighting_T, ParticleAccessor_T > > walberla::lbm_mesapd_coupling::psm::makePSMSweep ( const BlockDataID srcID,
const BlockDataID dstID,
const BlockDataID particleAndVolumeFractionFieldID,
const shared_ptr< StructuredBlockStorage > &  blockStorage,
const shared_ptr< ParticleAccessor_T > &  ac,
const ConstBlockDataID flagFieldID,
const Set< FlagUID > &  cellsToEvaluate,
const BlockDataID velocityFieldID 
)

◆ makePSMSweep() [9/10]

template<typename LatticeModel_T , typename FlagField_T , typename VelocityField_T , typename DensityField_T , int SolidCollision_T, int Weighting_T, typename ParticleAccessor_T >
shared_ptr< PSMSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::DensityVelocityCallback< VelocityField_T, DensityField_T >, SolidCollision_T, Weighting_T, ParticleAccessor_T > > walberla::lbm_mesapd_coupling::psm::makePSMSweep ( const BlockDataID srcID,
const BlockDataID dstID,
const BlockDataID particleAndVolumeFractionFieldID,
const shared_ptr< StructuredBlockStorage > &  blockStorage,
const shared_ptr< ParticleAccessor_T > &  ac,
const ConstBlockDataID flagFieldID,
const Set< FlagUID > &  cellsToEvaluate,
const BlockDataID velocityFieldID,
const BlockDataID densityFieldID 
)

◆ makePSMSweep() [10/10]

template<typename LatticeModel_T , typename Filter_T , typename DensityVelocityIn_T , typename DensityVelocityOut_T , int SolidCollision_T, int Weighting_T, typename ParticleAccessor_T >
shared_ptr< PSMSweep< LatticeModel_T, Filter_T, DensityVelocityIn_T, DensityVelocityOut_T, SolidCollision_T, Weighting_T, ParticleAccessor_T > > walberla::lbm_mesapd_coupling::psm::makePSMSweep ( const BlockDataID srcID,
const BlockDataID dstID,
const BlockDataID particleAndVolumeFractionFieldID,
const shared_ptr< StructuredBlockStorage > &  blockStorage,
const shared_ptr< ParticleAccessor_T > &  ac,
const Filter_T &  filter,
const DensityVelocityIn_T &  densityVelocityIn,
const DensityVelocityOut_T &  densityVelocityOut 
)