OpenLB 1.8.1
Loading...
Searching...
No Matches
olb::RTLBMdynamicsMcHardyRK< T, DESCRIPTOR, MOMENTA > Struct Template Reference

#include <advectionDiffusionDynamics.h>

+ Inheritance diagram for olb::RTLBMdynamicsMcHardyRK< T, DESCRIPTOR, MOMENTA >:
+ Collaboration diagram for olb::RTLBMdynamicsMcHardyRK< T, DESCRIPTOR, MOMENTA >:

Public Types

using MomentaF = typename MOMENTA::template type<DESCRIPTOR>
 
using EquilibriumF = typename equilibria::ZerothOrder::template type<DESCRIPTOR,MOMENTA>
 
using parameters = meta::list<Light::ABSORPTION,Light::SCATTERING,Light::ANISOMATRIX>
 
template<typename NEW_T >
using exchange_value_type = RTLBMdynamicsMcHardyRK<NEW_T,DESCRIPTOR,MOMENTA>
 
template<typename M >
using exchange_momenta = RTLBMdynamicsMcHardyRK<T,DESCRIPTOR,M>
 
- Public Types inherited from olb::dynamics::CustomCollision< T, DESCRIPTOR, momenta::BulkTuple >
using value_t
 
using descriptor_t
 
using MomentaF
 
- Public Types inherited from olb::Dynamics< T, DESCRIPTOR >
using value_t = T
 
using descriptor_t = DESCRIPTOR
 

Public Member Functions

std::type_index id () override
 Expose unique type-identifier for RTTI.
 
AbstractParameters< T, DESCRIPTOR > & getParameters (BlockLattice< T, DESCRIPTOR > &block) override
 Parameters access for legacy post processors.
 
template<typename CELL , typename PARAMETERS , typename V = typename CELL::value_t>
CellStatistic< V > collide (CELL &cell, PARAMETERS &parameters) any_platform
 
template<typename CELL , typename MATRIX , typename V = typename CELL::value_t>
void computeEquilibriumAniso (CELL &cell, Vector< V, DESCRIPTOR::q > &feq, MATRIX &anisoMatrix) const
 
template<typename CELL , typename VECTOR , typename V = typename CELL::value_t>
VECTOR doCollision (CELL &cell, VECTOR &feq, V absorption, V scattering) const
 
void computeEquilibrium (ConstCell< T, DESCRIPTOR > &cell, T irradiance, const T u[DESCRIPTOR::d], T fEq[DESCRIPTOR::q]) const override
 Return iPop equilibrium for given first and second momenta.
 
std::string getName () const override
 Return human-readable name.
 
 RTLBMdynamicsMcHardyRK (T latticeAbsorption, T latticeScattering, std::array< std::array< T, DESCRIPTOR::q >, DESCRIPTOR::q > &anisoMatrix)
 Constructor.
 
computeEquilibrium (int iPop, T rho, const T u[DESCRIPTOR::d], T uSqr) const override
 Compute equilibrium distribution function.
 
CellStatistic< T > collide (Cell< T, DESCRIPTOR > &cell, LatticeStatistics< T > &statistics) override
 Collision step.
 
getOmega () const
 Get local relaxation parameter of the dynamics.
 
void setOmega (T omega)
 Set local relaxation parameter of the dynamics.
 
template<typename M >
using exchange_momenta
 
- Public Member Functions inherited from olb::dynamics::CustomCollision< T, DESCRIPTOR, momenta::BulkTuple >
void initialize (Cell< T, DESCRIPTOR > &cell) override
 Initialize dynamics-specific data for cell.
 
computeRho (ConstCell< T, DESCRIPTOR > &cell) const override
 Compute particle density.
 
void computeU (ConstCell< T, DESCRIPTOR > &cell, T u[DESCRIPTOR::d]) const override
 Compute fluid velocity.
 
void computeJ (ConstCell< T, DESCRIPTOR > &cell, T j[DESCRIPTOR::d]) const override
 Compute fluid momentum.
 
void computeStress (ConstCell< T, DESCRIPTOR > &cell, T rho, const T u[DESCRIPTOR::d], T pi[util::TensorVal< DESCRIPTOR >::n]) const override
 Compute stress tensor.
 
void computeRhoU (ConstCell< T, DESCRIPTOR > &cell, T &rho, T u[DESCRIPTOR::d]) const override
 Compute fluid velocity and particle density.
 
void computeAllMomenta (ConstCell< T, DESCRIPTOR > &cell, T &rho, T u[DESCRIPTOR::d], T pi[util::TensorVal< DESCRIPTOR >::n]) const override
 Compute all momenta up to second order.
 
void defineRho (Cell< T, DESCRIPTOR > &cell, T rho) override
 Set particle density.
 
void defineU (Cell< T, DESCRIPTOR > &cell, const T u[DESCRIPTOR::d]) override
 Set fluid velocity.
 
void defineRhoU (Cell< T, DESCRIPTOR > &cell, T rho, const T u[DESCRIPTOR::d]) override
 Define fluid velocity and particle density.
 
void defineAllMomenta (Cell< T, DESCRIPTOR > &cell, T rho, const T u[DESCRIPTOR::d], const T pi[util::TensorVal< DESCRIPTOR >::n]) override
 Define all momenta up to second order.
 
void inverseShiftRhoU (ConstCell< T, DESCRIPTOR > &cell, T &rho, T u[DESCRIPTOR::d]) const override
 Calculate population momenta s.t. the physical momenta are reproduced by the computeRhoU.
 
- Public Member Functions inherited from olb::Dynamics< T, DESCRIPTOR >
virtual ~Dynamics () any_platform
 
virtual CellStatistic< T > collide (Cell< T, DESCRIPTOR > &cell)
 Perform purely-local collision step on Cell interface (legacy, to be deprecated)
 
void iniEquilibrium (Cell< T, DESCRIPTOR > &cell, T rho, const T u[DESCRIPTOR::d])
 Initialize to equilibrium distribution.
 
void iniRegularized (Cell< T, DESCRIPTOR > &cell, T rho, const T u[DESCRIPTOR::d], const T pi[util::TensorVal< DESCRIPTOR >::n])
 Initialize cell to equilibrium and non-equilibrum part.
 
- Public Member Functions inherited from olb::legacy::BasicDynamics< T, DESCRIPTOR, momenta::BulkTuple >

Detailed Description

template<typename T, typename DESCRIPTOR, typename MOMENTA = momenta::BulkTuple>
struct olb::RTLBMdynamicsMcHardyRK< T, DESCRIPTOR, MOMENTA >

Definition at line 71 of file rtlbmDynamics.h.

Member Typedef Documentation

◆ EquilibriumF

template<typename T , typename DESCRIPTOR , typename MOMENTA = momenta::BulkTuple>
using olb::RTLBMdynamicsMcHardyRK< T, DESCRIPTOR, MOMENTA >::EquilibriumF = typename equilibria::ZerothOrder::template type<DESCRIPTOR,MOMENTA>

Definition at line 843 of file advectionDiffusionDynamics.h.

◆ exchange_momenta

template<typename T , typename DESCRIPTOR , typename MOMENTA = momenta::BulkTuple>
template<typename M >
using olb::RTLBMdynamicsMcHardyRK< T, DESCRIPTOR, MOMENTA >::exchange_momenta = RTLBMdynamicsMcHardyRK<T,DESCRIPTOR,M>

Definition at line 850 of file advectionDiffusionDynamics.h.

◆ exchange_value_type

template<typename T , typename DESCRIPTOR , typename MOMENTA = momenta::BulkTuple>
template<typename NEW_T >
using olb::RTLBMdynamicsMcHardyRK< T, DESCRIPTOR, MOMENTA >::exchange_value_type = RTLBMdynamicsMcHardyRK<NEW_T,DESCRIPTOR,MOMENTA>

Definition at line 847 of file advectionDiffusionDynamics.h.

◆ MomentaF

template<typename T , typename DESCRIPTOR , typename MOMENTA = momenta::BulkTuple>
using olb::RTLBMdynamicsMcHardyRK< T, DESCRIPTOR, MOMENTA >::MomentaF = typename MOMENTA::template type<DESCRIPTOR>

Definition at line 842 of file advectionDiffusionDynamics.h.

◆ parameters

template<typename T , typename DESCRIPTOR , typename MOMENTA = momenta::BulkTuple>
using olb::RTLBMdynamicsMcHardyRK< T, DESCRIPTOR, MOMENTA >::parameters = meta::list<Light::ABSORPTION,Light::SCATTERING,Light::ANISOMATRIX>

Definition at line 844 of file advectionDiffusionDynamics.h.

Constructor & Destructor Documentation

◆ RTLBMdynamicsMcHardyRK()

template<typename T , typename DESCRIPTOR , typename MOMENTA >
olb::RTLBMdynamicsMcHardyRK< T, DESCRIPTOR, MOMENTA >::RTLBMdynamicsMcHardyRK ( T latticeAbsorption,
T latticeScattering,
std::array< std::array< T, DESCRIPTOR::q >, DESCRIPTOR::q > & anisoMatrix )

Constructor.

Definition at line 92 of file rtlbmDynamics.hh.

94 : legacy::BasicDynamics<T, DESCRIPTOR, MOMENTA>(), _absorption(latticeAbsorption), _scattering(latticeScattering), _anisoMatrix(anisoMatrix)
95{
96 this->getName() = "RTLBMdynamicsMcHardyRK";
97}
std::string getName() const override
Return human-readable name.

References olb::RTLBMdynamicsMcHardyRK< T, DESCRIPTOR, MOMENTA >::getName().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

Member Function Documentation

◆ collide() [1/2]

template<typename T , typename DESCRIPTOR , typename MOMENTA = momenta::BulkTuple>
template<typename CELL , typename PARAMETERS , typename V = typename CELL::value_t>
CellStatistic< V > olb::RTLBMdynamicsMcHardyRK< T, DESCRIPTOR, MOMENTA >::collide ( CELL & cell,
PARAMETERS & parameters )
inline

Definition at line 861 of file advectionDiffusionDynamics.h.

861 {
862 const V absorption = parameters.template get<Light::ABSORPTION>();
863 const V scattering = parameters.template get<Light::SCATTERING>();
864 auto anisoMatrix = parameters.template get<Light::ANISOMATRIX>();
865 V irradiance = MomentaF().computeRho(cell);
867 Vector<V, DESCRIPTOR::q> f_pre_collision;
868
869 // calculate f_pre_collision
870 for (int iPop = 0; iPop < DESCRIPTOR::q; ++iPop) {
871 f_pre_collision[iPop] = cell[iPop] + descriptors::t<V,DESCRIPTOR>(iPop);
872 }
873
874 // preparation for RK
875 for (int iPop = 0; iPop < DESCRIPTOR::q; ++iPop) {
876 cell[iPop] += descriptors::t<V,DESCRIPTOR>(iPop);
877 }
878
879 // Runge-Kutta 4
880 computeEquilibriumAniso(cell, feq, anisoMatrix);
881 auto k1 = doCollision(cell, feq, absorption, scattering);
882
883 for (int iPop = 0; iPop < DESCRIPTOR::q; ++iPop) {
884 cell[iPop] = f_pre_collision[iPop] + 0.5 * k1[iPop];
885 }
886
887 computeEquilibriumAniso(cell, feq, anisoMatrix);
888 auto k2 = doCollision(cell, feq, absorption, scattering);
889
890 for (int iPop = 0; iPop < DESCRIPTOR::q; ++iPop) {
891 cell[iPop] = f_pre_collision[iPop] + 0.5 * k2[iPop];
892 }
893
894 computeEquilibriumAniso(cell, feq, anisoMatrix);
895 auto k3 = doCollision(cell, feq, absorption, scattering);
896
897 for (int iPop = 0; iPop < DESCRIPTOR::q; ++iPop) {
898 cell[iPop] = f_pre_collision[iPop] + k3[iPop];
899 }
900
901 computeEquilibriumAniso(cell, feq, anisoMatrix);
902 auto k4 = doCollision(cell, feq, absorption, scattering);
903
904 for (int iPop = 0; iPop < DESCRIPTOR::q; ++iPop) {
905 cell[iPop] = f_pre_collision[iPop] + (1.0 / 6.0) * (k1[iPop] + 2 * k2[iPop] + 2 * k3[iPop] + k4[iPop])
907 }
908
909 return {irradiance, V(0)};
910 };
constexpr T t(unsigned iPop, tag::CUM) any_platform
Definition cum.h:108
Vector(T &&t, Ts &&... ts) -> Vector< std::remove_cvref_t< T >, 1+sizeof...(Ts)>
meta::list< Light::ABSORPTION, Light::SCATTERING, Light::ANISOMATRIX > parameters
void computeEquilibriumAniso(CELL &cell, Vector< V, DESCRIPTOR::q > &feq, MATRIX &anisoMatrix) const
VECTOR doCollision(CELL &cell, VECTOR &feq, V absorption, V scattering) const
typename MOMENTA::template type< DESCRIPTOR > MomentaF

References olb::RTLBMdynamicsMcHardyRK< T, DESCRIPTOR, MOMENTA >::computeEquilibriumAniso(), olb::RTLBMdynamicsMcHardyRK< T, DESCRIPTOR, MOMENTA >::doCollision(), and olb::descriptors::t().

+ Here is the call graph for this function:

◆ collide() [2/2]

template<typename T , typename DESCRIPTOR , typename MOMENTA >
CellStatistic< T > olb::RTLBMdynamicsMcHardyRK< T, DESCRIPTOR, MOMENTA >::collide ( Cell< T, DESCRIPTOR > & cell,
LatticeStatistics< T > & statistics )
override

Collision step.

Definition at line 129 of file rtlbmDynamics.hh.

130{
131 std::array<T,DESCRIPTOR::q> feq;
132 std::array<T,DESCRIPTOR::q> f_pre_collision;
133 // separate cell and precollision f_i
134 for ( int iPop = 0; iPop < DESCRIPTOR::q; ++iPop ) {
135 f_pre_collision[iPop] = cell[iPop] + descriptors::t<T,DESCRIPTOR>(iPop);
136 }
137
138 // shift only first collision und equilibrium and then at the very end
139 for ( int iPop = 0; iPop < DESCRIPTOR::q; ++iPop ) {
140 cell[iPop] += descriptors::t<T,DESCRIPTOR>(iPop);
141 }
142 computeEquilibriumAniso(cell,feq);
143 std::array<T,DESCRIPTOR::q> k1 = doCollision(cell,feq);
144 // update cell
145 for (int iPop = 0; iPop < DESCRIPTOR::q; ++iPop) {
146 cell[iPop] = f_pre_collision[iPop] + 0.5*k1[iPop];
147 }
148
149 computeEquilibriumAniso(cell,feq);
150 std::array<T,DESCRIPTOR::q> k2 = doCollision(cell,feq);
151 // update cell
152 for (int iPop = 0; iPop < DESCRIPTOR::q; ++iPop) {
153 cell[iPop] = f_pre_collision[iPop] + 0.5*k2[iPop];
154 }
155
156 computeEquilibriumAniso(cell,feq);
157 std::array<T,DESCRIPTOR::q> k3 = doCollision(cell,feq);
158 // update cell
159 for (int iPop = 0; iPop < DESCRIPTOR::q; ++iPop) {
160 cell[iPop] = f_pre_collision[iPop] + k3[iPop];
161 }
162
163 computeEquilibriumAniso(cell,feq);
164 std::array<T,DESCRIPTOR::q> k4 = doCollision(cell,feq);
165 // update cell
166 for (int iPop = 0; iPop < DESCRIPTOR::q; ++iPop) {
167 cell[iPop] = f_pre_collision[iPop] + 1/6.*(k1[iPop] + 2*k2[iPop] + 2*k3[iPop] + k4[iPop])
168 - descriptors::t<T,DESCRIPTOR>(iPop); // back shift for OpenLB
169 }
170 T temperature = lbm<DESCRIPTOR>::computeRho(cell);
171 statistics.incrementStats( temperature, T() );
172}
static V computeRho(CELL &cell) any_platform
Computation of density.
Definition lbm.h:267

References olb::lbm< DESCRIPTOR >::computeRho(), olb::LatticeStatistics< T >::incrementStats(), and olb::descriptors::t().

+ Here is the call graph for this function:

◆ computeEquilibrium() [1/2]

template<typename T , typename DESCRIPTOR , typename MOMENTA = momenta::BulkTuple>
void olb::RTLBMdynamicsMcHardyRK< T, DESCRIPTOR, MOMENTA >::computeEquilibrium ( ConstCell< T, DESCRIPTOR > & cell,
T rho,
const T u[DESCRIPTOR::d],
T fEq[DESCRIPTOR::q] ) const
inlineoverridevirtual

Return iPop equilibrium for given first and second momenta.

Implements olb::Dynamics< T, DESCRIPTOR >.

Definition at line 938 of file advectionDiffusionDynamics.h.

938 {
939 EquilibriumF().compute(cell, irradiance, u, fEq);
940 };
typename equilibria::ZerothOrder::template type< DESCRIPTOR, MOMENTA > EquilibriumF

◆ computeEquilibrium() [2/2]

template<typename T , typename DESCRIPTOR , typename MOMENTA >
T olb::RTLBMdynamicsMcHardyRK< T, DESCRIPTOR, MOMENTA >::computeEquilibrium ( int iPop,
T rho,
const T u[DESCRIPTOR::d],
T uSqr ) const
override

Compute equilibrium distribution function.

Definition at line 100 of file rtlbmDynamics.hh.

101{
103}

References olb::descriptors::t().

+ Here is the call graph for this function:

◆ computeEquilibriumAniso()

template<typename T , typename DESCRIPTOR , typename MOMENTA = momenta::BulkTuple>
template<typename CELL , typename MATRIX , typename V = typename CELL::value_t>
void olb::RTLBMdynamicsMcHardyRK< T, DESCRIPTOR, MOMENTA >::computeEquilibriumAniso ( CELL & cell,
Vector< V, DESCRIPTOR::q > & feq,
MATRIX & anisoMatrix ) const
inline

Definition at line 913 of file advectionDiffusionDynamics.h.

913 {
914// feq.fill(V());
915 for (int i = 0; i < DESCRIPTOR::q; ++i) {
916 feq[i] = V();
917 }
918 for (int iPop = 0; iPop < DESCRIPTOR::q; ++iPop) {
919 for (int jPop = 0; jPop < DESCRIPTOR::q; ++jPop) {
920 feq[iPop] += cell[jPop] * anisoMatrix[jPop+iPop];
921 }
922 feq[iPop] *= descriptors::t<V, DESCRIPTOR>(iPop);
923 }
924 };

References olb::descriptors::t().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ doCollision()

template<typename T , typename DESCRIPTOR , typename MOMENTA = momenta::BulkTuple>
template<typename CELL , typename VECTOR , typename V = typename CELL::value_t>
VECTOR olb::RTLBMdynamicsMcHardyRK< T, DESCRIPTOR, MOMENTA >::doCollision ( CELL & cell,
VECTOR & feq,
V absorption,
V scattering ) const
inline

Definition at line 927 of file advectionDiffusionDynamics.h.

927 {
928 VECTOR k;
929 for (int iPop = 0; iPop < DESCRIPTOR::q; ++iPop) {
931 V fFull = cell[iPop];
932 k[iPop] = -normC * (absorption + scattering) * fFull
933 + normC * scattering * feq[iPop];
934 }
935 return k;
936 };
constexpr int c(unsigned iPop, unsigned iDim) any_platform
Definition functions.h:83
T norm(const std::vector< T > &a)
l2 norm of a vector of arbitrary length

References olb::descriptors::c(), and olb::util::norm().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ exchange_momenta()

template<typename T , typename DESCRIPTOR , typename MOMENTA = momenta::BulkTuple>
template<typename M >
typedef RTLBMdynamicsMcHardyRK<T,DESCRIPTOR,M> olb::RTLBMdynamicsMcHardyRK< T, DESCRIPTOR, MOMENTA >::exchange_momenta

Definition at line 76 of file rtlbmDynamics.h.

◆ getName()

template<typename T , typename DESCRIPTOR , typename MOMENTA = momenta::BulkTuple>
std::string olb::RTLBMdynamicsMcHardyRK< T, DESCRIPTOR, MOMENTA >::getName ( ) const
inlineoverridevirtual

Return human-readable name.

Reimplemented from olb::Dynamics< T, DESCRIPTOR >.

Definition at line 942 of file advectionDiffusionDynamics.h.

942 {
943 return "RTLBMdynamicsMcHardyRK<" + MomentaF().getName() + ">";
944 };
+ Here is the caller graph for this function:

◆ getOmega()

template<typename T , typename DESCRIPTOR , typename MOMENTA >
T olb::RTLBMdynamicsMcHardyRK< T, DESCRIPTOR, MOMENTA >::getOmega ( ) const

Get local relaxation parameter of the dynamics.

Definition at line 175 of file rtlbmDynamics.hh.

176{
177 return -1;
178}

◆ getParameters()

template<typename T , typename DESCRIPTOR , typename MOMENTA = momenta::BulkTuple>
AbstractParameters< T, DESCRIPTOR > & olb::RTLBMdynamicsMcHardyRK< T, DESCRIPTOR, MOMENTA >::getParameters ( BlockLattice< T, DESCRIPTOR > & block)
inlineoverridevirtual

Parameters access for legacy post processors.

Implements olb::Dynamics< T, DESCRIPTOR >.

Definition at line 856 of file advectionDiffusionDynamics.h.

856 {
857 return block.template getData<OperatorParameters<RTLBMdynamicsMcHardyRK>>();
858 }

◆ id()

template<typename T , typename DESCRIPTOR , typename MOMENTA = momenta::BulkTuple>
std::type_index olb::RTLBMdynamicsMcHardyRK< T, DESCRIPTOR, MOMENTA >::id ( )
inlineoverridevirtual

Expose unique type-identifier for RTTI.

Implements olb::Dynamics< T, DESCRIPTOR >.

Definition at line 852 of file advectionDiffusionDynamics.h.

852 {
853 return typeid(RTLBMdynamicsMcHardyRK);
854 };
RTLBMdynamicsMcHardyRK(T latticeAbsorption, T latticeScattering, std::array< std::array< T, DESCRIPTOR::q >, DESCRIPTOR::q > &anisoMatrix)
Constructor.

References olb::RTLBMdynamicsMcHardyRK< T, DESCRIPTOR, MOMENTA >::RTLBMdynamicsMcHardyRK().

+ Here is the call graph for this function:

◆ setOmega()

template<typename T , typename DESCRIPTOR , typename MOMENTA >
void olb::RTLBMdynamicsMcHardyRK< T, DESCRIPTOR, MOMENTA >::setOmega ( T omega)

Set local relaxation parameter of the dynamics.

Definition at line 181 of file rtlbmDynamics.hh.

182{
183}

The documentation for this struct was generated from the following files: