Nektar++
Static Public Member Functions | Static Public Attributes | Protected Types | Protected Member Functions | Protected Attributes | Private Member Functions | Friends | List of all members
Nektar::SolverUtils::ForcingAbsorption Class Reference

#include <ForcingAbsorption.h>

Inheritance diagram for Nektar::SolverUtils::ForcingAbsorption:
[legend]

Static Public Member Functions

static SOLVER_UTILS_EXPORT ForcingSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields, const TiXmlElement *pForce)
 Creates an instance of this class. More...
 
- Static Public Member Functions inherited from Nektar::SolverUtils::Forcing
static SOLVER_UTILS_EXPORT std::vector< ForcingSharedPtrLoad (const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields=0)
 

Static Public Attributes

static std::string className
 Name of the class. More...
 

Protected Types

typedef bg::model::point< NekDouble, 3, bg::cs::cartesian > BPoint
 
typedef std::pair< BPoint, unsigned int > BPointPair
 
typedef bgi::rtree< BPointPair, bgi::rstar< 16 > > BRTree
 

Protected Member Functions

void CalculateForcing (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, const NekDouble &time)
 
virtual SOLVER_UTILS_EXPORT void v_InitObject (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields, const TiXmlElement *pForce) override
 
virtual SOLVER_UTILS_EXPORT void v_Apply (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble &time) override
 
virtual SOLVER_UTILS_EXPORT void v_ApplyCoeff (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble &time) override
 
- Protected Member Functions inherited from Nektar::SolverUtils::Forcing
SOLVER_UTILS_EXPORT Forcing (const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation)
 Constructor. More...
 
virtual SOLVER_UTILS_EXPORT void v_InitObject (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields, const TiXmlElement *pForce)=0
 
virtual SOLVER_UTILS_EXPORT void v_Apply (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble &time)=0
 
virtual SOLVER_UTILS_EXPORT void v_PreApply (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble &time)
 
virtual SOLVER_UTILS_EXPORT void v_ApplyCoeff (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble &time)
 
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const LibUtilities::SessionReaderSharedPtr &pSession, std::string pName, bool pCache=false)
 Get a SessionFunction by name. More...
 
SOLVER_UTILS_EXPORT void EvaluateTimeFunction (LibUtilities::SessionReaderSharedPtr pSession, std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, NekDouble pTime=NekDouble(0))
 
SOLVER_UTILS_EXPORT void EvaluateTimeFunction (const NekDouble pTime, const LibUtilities::EquationSharedPtr &pEqn, Array< OneD, NekDouble > &pArray)
 

Protected Attributes

bool m_hasRefFlow
 
bool m_hasRefFlowTime
 
Array< OneD, Array< OneD, NekDouble > > m_Absorption
 
Array< OneD, Array< OneD, NekDouble > > m_Refflow
 
std::string m_funcNameTime
 
std::vector< unsigned int > m_bRegions
 
std::shared_ptr< BRTreem_rtree
 
- Protected Attributes inherited from Nektar::SolverUtils::Forcing
LibUtilities::SessionReaderSharedPtr m_session
 Session reader. More...
 
const std::weak_ptr< EquationSystemm_equ
 Weak pointer to equation system using this forcing. More...
 
Array< OneD, Array< OneD, NekDouble > > m_Forcing
 Evaluated forcing function. More...
 
int m_NumVariable
 Number of variables. More...
 
std::map< std::string, SolverUtils::SessionFunctionSharedPtrm_sessionFunctions
 Map of known SessionFunctions. More...
 

Private Member Functions

 ForcingAbsorption (const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation)
 
virtual ~ForcingAbsorption (void)
 
void CalcAbsorption (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const TiXmlElement *pForce)
 

Friends

class MemoryManager< ForcingAbsorption >
 

Additional Inherited Members

- Public Member Functions inherited from Nektar::SolverUtils::Forcing
virtual SOLVER_UTILS_EXPORT ~Forcing ()
 
SOLVER_UTILS_EXPORT void InitObject (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields, const TiXmlElement *pForce)
 Initialise the forcing object. More...
 
SOLVER_UTILS_EXPORT void Apply (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble &time)
 Apply the forcing. More...
 
SOLVER_UTILS_EXPORT void PreApply (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble &time)
 Change the advection velocity before applying the forcing. For example, subtracting the frame velocity from the advection velocity in the MovingRefercenceFrame. More...
 
SOLVER_UTILS_EXPORT void ApplyCoeff (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble &time)
 Apply the forcing. More...
 
SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & GetForces ()
 
SOLVER_UTILS_EXPORT Array< OneD, Array< OneD, NekDouble > > & UpdateForces ()
 

Detailed Description

Definition at line 57 of file ForcingAbsorption.h.

Member Typedef Documentation

◆ BPoint

typedef bg::model::point<NekDouble, 3, bg::cs::cartesian> Nektar::SolverUtils::ForcingAbsorption::BPoint
protected

Definition at line 80 of file ForcingAbsorption.h.

◆ BPointPair

typedef std::pair<BPoint, unsigned int> Nektar::SolverUtils::ForcingAbsorption::BPointPair
protected

Definition at line 81 of file ForcingAbsorption.h.

◆ BRTree

typedef bgi::rtree<BPointPair, bgi::rstar<16> > Nektar::SolverUtils::ForcingAbsorption::BRTree
protected

Definition at line 82 of file ForcingAbsorption.h.

Constructor & Destructor Documentation

◆ ForcingAbsorption()

Nektar::SolverUtils::ForcingAbsorption::ForcingAbsorption ( const LibUtilities::SessionReaderSharedPtr pSession,
const std::weak_ptr< EquationSystem > &  pEquation 
)
private

Definition at line 53 of file ForcingAbsorption.cpp.

56 : Forcing(pSession, pEquation), m_hasRefFlow(false), m_hasRefFlowTime(false)
57{
58}
SOLVER_UTILS_EXPORT Forcing(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation)
Constructor.
Definition: Forcing.cpp:50

◆ ~ForcingAbsorption()

virtual Nektar::SolverUtils::ForcingAbsorption::~ForcingAbsorption ( void  )
inlineprivatevirtual

Definition at line 117 of file ForcingAbsorption.h.

117{};

Member Function Documentation

◆ CalcAbsorption()

void Nektar::SolverUtils::ForcingAbsorption::CalcAbsorption ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields,
const TiXmlElement *  pForce 
)
private

Definition at line 111 of file ForcingAbsorption.cpp.

115{
116 const TiXmlElement *funcNameElmt = pForce->FirstChildElement("COEFF");
117 ASSERTL0(funcNameElmt,
118 "Requires COEFF tag, specifying function "
119 "name which prescribes absorption layer coefficient.");
120 string funcName = funcNameElmt->GetText();
121 ASSERTL0(m_session->DefinesFunction(funcName),
122 "Function '" + funcName + "' not defined.");
123
124 int npts = pFields[0]->GetTotPoints();
125
126 m_Absorption = Array<OneD, Array<OneD, NekDouble>>(m_NumVariable);
127 for (int i = 0; i < m_NumVariable; ++i)
128 {
129 m_Absorption[i] = Array<OneD, NekDouble>(npts, 0.0);
130 }
131
132 funcNameElmt = pForce->FirstChildElement("BOUNDARYREGIONS");
133 if (funcNameElmt)
134 {
135 ASSERTL0(
136 ParseUtils::GenerateVector(funcNameElmt->GetText(), m_bRegions),
137 "Unable to process list of BOUNDARYREGIONS in Absorption "
138 "Forcing: " +
139 std::string(funcNameElmt->GetText()));
140
141 // alter m_bRegions so that it contains the boundaryRegions of this rank
142 std::vector<unsigned int> localBRegions;
143 SpatialDomains::BoundaryConditions bcs(m_session,
144 pFields[0]->GetGraph());
146 bcs.GetBoundaryRegions();
147 SpatialDomains::BoundaryRegionCollection::iterator it1;
148 int n = 0;
149 for (it1 = regions.begin(); it1 != regions.end(); ++it1)
150 {
151 if (std::find(m_bRegions.begin(), m_bRegions.end(), it1->first) !=
152 m_bRegions.end())
153 {
154 localBRegions.push_back(n);
155 }
156 n++;
157 }
158 m_bRegions = localBRegions;
159
160 if (m_bRegions.size() == 0)
161 {
162 return;
163 }
164
165 std::vector<Array<OneD, const NekDouble>> points;
166
167 Array<OneD, Array<OneD, NekDouble>> x(3);
168 for (int i = 0; i < 3; i++)
169 {
170 x[i] = Array<OneD, NekDouble>(npts, 0.0);
171 }
172 pFields[0]->GetCoords(x[0], x[1], x[2]);
173 for (int i = 0; i < 3; i++)
174 {
175 points.push_back(x[i]);
176 }
177
178 Array<OneD, NekDouble> t(npts, 0.0);
179 points.push_back(t);
180
181 Array<OneD, NekDouble> r(npts, 0.0);
182 std::vector<unsigned int>::iterator it;
183 std::vector<BPointPair> inPoints;
184 Array<OneD, Array<OneD, NekDouble>> b(3);
185 for (it = m_bRegions.begin(); it != m_bRegions.end(); ++it)
186 {
187 int bpts = pFields[0]->GetBndCondExpansions()[*it]->GetNpoints();
188 for (int i = 0; i < 3; i++)
189 {
190 b[i] = Array<OneD, NekDouble>(bpts, 0.0);
191 }
192 pFields[0]->GetBndCondExpansions()[*it]->GetCoords(b[0], b[1],
193 b[2]);
194 for (int i = 0;
195 i < pFields[0]->GetBndCondExpansions()[*it]->GetNpoints(); ++i)
196 {
197 inPoints.push_back(
198 BPointPair(BPoint(b[0][i], b[1][i], b[2][i]), i));
199 }
200 }
202 m_rtree->insert(inPoints.begin(), inPoints.end());
203
204 for (int i = 0; i < npts; ++i)
205 {
206 std::vector<BPointPair> result;
207 BPoint sPoint(x[0][i], x[1][i], x[2][i]);
208 m_rtree->query(bgi::nearest(sPoint, 1), std::back_inserter(result));
209 r[i] = bg::distance(sPoint, result[0].first);
210 }
211 points.push_back(r);
212
213 std::string s_FieldStr;
214 for (int i = 0; i < m_NumVariable; ++i)
215 {
216 s_FieldStr = m_session->GetVariable(i);
217 ASSERTL0(m_session->DefinesFunction(funcName, s_FieldStr),
218 "Variable '" + s_FieldStr + "' not defined.");
219
221 m_session->GetFunction(funcName, s_FieldStr);
222 ASSERTL0(ffunc->GetVlist() == "x y z t r",
223 "EVARS of " + funcName + " must be 'r'");
224
225 ffunc->Evaluate(points, m_Absorption[i]);
226 }
227 }
228 else
229 {
230 for (int i = 0; i < m_NumVariable; ++i)
231 {
232 std::string s_FieldStr = m_session->GetVariable(i);
233 GetFunction(pFields, m_session, funcName)
234 ->Evaluate(s_FieldStr, m_Absorption[i]);
235 }
236 }
237}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
static bool GenerateVector(const std::string &str, std::vector< T > &out)
Takes a comma-separated string and converts it to entries in a vector.
Definition: ParseUtils.cpp:131
bg::model::point< NekDouble, 3, bg::cs::cartesian > BPoint
std::pair< BPoint, unsigned int > BPointPair
std::vector< unsigned int > m_bRegions
Array< OneD, Array< OneD, NekDouble > > m_Absorption
int m_NumVariable
Number of variables.
Definition: Forcing.h:123
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const LibUtilities::SessionReaderSharedPtr &pSession, std::string pName, bool pCache=false)
Get a SessionFunction by name.
Definition: Forcing.cpp:194
LibUtilities::SessionReaderSharedPtr m_session
Session reader.
Definition: Forcing.h:117
std::shared_ptr< Equation > EquationSharedPtr
Definition: Equation.h:129
std::map< int, BoundaryRegionShPtr > BoundaryRegionCollection
Definition: Conditions.h:210
InputIterator find(InputIterator first, InputIterator last, InputIterator startingpoint, const EqualityComparable &value)
Definition: StdRegions.hpp:453

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::StdRegions::find(), Nektar::ParseUtils::GenerateVector(), Nektar::SpatialDomains::BoundaryConditions::GetBoundaryRegions(), Nektar::SolverUtils::Forcing::GetFunction(), m_Absorption, m_bRegions, Nektar::SolverUtils::Forcing::m_NumVariable, m_rtree, and Nektar::SolverUtils::Forcing::m_session.

Referenced by v_InitObject().

◆ CalculateForcing()

void Nektar::SolverUtils::ForcingAbsorption::CalculateForcing ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble > > &  inarray,
const NekDouble time 
)
protected

Definition at line 268 of file ForcingAbsorption.cpp.

271{
272 boost::ignore_unused(fields);
273 int nq = m_Forcing[0].size();
274
275 std::string s_FieldStr;
276 Array<OneD, NekDouble> TimeScale(1);
277 Array<OneD, Array<OneD, NekDouble>> RefflowScaled(m_NumVariable);
278
279 if (m_hasRefFlow)
280 {
281 for (int i = 0; i < m_NumVariable; i++)
282 {
283 RefflowScaled[i] = Array<OneD, NekDouble>(nq);
285 {
286 s_FieldStr = m_session->GetVariable(i);
287
288 std::string s_FieldStr = m_session->GetVariable(i);
290 ->Evaluate(s_FieldStr, m_Refflow[i], time);
291 Vmath::Vcopy(nq, m_Refflow[i], 1, RefflowScaled[i], 1);
292 }
293 else
294 {
295 Vmath::Vcopy(nq, m_Refflow[i], 1, RefflowScaled[i], 1);
296 }
297
298 Vmath::Vsub(nq, inarray[i], 1, RefflowScaled[i], 1, m_Forcing[i],
299 1);
300 Vmath::Vmul(nq, m_Absorption[i], 1, m_Forcing[i], 1, m_Forcing[i],
301 1);
302 }
303 }
304 else
305 {
306 for (int i = 0; i < m_NumVariable; i++)
307 {
308 Vmath::Vmul(nq, m_Absorption[i], 1, inarray[i], 1, m_Forcing[i], 1);
309 }
310 }
311}
Array< OneD, Array< OneD, NekDouble > > m_Refflow
Array< OneD, Array< OneD, NekDouble > > m_Forcing
Evaluated forcing function.
Definition: Forcing.h:121
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:207
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:414

References Nektar::SolverUtils::Forcing::GetFunction(), m_Absorption, Nektar::SolverUtils::Forcing::m_Forcing, m_funcNameTime, m_hasRefFlow, m_hasRefFlowTime, Nektar::SolverUtils::Forcing::m_NumVariable, m_Refflow, Nektar::SolverUtils::Forcing::m_session, Vmath::Vcopy(), Vmath::Vmul(), and Vmath::Vsub().

Referenced by v_Apply(), and v_ApplyCoeff().

◆ create()

static SOLVER_UTILS_EXPORT ForcingSharedPtr Nektar::SolverUtils::ForcingAbsorption::create ( const LibUtilities::SessionReaderSharedPtr pSession,
const std::weak_ptr< EquationSystem > &  pEquation,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields,
const unsigned int &  pNumForcingFields,
const TiXmlElement *  pForce 
)
inlinestatic

Creates an instance of this class.

Definition at line 63 of file ForcingAbsorption.h.

68 {
71 pEquation);
72 p->InitObject(pFields, pNumForcingFields, pForce);
73 return p;
74 }
SOLVER_UTILS_EXPORT typedef std::shared_ptr< Forcing > ForcingSharedPtr
A shared pointer to an EquationSystem object.
Definition: Forcing.h:55

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::SolverUtils::ForcingSharedPtr, and CellMLToNektar.cellml_metadata::p.

◆ v_Apply()

void Nektar::SolverUtils::ForcingAbsorption::v_Apply ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble time 
)
overrideprotectedvirtual

Implements Nektar::SolverUtils::Forcing.

Definition at line 239 of file ForcingAbsorption.cpp.

243{
244 int nq = m_Forcing[0].size();
245 CalculateForcing(fields, inarray, time);
246 for (int i = 0; i < m_NumVariable; ++i)
247 {
248 Vmath::Vadd(nq, m_Forcing[i], 1, outarray[i], 1, outarray[i], 1);
249 }
250}
void CalculateForcing(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, const NekDouble &time)
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:354

References CalculateForcing(), Nektar::SolverUtils::Forcing::m_Forcing, Nektar::SolverUtils::Forcing::m_NumVariable, and Vmath::Vadd().

◆ v_ApplyCoeff()

void Nektar::SolverUtils::ForcingAbsorption::v_ApplyCoeff ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble time 
)
overrideprotectedvirtual

Reimplemented from Nektar::SolverUtils::Forcing.

Definition at line 252 of file ForcingAbsorption.cpp.

256{
257 int ncoeff = outarray[m_NumVariable - 1].size();
258 Array<OneD, NekDouble> tmp(ncoeff, 0.0);
259 CalculateForcing(fields, inarray, time);
260
261 for (int i = 0; i < m_NumVariable; ++i)
262 {
263 fields[i]->FwdTrans(m_Forcing[i], tmp);
264 Vmath::Vadd(ncoeff, tmp, 1, outarray[i], 1, outarray[i], 1);
265 }
266}

References CalculateForcing(), Nektar::SolverUtils::Forcing::m_Forcing, Nektar::SolverUtils::Forcing::m_NumVariable, and Vmath::Vadd().

◆ v_InitObject()

void Nektar::SolverUtils::ForcingAbsorption::v_InitObject ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields,
const unsigned int &  pNumForcingFields,
const TiXmlElement *  pForce 
)
overrideprotectedvirtual

Implements Nektar::SolverUtils::Forcing.

Definition at line 60 of file ForcingAbsorption.cpp.

63{
64 m_NumVariable = pNumForcingFields;
65 int npts = pFields[0]->GetTotPoints();
66
67 CalcAbsorption(pFields, pForce);
68
69 m_Forcing = Array<OneD, Array<OneD, NekDouble>>(m_NumVariable);
70 for (int i = 0; i < m_NumVariable; ++i)
71 {
72 m_Forcing[i] = Array<OneD, NekDouble>(npts, 0.0);
73 }
74
75 const TiXmlElement *funcNameElmt = pForce->FirstChildElement("REFFLOW");
76 if (funcNameElmt)
77 {
78 string funcName = funcNameElmt->GetText();
79 ASSERTL0(m_session->DefinesFunction(funcName),
80 "Function '" + funcName + "' not defined.");
81 m_Refflow = Array<OneD, Array<OneD, NekDouble>>(m_NumVariable);
82 for (int i = 0; i < m_NumVariable; ++i)
83 {
84 std::string s_FieldStr = m_session->GetVariable(i);
85 ASSERTL0(m_session->DefinesFunction(funcName, s_FieldStr),
86 "Variable '" + s_FieldStr + "' not defined.");
87 m_Refflow[i] = Array<OneD, NekDouble>(npts, 0.0);
88 GetFunction(pFields, m_session, funcName)
89 ->Evaluate(s_FieldStr, m_Refflow[i]);
90 }
91 m_hasRefFlow = true;
92 }
93
94 funcNameElmt = pForce->FirstChildElement("REFFLOWTIME");
95 if (funcNameElmt)
96 {
97 m_funcNameTime = funcNameElmt->GetText();
98 m_hasRefFlowTime = true;
99 m_hasRefFlow = true;
100 m_Refflow = Array<OneD, Array<OneD, NekDouble>>(m_NumVariable);
101 for (int i = 0; i < m_NumVariable; ++i)
102 {
103 std::string s_FieldStr = m_session->GetVariable(i);
104 ASSERTL0(m_session->DefinesFunction(m_funcNameTime, s_FieldStr),
105 "Variable '" + s_FieldStr + "' not defined.");
106 m_Refflow[i] = Array<OneD, NekDouble>(npts, 0.0);
107 }
108 }
109}
void CalcAbsorption(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const TiXmlElement *pForce)

References ASSERTL0, CalcAbsorption(), Nektar::SolverUtils::Forcing::GetFunction(), Nektar::SolverUtils::Forcing::m_Forcing, m_funcNameTime, m_hasRefFlow, m_hasRefFlowTime, Nektar::SolverUtils::Forcing::m_NumVariable, m_Refflow, and Nektar::SolverUtils::Forcing::m_session.

Friends And Related Function Documentation

◆ MemoryManager< ForcingAbsorption >

friend class MemoryManager< ForcingAbsorption >
friend

Definition at line 1 of file ForcingAbsorption.h.

Member Data Documentation

◆ className

std::string Nektar::SolverUtils::ForcingAbsorption::className
static
Initial value:
=
"Absorption", ForcingAbsorption::create, "Forcing Absorption")
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
static SOLVER_UTILS_EXPORT ForcingSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields, const TiXmlElement *pForce)
Creates an instance of this class.
ForcingFactory & GetForcingFactory()
Declaration of the forcing factory singleton.
Definition: Forcing.cpp:44

Name of the class.

Definition at line 77 of file ForcingAbsorption.h.

◆ m_Absorption

Array<OneD, Array<OneD, NekDouble> > Nektar::SolverUtils::ForcingAbsorption::m_Absorption
protected

Definition at line 86 of file ForcingAbsorption.h.

Referenced by CalcAbsorption(), and CalculateForcing().

◆ m_bRegions

std::vector<unsigned int> Nektar::SolverUtils::ForcingAbsorption::m_bRegions
protected

Definition at line 89 of file ForcingAbsorption.h.

Referenced by CalcAbsorption().

◆ m_funcNameTime

std::string Nektar::SolverUtils::ForcingAbsorption::m_funcNameTime
protected

Definition at line 88 of file ForcingAbsorption.h.

Referenced by CalculateForcing(), and v_InitObject().

◆ m_hasRefFlow

bool Nektar::SolverUtils::ForcingAbsorption::m_hasRefFlow
protected

Definition at line 84 of file ForcingAbsorption.h.

Referenced by CalculateForcing(), and v_InitObject().

◆ m_hasRefFlowTime

bool Nektar::SolverUtils::ForcingAbsorption::m_hasRefFlowTime
protected

Definition at line 85 of file ForcingAbsorption.h.

Referenced by CalculateForcing(), and v_InitObject().

◆ m_Refflow

Array<OneD, Array<OneD, NekDouble> > Nektar::SolverUtils::ForcingAbsorption::m_Refflow
protected

Definition at line 87 of file ForcingAbsorption.h.

Referenced by CalculateForcing(), and v_InitObject().

◆ m_rtree

std::shared_ptr<BRTree> Nektar::SolverUtils::ForcingAbsorption::m_rtree
protected

Definition at line 90 of file ForcingAbsorption.h.

Referenced by CalcAbsorption().