template<typename LatticeModel_T, typename FlagField_T, typename ParticleAccessor_T>
class walberla::lbm_mesapd_coupling::CurvedLinear< LatticeModel_T, FlagField_T, ParticleAccessor_T >
Linear boundary handling for moving particles.
This boundary condition implements the CLI scheme from Ginzburg et al. - "Two-Relaxation-Time Lattice Boltzmann Scheme: About Parametrization, Velocity, Pressure and Mixed Boundary Conditions" (2008) References to equations and tables in the code documentation are contained in this paper.
It uses one additional cell in inverse direction of the boundary handling for a linear interpolation.
In case not enough cells are available, one fall back solution is available: If no additional cell is available, the scheme is replaced by the bounce back scheme (see SimpleBB.h).
|
| CurvedLinear (const BoundaryUID &boundaryUID, const FlagUID &uid, PDFField_T *const pdfField, const FlagField_T *const flagField, ParticleField_T *const particleField, const shared_ptr< ParticleAccessor_T > &ac, const flag_t domain, const StructuredBlockStorage &blockStorage, const IBlock &block, std::function< real_t(const Vector3< real_t > &)> hydrostaticDensityFct=nullptr) |
|
void | pushFlags (std::vector< FlagUID > &uids) const |
|
void | beforeBoundaryTreatment () const |
|
void | afterBoundaryTreatment () const |
|
template<typename Buffer_T > |
void | packCell (Buffer_T &, const cell_idx_t, const cell_idx_t, const cell_idx_t) const |
|
template<typename Buffer_T > |
void | registerCell (Buffer_T &, const flag_t, const cell_idx_t, const cell_idx_t, const cell_idx_t) |
|
void | registerCell (const flag_t, const cell_idx_t, const cell_idx_t, const cell_idx_t, const BoundaryConfiguration &) |
|
void | registerCells (const flag_t, const CellInterval &, const BoundaryConfiguration &) |
|
template<typename CellIterator > |
void | registerCells (const flag_t, const CellIterator &, const CellIterator &, const BoundaryConfiguration &) |
|
void | unregisterCell (const flag_t, const cell_idx_t, const cell_idx_t, const cell_idx_t) const |
|
void | treatDirection (const cell_idx_t x, const cell_idx_t y, const cell_idx_t z, const stencil::Direction dir, const cell_idx_t nx, const cell_idx_t ny, const cell_idx_t nz, const flag_t mask) |
|
| Boundary (const BoundaryUID &boundaryUID) |
|
void | setMask (const flag_t mask) |
|
flag_t | getMask () const |
|
const BoundaryUID & | getUID () const |
|