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

Base class for lb sponge layers. More...

#include <lbsrctermsponge.h>

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

Classes

struct  SpongeCell
 

Public Member Functions

 LbSrcTerm_sponge (LbSolverDxQy< nDim, nDist, SysEqn > *p_solver)
 
void init () override
 Initialize properties common by all lb sponge terms. More...
 
- 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

MInt a_spongeCellId (const MInt cellId)
 Get index of cellId in m_spongeCells, returns -1 if unavailable. More...
 
void readProperties () override
 Reading basics properties common for all lb sponge terms. More...
 
virtual void readProperties ()=0
 

Protected Attributes

LbSolverDxQy< nDim, nDist, SysEqn > * m_solver
 
MBool m_isActive = false
 
std::vector< SpongeCellm_spongeCells
 
std::unordered_map< MInt, MIntm_mapCellId2SpongeCellId
 
MFloat m_spongeLayerThickness = 0.0
 
std::array< MFloat, 2 *nDimm_spongeThicknessFactor {}
 
MFloat m_spongeSigma = 0.0
 
MFloat m_spongeBeta = 2.0
 
MFloat m_trgRho = 1.0
 
std::array< MFloat, nDimm_trgU {}
 

Private Member Functions

void updateMapCellId2PmlCellId ()
 Update m_mapCellIds2PmlCellId. 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>
class maia::lb::LbSrcTerm_sponge< nDim, nDist, SysEqn >
Author
Miro Gondrum
Date
22.12.2023

Definition at line 22 of file lbsrctermsponge.h.

Constructor & Destructor Documentation

◆ LbSrcTerm_sponge()

template<MInt nDim, MInt nDist, class SysEqn >
maia::lb::LbSrcTerm_sponge< nDim, nDist, SysEqn >::LbSrcTerm_sponge ( LbSolverDxQy< nDim, nDist, SysEqn > *  p_solver)
inline

Definition at line 62 of file lbsrctermsponge.h.

62: m_solver(p_solver) { readProperties(); };
LbSolverDxQy< nDim, nDist, SysEqn > * m_solver
void readProperties() override
Reading basics properties common for all lb sponge terms.

Member Function Documentation

◆ a_spongeCellId()

template<MInt nDim, MInt nDist, class SysEqn >
MInt maia::lb::LbSrcTerm_sponge< nDim, nDist, SysEqn >::a_spongeCellId ( const MInt  cellId)
inlineprotected
Author
Miro Gondrum
Date
22.12.2023

Definition at line 50 of file lbsrctermsponge.h.

50 {
51 MInt spongeCellId = -1;
52 auto search = m_mapCellId2SpongeCellId.find(cellId);
53 if(search != m_mapCellId2SpongeCellId.end()) {
54 spongeCellId = search->second;
55 }
56 return spongeCellId;
57 };
std::unordered_map< MInt, MInt > m_mapCellId2SpongeCellId
int32_t MInt
Definition: maiatypes.h:62

◆ init()

template<MInt nDim, MInt nDist, class SysEqn >
void maia::lb::LbSrcTerm_sponge< nDim, nDist, SysEqn >::init
overridevirtual
Author
Miro Gondrum
Date
22.12.2023

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

Definition at line 135 of file lbsrctermsponge.h.

135 {
136 TRACE();
137 m_isActive = (m_spongeLayerThickness > 0.0) && m_solver->isActive();
138 if(!m_isActive) return;
139
140 const MFloat epsilon = m_solver->c_cellLengthAtLevel(m_solver->maxLevel()) / 1000.0;
141 std::array<MFloat, 2 * nDim> domainBoundaries{};
142 {
143 const MFloat halfCellWidth = F1B2 * m_solver->grid().cellLengthAtCell(0);
144 for(MInt i = 0; i < nDim; i++) {
145 domainBoundaries[2 * i + 0] = m_solver->a_coordinate(0, i) - halfCellWidth;
146 domainBoundaries[2 * i + 1] = m_solver->a_coordinate(0, i) + halfCellWidth;
147 }
148 }
149
150 for(MInt cellId = 0; cellId < m_solver->m_cells.size(); cellId++) {
151 if(!m_solver->c_isLeafCell(cellId)) continue;
152 const MFloat halfCellWidth = F1B2 * m_solver->grid().cellLengthAtCell(cellId);
153 for(MInt i = 0; i < nDim; i++) {
154 domainBoundaries[2 * i + 0] =
155 std::min(domainBoundaries[2 * i + 0], m_solver->a_coordinate(cellId, i) - halfCellWidth);
156 domainBoundaries[2 * i + 1] =
157 std::max(domainBoundaries[2 * i + 1], m_solver->a_coordinate(cellId, i) + halfCellWidth);
158 }
159 }
160 // exchange domain boundaries
161 std::array<MFloat, nDim> tmpMin{};
162 std::array<MFloat, nDim> tmpMax{};
163 for(MInt i = 0; i < nDim; i++) {
164 tmpMin[i] = domainBoundaries[2 * i + 0];
165 tmpMax[i] = domainBoundaries[2 * i + 1];
166 }
167 MPI_Allreduce(MPI_IN_PLACE, tmpMin.data(), nDim, maia::type_traits<MFloat>::mpiType(), MPI_MIN, m_solver->mpiComm(),
168 AT_, "MPI_IN_PLACE", "tmpMin");
169 MPI_Allreduce(MPI_IN_PLACE, tmpMax.data(), nDim, maia::type_traits<MFloat>::mpiType(), MPI_MAX, m_solver->mpiComm(),
170 AT_, "MPI_IN_PLACE", "tmpMax");
171 for(MInt i = 0; i < nDim; i++) {
172 domainBoundaries[2 * i + 0] = tmpMin[i];
173 domainBoundaries[2 * i + 1] = tmpMax[i];
174 }
175
176 // fill sponge coordinate
177 std::array<MFloat, 2 * nDim> spongeCoord{};
178 for(MInt i = 0; i < nDim; i++) {
179 spongeCoord[2 * i + 0] = domainBoundaries[2 * i + 0] + m_spongeLayerThickness * m_spongeThicknessFactor[2 * i + 0];
180 spongeCoord[2 * i + 1] = domainBoundaries[2 * i + 1] - m_spongeLayerThickness * m_spongeThicknessFactor[2 * i + 1];
181 }
182 std::array<MFloat, 2 * nDim> F1BspongeCoord{};
183 for(MInt i = 0; i < nDim; i++) {
184 F1BspongeCoord[2 * i + 0] = 1.0 / std::max(epsilon, m_spongeLayerThickness * m_spongeThicknessFactor[2 * i + 0]);
185 F1BspongeCoord[2 * i + 1] = 1.0 / std::max(epsilon, m_spongeLayerThickness * m_spongeThicknessFactor[2 * i + 1]);
186 }
187
188 // search sponge ids
189 for(MInt cellId = 0; cellId < m_solver->m_cells.size(); cellId++) {
190 if(m_solver->a_isHalo(cellId)) continue;
191
192 std::array<MFloat, nDim> dist{};
193 for(MInt i = 0; i < nDim; i++) {
194 const MUint distIdNegSide = 2 * i + 0;
195 const MUint distIdPosSide = 2 * i + 1;
196 const MFloat distNegSide =
197 (spongeCoord[distIdNegSide] - m_solver->a_coordinate(cellId, i)) * F1BspongeCoord[distIdNegSide];
198 const MFloat distPosSide =
199 (m_solver->a_coordinate(cellId, i) - spongeCoord[distIdPosSide]) * F1BspongeCoord[distIdPosSide];
200 dist[i] = std::max(distNegSide, distPosSide);
201 }
202
203 // choose the farthest sponge region
204 const MFloat relDist = *std::max_element(dist.begin(), dist.end());
205
206 if(relDist <= 0.0) continue;
207
208 const MFloat nodeDist = std::pow(relDist, m_spongeBeta);
209 const MFloat spongeFactor = m_spongeSigma * nodeDist;
210
211 constexpr MFloat epsSpongeFactor = 1e-15;
212 if(spongeFactor > epsSpongeFactor) {
213 auto& cell = m_spongeCells.emplace_back();
214 cell.cellId = cellId;
215 cell.spongeFactor = spongeFactor;
216 }
217 }
218 // Update map cellId->spongeCellId
220}
std::array< MFloat, 2 *nDim > m_spongeThicknessFactor
void updateMapCellId2PmlCellId()
Update m_mapCellIds2PmlCellId.
std::vector< SpongeCell > m_spongeCells
uint32_t MUint
Definition: maiatypes.h:63
double MFloat
Definition: maiatypes.h:52
int MPI_Allreduce(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, const MString &name, const MString &sndvarname, const MString &rcvvarname)
same as MPI_Allreduce
void const MInt cellId
Definition: collector.h:239
MFloat dist(const Point< DIM > &p, const Point< DIM > &q)
Definition: pointbox.h:54

◆ readProperties()

template<MInt nDim, MInt nDist, class SysEqn >
void maia::lb::LbSrcTerm_sponge< nDim, nDist, SysEqn >::readProperties
overrideprotectedvirtual
Author
Miro Gondrum
Date
22.12.2023

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

Definition at line 86 of file lbsrctermsponge.h.

86 {
87 TRACE();
88 const MInt solverId = m_solver->m_solverId;
89
90 std::stringstream ss;
91 ss << "--------------------------------------------------------------------------------" << std::endl;
92 ss << " Sponge properties" << std::endl;
93 ss << "--------------------------------------------------------------------------------" << std::endl;
94 ss << " [Sponge region]" << std::endl;
95
96 // sponge information
97 m_spongeSigma = Context::getSolverProperty<MFloat>("spongeSigma", solverId, AT_, &m_spongeSigma);
99 Context::getSolverProperty<MFloat>("spongeLayerThickness", solverId, AT_, &m_spongeLayerThickness);
100 ss << " spongeSigma: " << m_spongeSigma << std::endl;
101 ss << " spongeLayerThickness: " << m_spongeLayerThickness << std::endl;
102 ss << " spongeThicknessFactor: [ ";
103 for(MInt i = 0; i < 2 * nDim; i++) {
104 m_spongeThicknessFactor[i] = Context::getSolverProperty<MFloat>("spongeThicknessFactor", solverId, AT_, i);
105 ss << m_spongeThicknessFactor[i] << " ";
106 }
107 ss << "]" << std::endl;
108 ss << " [Sponge shape]" << std::endl;
109 m_spongeBeta = Context::getSolverProperty<MFloat>("spongeBeta", solverId, AT_, &m_spongeBeta);
110 ss << " spongeBeta: " << m_spongeBeta << std::endl;
111 ss << " [Sponge target state]" << std::endl;
112 m_trgRho = (m_solver->m_densityFluctuations) ? 0.0 : 1.0;
113 m_trgRho = Context::getSolverProperty<MFloat>("spongeTrgRho", solverId, AT_, &m_trgRho);
114 ss << " spongeTrgRho: " << m_trgRho << std::endl;
115 m_trgU.fill(0.0);
116 for(MInt i = 0; i < nDim; i++) {
117 m_trgU[i] = Context::getSolverProperty<MFloat>("spongeTrgU", solverId, AT_, &m_trgU[i], i) * m_solver->m_Ma * LBCS;
118 }
119 ss << " spongeTrgU: ";
120 for(MInt i = 0; i < nDim; i++) {
121 ss << m_trgU[i] << " ";
122 }
123 ss << std::endl;
124 ss << "--------------------------------------------------------------------------------" << std::endl;
125
126 if(m_solver->domainId() == 0) std::cout << ss.str();
127 m_log << ss.str();
128}
std::array< MFloat, nDim > m_trgU
InfoOutFile m_log

◆ updateMapCellId2PmlCellId()

template<MInt nDim, MInt nDist, class SysEqn >
void maia::lb::LbSrcTerm_sponge< nDim, nDist, SysEqn >::updateMapCellId2PmlCellId
private
Author
Miro Gondrum
Date
22.12.2023

Definition at line 72 of file lbsrctermsponge.h.

72 {
73 const MInt noSpongeCells = m_spongeCells.size();
75 m_mapCellId2SpongeCellId.reserve(noSpongeCells);
76 for(MInt spongeCellId = 0; spongeCellId < noSpongeCells; spongeCellId++) {
77 m_mapCellId2SpongeCellId.emplace(m_spongeCells[spongeCellId].cellId, spongeCellId);
78 }
79}

Member Data Documentation

◆ m_isActive

template<MInt nDim, MInt nDist, class SysEqn >
MBool maia::lb::LbSrcTerm_sponge< nDim, nDist, SysEqn >::m_isActive = false
protected

Definition at line 29 of file lbsrctermsponge.h.

◆ m_mapCellId2SpongeCellId

template<MInt nDim, MInt nDist, class SysEqn >
std::unordered_map<MInt, MInt> maia::lb::LbSrcTerm_sponge< nDim, nDist, SysEqn >::m_mapCellId2SpongeCellId
protected

Definition at line 37 of file lbsrctermsponge.h.

◆ m_solver

template<MInt nDim, MInt nDist, class SysEqn >
LbSolverDxQy<nDim, nDist, SysEqn>* maia::lb::LbSrcTerm_sponge< nDim, nDist, SysEqn >::m_solver
protected

Definition at line 27 of file lbsrctermsponge.h.

◆ m_spongeBeta

template<MInt nDim, MInt nDist, class SysEqn >
MFloat maia::lb::LbSrcTerm_sponge< nDim, nDist, SysEqn >::m_spongeBeta = 2.0
protected

Definition at line 42 of file lbsrctermsponge.h.

◆ m_spongeCells

template<MInt nDim, MInt nDist, class SysEqn >
std::vector<SpongeCell> maia::lb::LbSrcTerm_sponge< nDim, nDist, SysEqn >::m_spongeCells
protected

Definition at line 36 of file lbsrctermsponge.h.

◆ m_spongeLayerThickness

template<MInt nDim, MInt nDist, class SysEqn >
MFloat maia::lb::LbSrcTerm_sponge< nDim, nDist, SysEqn >::m_spongeLayerThickness = 0.0
protected

Definition at line 39 of file lbsrctermsponge.h.

◆ m_spongeSigma

template<MInt nDim, MInt nDist, class SysEqn >
MFloat maia::lb::LbSrcTerm_sponge< nDim, nDist, SysEqn >::m_spongeSigma = 0.0
protected

Definition at line 41 of file lbsrctermsponge.h.

◆ m_spongeThicknessFactor

template<MInt nDim, MInt nDist, class SysEqn >
std::array<MFloat, 2 * nDim> maia::lb::LbSrcTerm_sponge< nDim, nDist, SysEqn >::m_spongeThicknessFactor {}
protected

Definition at line 40 of file lbsrctermsponge.h.

◆ m_trgRho

template<MInt nDim, MInt nDist, class SysEqn >
MFloat maia::lb::LbSrcTerm_sponge< nDim, nDist, SysEqn >::m_trgRho = 1.0
protected

Definition at line 43 of file lbsrctermsponge.h.

◆ m_trgU

template<MInt nDim, MInt nDist, class SysEqn >
std::array<MFloat, nDim> maia::lb::LbSrcTerm_sponge< nDim, nDist, SysEqn >::m_trgU {}
protected

Definition at line 44 of file lbsrctermsponge.h.


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