walberla::pe_coupling::discrete_particle_methods Namespace Reference

Classes

class  AddedMassForceEvaluator
 Evaluator of the added (virtual) mass force based on the given added mass correlation. More...
 
class  AveragedInteractionForceFieldToForceFieldAdder
 Adds an accumulated force field value to another force field, applying continuous averaging. More...
 
class  BodyVelocityInitializer
 Initializes the bodies with the velocity of the surrounding fluid. More...
 
class  BodyVelocityTimeDerivativeEvaluator
 Evaluates the time derivative of the body velocity. More...
 
class  CombinedReductionFieldCommunication
 Special communication routine to += the values in the ghost layers and update the ghost layers. More...
 
class  EffectiveViscosityFieldEvaluator
 Evaluates the (local) effective viscosity. More...
 
class  ExternalForceToForceFieldAdder
 Applies a given constant external force to a force field. More...
 
class  FieldDataSwapper
 Swaps the data of two identical fields. More...
 
class  ForceFieldResetter
 Resets the values currently stored in the force field to (0,0,0). More...
 
class  GNSExternalForceToForceFieldAdder
 Applies a given constant external force to a force field. More...
 
class  GNSPressureFieldEvaluator
 Evaluates the pressure field based on the given (GNS) LBM field. More...
 
class  GNSSmagorinskyLESField
 Adjusts locally the fluid viscosity based on a LES-Smagorinsky turbulence model when using GNS-LBM. More...
 
class  GNSSweep
 Sweep for the LBM formulation of the generalized Navier Stokes equations. More...
 
class  GNSVelocityFieldEvaluator
 Evaluator for the fluid-phase velocity, given a (GNS) PDf field. More...
 
class  InteractionForceEvaluator
 Evaluator of the two most important interaction forces: drag and pressure gradient. More...
 
class  LiftForceEvaluator
 Evaluator of the lift force based on the given lift correlation. More...
 
class  LubricationForceEvaluator
 Evaluates the lubrication force between for a sphere-sphere or a sphere-plane pair. More...
 
class  PressureFieldEvaluator
 Evaluator of the pressure field, given a PDF field. More...
 
class  PressureGradientFieldEvaluator
 Evaluator of the pressure gradient field, given a pressure field. More...
 
class  SolidVolumeFractionFieldEvaluator
 Evaluator of the solid volume fraction field. More...
 
class  StressTensorGradientFieldEvaluator
 Evaluator of the stress tensor gradient field, given a velocity field. More...
 
class  VelocityCurlFieldEvaluator
 Evaluator of the curl of a given velocity field. More...
 
class  VelocityFieldEvaluator
 Evaluator of the velocity field, given a PDF field. More...
 
class  VelocityGradientFieldEvaluator
 Evaluator of the velocity gradient field, given a velocity field. More...
 
class  VelocityTotalTimeDerivativeFieldEvaluator
 Evaluator of the total derivative of the fluid velocity. More...
 

Functions

Vector3< real_taddedMassForceFinn (const Vector3< real_t > &timeDerivativeFluidVel, const Vector3< real_t > &timeDerivativeBodyVel, const real_t &bodyVolume, const real_t &fluidDensity)
 Correlation functions for the added mass force. More...
 
Vector3< real_tnoAddedMassForce (const Vector3< real_t > &, const Vector3< real_t > &, const real_t &, const real_t &)
 
real_t dragCoeffSchillerNaumann (real_t reynoldsNumber)
 Various correlation functions for the drag force. More...
 
real_t dragCoeffStokes (real_t fluidVolumeFraction, real_t diameter, real_t fluidDynamicViscosity)
 
Vector3< real_tdragForceStokes (const Vector3< real_t > &fluidVel, const Vector3< real_t > &particleVel, real_t solidVolumeFraction, real_t diameter, real_t fluidDynamicViscosity, real_t)
 
Vector3< real_tdragForceErgunWenYu (const Vector3< real_t > &fluidVel, const Vector3< real_t > &particleVel, real_t solidVolumeFraction, real_t diameter, real_t fluidDynamicViscosity, real_t fluidDensity)
 
Vector3< real_tdragForceTang (const Vector3< real_t > &fluidVel, const Vector3< real_t > &particleVel, real_t solidVolumeFraction, real_t diameter, real_t fluidDynamicViscosity, real_t fluidDensity)
 
Vector3< real_tdragForceFelice (const Vector3< real_t > &fluidVel, const Vector3< real_t > &particleVel, real_t solidVolumeFraction, real_t diameter, real_t fluidDynamicViscosity, real_t fluidDensity)
 
Vector3< real_tdragForceTenneti (const Vector3< real_t > &fluidVel, const Vector3< real_t > &particleVel, real_t solidVolumeFraction, real_t diameter, real_t fluidDynamicViscosity, real_t fluidDensity)
 
Vector3< real_tnoDragForce (const Vector3< real_t > &, const Vector3< real_t > &, real_t, real_t, real_t, real_t)
 
Vector3< real_tliftForceSaffman (const Vector3< real_t > &fluidVel, const Vector3< real_t > &curlFluidVel, const Vector3< real_t > &particleVel, real_t diameter, real_t fluidDynamicViscosity, real_t fluidDensity)
 Correlation functions for the lift force. More...
 
Vector3< real_tnoLiftForce (const Vector3< real_t > &, const Vector3< real_t > &, const Vector3< real_t > &, real_t, real_t, real_t)
 
real_t calculateUnchangedEffectiveViscosity (real_t fluidViscosity, real_t)
 
real_t calculateRescaledEffectiveViscosity (real_t fluidViscosity, real_t porosity)
 
real_t calculateEilersEffectiveViscosity (real_t fluidViscosity, real_t porosity)
 
template<typename LatticeModel_T , typename Filter_T , typename DensityVelocityIn_T , typename DensityVelocityOut_T >
shared_ptr< GNSSweep< LatticeModel_T, Filter_T, DensityVelocityIn_T, DensityVelocityOut_T > > makeGNSSweep (const BlockDataID &pdfFieldID, const BlockDataID &solidVolumeFractionFieldID, 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< GNSSweep< LatticeModel_T, Filter_T, DensityVelocityIn_T, DensityVelocityOut_T > > makeGNSSweep (const BlockDataID &srcID, const BlockDataID &dstID, const BlockDataID &solidVolumeFractionFieldID, const Filter_T &filter, const DensityVelocityIn_T &densityVelocityIn, const DensityVelocityOut_T &densityVelocityOut)
 
template<typename LatticeModel_T >
shared_ptr< GNSSweep< LatticeModel_T, walberla::field::DefaultEvaluationFilter, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::DefaultDensityVelocityCallback > > makeGNSSweep (const BlockDataID &pdfFieldID, const BlockDataID &solidVolumeFractionFieldID)
 
template<typename LatticeModel_T >
shared_ptr< GNSSweep< LatticeModel_T, walberla::field::DefaultEvaluationFilter, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::DefaultDensityVelocityCallback > > makeGNSSweep (const BlockDataID &srcID, const BlockDataID &dstID, const BlockDataID &solidVolumeFractionFieldID)
 
template<typename LatticeModel_T , typename FlagField_T >
shared_ptr< GNSSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::DefaultDensityVelocityCallback > > makeGNSSweep (const BlockDataID &pdfFieldID, const BlockDataID &solidVolumeFractionFieldID, const ConstBlockDataID &flagFieldID, const Set< FlagUID > &cellsToEvaluate)
 
template<typename LatticeModel_T , typename FlagField_T >
shared_ptr< GNSSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::DefaultDensityVelocityCallback > > makeGNSSweep (const BlockDataID &srcID, const BlockDataID &dstID, const BlockDataID &solidVolumeFractionFieldID, const ConstBlockDataID &flagFieldID, const Set< FlagUID > &cellsToEvaluate)
 
template<typename LatticeModel_T , typename FlagField_T , typename VelocityField_T >
shared_ptr< GNSSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::VelocityCallback< VelocityField_T > > > makeGNSSweep (const BlockDataID &pdfFieldID, const BlockDataID &solidVolumeFractionFieldID, const ConstBlockDataID &flagFieldID, const Set< FlagUID > &cellsToEvaluate, const BlockDataID &velocityFieldID)
 
template<typename LatticeModel_T , typename FlagField_T , typename VelocityField_T >
shared_ptr< GNSSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::VelocityCallback< VelocityField_T > > > makeGNSSweep (const BlockDataID &srcID, const BlockDataID &dstID, const BlockDataID &solidVolumeFractionFieldID, 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< GNSSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::DensityVelocityCallback< VelocityField_T, DensityField_T > > > makeGNSSweep (const BlockDataID &pdfFieldID, const BlockDataID &solidVolumeFractionFieldID, 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< GNSSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter< FlagField_T >, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::DensityVelocityCallback< VelocityField_T, DensityField_T > > > makeGNSSweep (const BlockDataID &srcID, const BlockDataID &dstID, const BlockDataID &solidVolumeFractionFieldID, const ConstBlockDataID &flagFieldID, const Set< FlagUID > &cellsToEvaluate, const BlockDataID &velocityFieldID, const BlockDataID &densityFieldID)
 

Variables

const real_t thresholdAbsoluteVelocityDifference = real_t(1e-10)
 

Function Documentation

◆ addedMassForceFinn()

Vector3<real_t> walberla::pe_coupling::discrete_particle_methods::addedMassForceFinn ( const Vector3< real_t > &  timeDerivativeFluidVel,
const Vector3< real_t > &  timeDerivativeBodyVel,
const real_t bodyVolume,
const real_t fluidDensity 
)

Correlation functions for the added mass force.

To be compatible with the interface of AddedMassForceEvaluator, all functions use the following signature: Vector3<real_t> ( const Vector3<real_t> & timeDerivativeFluidVel, const Vector3<real_t> & timeDerivativeBodyVel, const real_t & bodyVolume, const real_t & fluidDensity )

◆ calculateEilersEffectiveViscosity()

real_t walberla::pe_coupling::discrete_particle_methods::calculateEilersEffectiveViscosity ( real_t  fluidViscosity,
real_t  porosity 
)

◆ calculateRescaledEffectiveViscosity()

real_t walberla::pe_coupling::discrete_particle_methods::calculateRescaledEffectiveViscosity ( real_t  fluidViscosity,
real_t  porosity 
)

◆ calculateUnchangedEffectiveViscosity()

real_t walberla::pe_coupling::discrete_particle_methods::calculateUnchangedEffectiveViscosity ( real_t  fluidViscosity,
real_t   
)

◆ dragCoeffSchillerNaumann()

real_t walberla::pe_coupling::discrete_particle_methods::dragCoeffSchillerNaumann ( real_t  reynoldsNumber)

Various correlation functions for the drag force.

These functions calculate the drag force for fluid-particle interactions based on different empirical correlations from literature. Always be aware that those empirical formulas were obtained by different setups and are thus not generally applicable!

To be compatible with the interface of InteractionForceEvaluator, all functions use the following signature: Vector3<real_t> ( const Vector3<real_t> & fluidVel, const Vector3<real_t> & particleVel, real_t solidVolumeFraction, real_t diameter, real_t fluidDynamicViscosity, real_t fluidDensity )

◆ dragCoeffStokes()

real_t walberla::pe_coupling::discrete_particle_methods::dragCoeffStokes ( real_t  fluidVolumeFraction,
real_t  diameter,
real_t  fluidDynamicViscosity 
)

◆ dragForceErgunWenYu()

Vector3<real_t> walberla::pe_coupling::discrete_particle_methods::dragForceErgunWenYu ( const Vector3< real_t > &  fluidVel,
const Vector3< real_t > &  particleVel,
real_t  solidVolumeFraction,
real_t  diameter,
real_t  fluidDynamicViscosity,
real_t  fluidDensity 
)

◆ dragForceFelice()

Vector3<real_t> walberla::pe_coupling::discrete_particle_methods::dragForceFelice ( const Vector3< real_t > &  fluidVel,
const Vector3< real_t > &  particleVel,
real_t  solidVolumeFraction,
real_t  diameter,
real_t  fluidDynamicViscosity,
real_t  fluidDensity 
)

◆ dragForceStokes()

Vector3<real_t> walberla::pe_coupling::discrete_particle_methods::dragForceStokes ( const Vector3< real_t > &  fluidVel,
const Vector3< real_t > &  particleVel,
real_t  solidVolumeFraction,
real_t  diameter,
real_t  fluidDynamicViscosity,
real_t   
)

◆ dragForceTang()

Vector3<real_t> walberla::pe_coupling::discrete_particle_methods::dragForceTang ( const Vector3< real_t > &  fluidVel,
const Vector3< real_t > &  particleVel,
real_t  solidVolumeFraction,
real_t  diameter,
real_t  fluidDynamicViscosity,
real_t  fluidDensity 
)

◆ dragForceTenneti()

Vector3<real_t> walberla::pe_coupling::discrete_particle_methods::dragForceTenneti ( const Vector3< real_t > &  fluidVel,
const Vector3< real_t > &  particleVel,
real_t  solidVolumeFraction,
real_t  diameter,
real_t  fluidDynamicViscosity,
real_t  fluidDensity 
)

◆ liftForceSaffman()

Vector3<real_t> walberla::pe_coupling::discrete_particle_methods::liftForceSaffman ( const Vector3< real_t > &  fluidVel,
const Vector3< real_t > &  curlFluidVel,
const Vector3< real_t > &  particleVel,
real_t  diameter,
real_t  fluidDynamicViscosity,
real_t  fluidDensity 
)

Correlation functions for the lift force.

To be compatible with the interface of LiftForceEvaluator, all functions use the following signature: Vector3<real_t> ( const Vector3<real_t> & fluidVel, const Vector3<real_t> & curlFluidVel, const Vector3<real_t> & particleVel, real_t diameter, real_t fluidDynamicViscosity, real_t fluidDensity )

◆ makeGNSSweep() [1/10]

template<typename LatticeModel_T >
shared_ptr< GNSSweep< LatticeModel_T, walberla::field::DefaultEvaluationFilter, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::DefaultDensityVelocityCallback > > walberla::pe_coupling::discrete_particle_methods::makeGNSSweep ( const BlockDataID pdfFieldID,
const BlockDataID solidVolumeFractionFieldID 
)

◆ makeGNSSweep() [2/10]

template<typename LatticeModel_T , typename FlagField_T >
shared_ptr< GNSSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter<FlagField_T>, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::DefaultDensityVelocityCallback > > walberla::pe_coupling::discrete_particle_methods::makeGNSSweep ( const BlockDataID pdfFieldID,
const BlockDataID solidVolumeFractionFieldID,
const ConstBlockDataID flagFieldID,
const Set< FlagUID > &  cellsToEvaluate 
)

◆ makeGNSSweep() [3/10]

template<typename LatticeModel_T , typename FlagField_T , typename VelocityField_T >
shared_ptr< GNSSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter<FlagField_T>, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::VelocityCallback<VelocityField_T> > > walberla::pe_coupling::discrete_particle_methods::makeGNSSweep ( const BlockDataID pdfFieldID,
const BlockDataID solidVolumeFractionFieldID,
const ConstBlockDataID flagFieldID,
const Set< FlagUID > &  cellsToEvaluate,
const BlockDataID velocityFieldID 
)

◆ makeGNSSweep() [4/10]

template<typename LatticeModel_T , typename FlagField_T , typename VelocityField_T , typename DensityField_T >
shared_ptr< GNSSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter<FlagField_T>, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::DensityVelocityCallback<VelocityField_T,DensityField_T> > > walberla::pe_coupling::discrete_particle_methods::makeGNSSweep ( const BlockDataID pdfFieldID,
const BlockDataID solidVolumeFractionFieldID,
const ConstBlockDataID flagFieldID,
const Set< FlagUID > &  cellsToEvaluate,
const BlockDataID velocityFieldID,
const BlockDataID densityFieldID 
)

◆ makeGNSSweep() [5/10]

template<typename LatticeModel_T , typename Filter_T , typename DensityVelocityIn_T , typename DensityVelocityOut_T >
shared_ptr< GNSSweep< LatticeModel_T, Filter_T, DensityVelocityIn_T, DensityVelocityOut_T > > walberla::pe_coupling::discrete_particle_methods::makeGNSSweep ( const BlockDataID pdfFieldID,
const BlockDataID solidVolumeFractionFieldID,
const Filter_T &  filter,
const DensityVelocityIn_T &  densityVelocityIn,
const DensityVelocityOut_T &  densityVelocityOut 
)

◆ makeGNSSweep() [6/10]

template<typename LatticeModel_T >
shared_ptr< GNSSweep< LatticeModel_T, walberla::field::DefaultEvaluationFilter, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::DefaultDensityVelocityCallback > > walberla::pe_coupling::discrete_particle_methods::makeGNSSweep ( const BlockDataID srcID,
const BlockDataID dstID,
const BlockDataID solidVolumeFractionFieldID 
)

◆ makeGNSSweep() [7/10]

template<typename LatticeModel_T , typename FlagField_T >
shared_ptr< GNSSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter<FlagField_T>, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::DefaultDensityVelocityCallback > > walberla::pe_coupling::discrete_particle_methods::makeGNSSweep ( const BlockDataID srcID,
const BlockDataID dstID,
const BlockDataID solidVolumeFractionFieldID,
const ConstBlockDataID flagFieldID,
const Set< FlagUID > &  cellsToEvaluate 
)

◆ makeGNSSweep() [8/10]

template<typename LatticeModel_T , typename FlagField_T , typename VelocityField_T >
shared_ptr< GNSSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter<FlagField_T>, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::VelocityCallback<VelocityField_T> > > walberla::pe_coupling::discrete_particle_methods::makeGNSSweep ( const BlockDataID srcID,
const BlockDataID dstID,
const BlockDataID solidVolumeFractionFieldID,
const ConstBlockDataID flagFieldID,
const Set< FlagUID > &  cellsToEvaluate,
const BlockDataID velocityFieldID 
)

◆ makeGNSSweep() [9/10]

template<typename LatticeModel_T , typename FlagField_T , typename VelocityField_T , typename DensityField_T >
shared_ptr< GNSSweep< LatticeModel_T, walberla::field::FlagFieldEvaluationFilter<FlagField_T>, lbm::DefaultDensityEquilibriumVelocityCalculation, lbm::DensityVelocityCallback<VelocityField_T,DensityField_T> > > walberla::pe_coupling::discrete_particle_methods::makeGNSSweep ( const BlockDataID srcID,
const BlockDataID dstID,
const BlockDataID solidVolumeFractionFieldID,
const ConstBlockDataID flagFieldID,
const Set< FlagUID > &  cellsToEvaluate,
const BlockDataID velocityFieldID,
const BlockDataID densityFieldID 
)

◆ makeGNSSweep() [10/10]

template<typename LatticeModel_T , typename Filter_T , typename DensityVelocityIn_T , typename DensityVelocityOut_T >
shared_ptr< GNSSweep< LatticeModel_T, Filter_T, DensityVelocityIn_T, DensityVelocityOut_T > > walberla::pe_coupling::discrete_particle_methods::makeGNSSweep ( const BlockDataID srcID,
const BlockDataID dstID,
const BlockDataID solidVolumeFractionFieldID,
const Filter_T &  filter,
const DensityVelocityIn_T &  densityVelocityIn,
const DensityVelocityOut_T &  densityVelocityOut 
)

◆ noAddedMassForce()

Vector3<real_t> walberla::pe_coupling::discrete_particle_methods::noAddedMassForce ( const Vector3< real_t > &  ,
const Vector3< real_t > &  ,
const real_t ,
const real_t  
)

◆ noDragForce()

Vector3<real_t> walberla::pe_coupling::discrete_particle_methods::noDragForce ( const Vector3< real_t > &  ,
const Vector3< real_t > &  ,
real_t  ,
real_t  ,
real_t  ,
real_t   
)

◆ noLiftForce()

Vector3<real_t> walberla::pe_coupling::discrete_particle_methods::noLiftForce ( const Vector3< real_t > &  ,
const Vector3< real_t > &  ,
const Vector3< real_t > &  ,
real_t  ,
real_t  ,
real_t   
)

Variable Documentation

◆ thresholdAbsoluteVelocityDifference

const real_t walberla::pe_coupling::discrete_particle_methods::thresholdAbsoluteVelocityDifference = real_t(1e-10)