MAIA bb96820c
Multiphysics at AIA
Loading...
Searching...
No Matches
maia::lb::LbSrcTerm_smagorinsky< nDim, nDist, SysEqn, bubble > Class Template Reference

Class to handle sub-grid scale model (Smagorinsky) More...

Inheritance diagram for maia::lb::LbSrcTerm_smagorinsky< nDim, nDist, SysEqn, bubble >:
[legend]
Collaboration diagram for maia::lb::LbSrcTerm_smagorinsky< nDim, nDist, SysEqn, bubble >:
[legend]

Public Member Functions

 LbSrcTerm_smagorinsky (LbSolverDxQy< nDim, nDist, SysEqn > *p_solver)
 
void init () override
 
void apply_preCollision () override
 
void apply_postCollision () override
 
void apply_postPropagation () override
 
- Public Member Functions inherited from maia::lb::LbSrcTerm< nDim, nDist, SysEqn >
virtual void init ()=0
 
virtual void apply_preCollision ()=0
 
virtual void apply_postCollision ()=0
 
virtual void apply_postPropagation ()=0
 
virtual ~LbSrcTerm ()=default
 

Protected Member Functions

void readProperties () override
 
virtual void readProperties ()=0
 

Protected Attributes

LbSolverDxQy< nDim, nDist, SysEqn > * m_solver
 
MFloat m_Cs = 0.1
 Smagorinsky constant. More...
 
const MFloat m_deltaX = 1.0
 Filter width. More...
 

Additional Inherited Members

- Public Types inherited from maia::lb::LbSrcTerm< nDim, nDist, SysEqn >
using SysEqn = SysEqn
 
- Static Public Attributes inherited from maia::lb::LbSrcTerm< nDim, nDist, SysEqn >
static constexpr MInt nDim
 
static constexpr MInt nDist
 

Detailed Description

template<MInt nDim, MInt nDist, class SysEqn, MBool bubble = false>
class maia::lb::LbSrcTerm_smagorinsky< nDim, nDist, SysEqn, bubble >
Author
Miro Gondrum
Date
01.02.2022

Definition at line 23 of file lbsrctermsgs.cpp.

Constructor & Destructor Documentation

◆ LbSrcTerm_smagorinsky()

template<MInt nDim, MInt nDist, class SysEqn , MBool bubble = false>
maia::lb::LbSrcTerm_smagorinsky< nDim, nDist, SysEqn, bubble >::LbSrcTerm_smagorinsky ( LbSolverDxQy< nDim, nDist, SysEqn > *  p_solver)
inline

Definition at line 34 of file lbsrctermsgs.cpp.

34: m_solver(p_solver) { readProperties(); };
LbSolverDxQy< nDim, nDist, SysEqn > * m_solver

Member Function Documentation

◆ apply_postCollision()

template<MInt nDim, MInt nDist, class SysEqn , MBool bubble = false>
void maia::lb::LbSrcTerm_smagorinsky< nDim, nDist, SysEqn, bubble >::apply_postCollision ( )
inlineoverridevirtual

Implements maia::lb::LbSrcTerm< nDim, nDist, SysEqn >.

Definition at line 38 of file lbsrctermsgs.cpp.

38{/*NOT USED*/};

◆ apply_postPropagation()

template<MInt nDim, MInt nDist, class SysEqn , MBool bubble = false>
void maia::lb::LbSrcTerm_smagorinsky< nDim, nDist, SysEqn, bubble >::apply_postPropagation ( )
inlineoverridevirtual

Implements maia::lb::LbSrcTerm< nDim, nDist, SysEqn >.

Definition at line 39 of file lbsrctermsgs.cpp.

39{/*NOT USED*/};

◆ apply_preCollision()

template<MInt nDim, MInt nDist, class SysEqn , MBool bubble>
void maia::lb::LbSrcTerm_smagorinsky< nDim, nDist, SysEqn, bubble >::apply_preCollision
overridevirtual

Implements maia::lb::LbSrcTerm< nDim, nDist, SysEqn >.

Definition at line 53 of file lbsrctermsgs.cpp.

53 {
54 TRACE();
55 LbSolverDxQy<nDim, nDist, SysEqn>* const s = m_solver; // alias for readability
56
57 s->m_nu = s->m_Ma * LBCS / s->m_Re * s->m_referenceLength;
58 maia::parallelFor<false>(0, s->m_currentMaxNoCells, [=](MInt index) {
59 const MInt pCellId = s->m_activeCellList[index];
60 const MInt lvlDiff = s->maxLevel() - s->a_level(pCellId);
61 if((globalTimeStep - 1) % IPOW2(lvlDiff) != 0) return;
62 // TODO labels:LB daniell: This is critical for setting a_oldNu ?! But without
63 // test case 3D_LB-FV_Euler-Euler_bubble_sphere fails
64 s->a_nu(pCellId) = s->m_nu;
65 // Ref: Hou, Sterling, Chen, and Doolen 1994
66 // Also worth reading: Siebler, Krafczyk, Freudiger, and Geier 2011
67 constexpr MInt nDimSqr = nDim * nDim;
68 MFloat c[nDimSqr]{};
69 const MFloat rho = s->a_variable(pCellId, s->PV->RHO);
70#ifndef WAR_NVHPC_PSTL
71 const MFloat* const u = &s->a_variable(pCellId, s->PV->U);
72 const MFloat* const dist = &s->a_oldDistribution(pCellId, 0);
73 s->sysEqn().calcMomentumFlux(rho, u, dist, c);
74#else
75 MFloat u[nDim] = {F0};
76 MFloat dist[nDist] = {F0};
77 for(MInt d = 0; d < nDim; d++) {
78 u[d] = s->a_variable(pCellId, d);
79 }
80 for(MInt d = 0; d < nDist; d++) {
81 dist[d] = s->a_oldDistribution(pCellId, d);
82 }
83 s->sysEqn().calcMomentumFlux(rho, u, dist, c);
84#endif
85 const MFloat Q = sqrt(2.0 * std::inner_product(&c[0], &c[nDimSqr], &c[0], .0));
86 MFloat nu_BIT = 0.0;
87 if constexpr(bubble) {
88 if(s->m_isEELiquid) {
89 // BIT (bubble induced turbulence) see Mohammadi 19 / Sato 81 / Milelli 01
90 const MFloat u_d = sqrt(POW2(s->a_uOtherPhase(pCellId, 0) - u[0]) + POW2(s->a_uOtherPhase(pCellId, 1) - u[1])
91 + POW2(s->a_uOtherPhase(pCellId, 2) - u[2]));
92 nu_BIT = m_Cs * s->a_alphaGasLim(pCellId) * u_d;
93 }
94 }
95 // The stress tensor S_mean depends on tau and tau depends on S_mean.
96 // This leads to a quadratic equation, that is solved by this.
97 MFloat tauBase = F1B2 + 3.0 * (s->m_nu + nu_BIT) * FFPOW2(s->maxLevel() - s->a_level(pCellId));
98 const MFloat tau = F1B2
99 * (tauBase
100 + sqrt(tauBase * tauBase
101 + 2.0 * m_Cs * m_Cs * m_deltaX * m_deltaX * (F1BCSsq * F1BCSsq) * Q
102 * FPOW2(s->maxLevel() - s->a_level(pCellId))));
103
104 // eddy viscosity
105 if constexpr(bubble) {
106 if(s->m_cells.saveOldNu()) {
107 s->a_oldNu(pCellId) = s->a_nu(pCellId);
108 }
109 }
110 s->a_nu(pCellId) = (2.0 * tau - 1.0) * CSsq * FPOW2(s->maxLevel() - s->a_level(pCellId)) / 2.0;
111 if constexpr(bubble) {
112 if(s->m_isEELiquid) {
113 if(s->m_cells.saveOldNu()) {
114 s->a_oldNuT(pCellId) = s->a_nuT(pCellId);
115 }
116 s->a_nuT(pCellId) = s->a_nu(pCellId) - s->m_nu - nu_BIT;
117 }
118 }
119 });
120}
This class represents all LB models.
Definition: lbsolverdxqy.h:29
constexpr SysEqn sysEqn() const
Definition: lbsolverdxqy.h:194
MFloat m_Cs
Smagorinsky constant.
const MFloat m_deltaX
Filter width.
constexpr Real POW2(const Real x)
Definition: functions.h:119
constexpr MFloat FPOW2(MInt x)
constexpr MFloat FFPOW2(MInt x)
int32_t MInt
Definition: maiatypes.h:62
double MFloat
Definition: maiatypes.h:52
MFloat dist(const Point< DIM > &p, const Point< DIM > &q)
Definition: pointbox.h:54

◆ init()

template<MInt nDim, MInt nDist, class SysEqn , MBool bubble = false>
void maia::lb::LbSrcTerm_smagorinsky< nDim, nDist, SysEqn, bubble >::init ( )
inlineoverridevirtual

Implements maia::lb::LbSrcTerm< nDim, nDist, SysEqn >.

Definition at line 36 of file lbsrctermsgs.cpp.

36{/*NOT USED*/};

◆ readProperties()

template<MInt nDim, MInt nDist, class SysEqn , MBool bubble>
void maia::lb::LbSrcTerm_smagorinsky< nDim, nDist, SysEqn, bubble >::readProperties
overrideprotectedvirtual

Implements maia::lb::LbSrcTerm< nDim, nDist, SysEqn >.

Definition at line 44 of file lbsrctermsgs.cpp.

44 {
45 TRACE();
46 const MInt solverId = m_solver->m_solverId;
47 if(Context::propertyExists("smagorinskyConstant")) {
48 m_Cs = Context::getSolverProperty<MFloat>("smagorinskyConstant", solverId, AT_);
49 }
50}
static MBool propertyExists(const MString &name, MInt solver=m_noSolvers)
This function checks if a property exists in general.
Definition: context.cpp:494

Member Data Documentation

◆ m_Cs

template<MInt nDim, MInt nDist, class SysEqn , MBool bubble = false>
MFloat maia::lb::LbSrcTerm_smagorinsky< nDim, nDist, SysEqn, bubble >::m_Cs = 0.1
protected

Definition at line 26 of file lbsrctermsgs.cpp.

◆ m_deltaX

template<MInt nDim, MInt nDist, class SysEqn , MBool bubble = false>
const MFloat maia::lb::LbSrcTerm_smagorinsky< nDim, nDist, SysEqn, bubble >::m_deltaX = 1.0
protected

Definition at line 27 of file lbsrctermsgs.cpp.

◆ m_solver

template<MInt nDim, MInt nDist, class SysEqn , MBool bubble = false>
LbSolverDxQy<nDim, nDist, SysEqn>* maia::lb::LbSrcTerm_smagorinsky< nDim, nDist, SysEqn, bubble >::m_solver
protected

Definition at line 25 of file lbsrctermsgs.cpp.


The documentation for this class was generated from the following file: