Nektar++
Static Public Member Functions | Static Public Attributes | Protected Types | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes | 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)
 
SOLVER_UTILS_EXPORT void v_InitObject (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields, const TiXmlElement *pForce) override
 
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
 
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)
 
 ~ForcingAbsorption (void) override
 
void CalcAbsorption (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const TiXmlElement *pForce)
 

Private Attributes

bool m_homogeneous
 

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 55 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 78 of file ForcingAbsorption.h.

◆ BPointPair

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

Definition at line 79 of file ForcingAbsorption.h.

◆ BRTree

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

Definition at line 80 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 49 of file ForcingAbsorption.cpp.

52 : Forcing(pSession, pEquation), m_hasRefFlow(false), m_hasRefFlowTime(false)
53{
54}
SOLVER_UTILS_EXPORT Forcing(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation)
Constructor.
Definition: Forcing.cpp:48

◆ ~ForcingAbsorption()

Nektar::SolverUtils::ForcingAbsorption::~ForcingAbsorption ( void  )
inlineoverrideprivate

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 112 of file ForcingAbsorption.cpp.

116{
117 const TiXmlElement *funcNameElmt = pForce->FirstChildElement("COEFF");
118 ASSERTL0(funcNameElmt,
119 "Requires COEFF tag, specifying function "
120 "name which prescribes absorption layer coefficient.");
121 string funcName = funcNameElmt->GetText();
122 ASSERTL0(m_session->DefinesFunction(funcName),
123 "Function '" + funcName + "' not defined.");
124
125 int npts = pFields[0]->GetTotPoints();
126
127 m_Absorption = Array<OneD, Array<OneD, NekDouble>>(m_NumVariable);
128 for (int i = 0; i < m_NumVariable; ++i)
129 {
130 m_Absorption[i] = Array<OneD, NekDouble>(npts, 0.0);
131 }
132
133 funcNameElmt = pForce->FirstChildElement("BOUNDARYREGIONS");
134 if (funcNameElmt)
135 {
136 ASSERTL0(
137 ParseUtils::GenerateVector(funcNameElmt->GetText(), m_bRegions),
138 "Unable to process list of BOUNDARYREGIONS in Absorption "
139 "Forcing: " +
140 std::string(funcNameElmt->GetText()));
141
142 // alter m_bRegions so that it contains the boundaryRegions of this rank
143 std::vector<unsigned int> localBRegions;
144 SpatialDomains::BoundaryConditions bcs(m_session,
145 pFields[0]->GetGraph());
147 bcs.GetBoundaryRegions();
148 SpatialDomains::BoundaryRegionCollection::iterator it1;
149 int n = 0;
150 for (it1 = regions.begin(); it1 != regions.end(); ++it1)
151 {
152 if (std::find(m_bRegions.begin(), m_bRegions.end(), it1->first) !=
153 m_bRegions.end())
154 {
155 localBRegions.push_back(n);
156 }
157 n++;
158 }
159 m_bRegions = localBRegions;
160
161 if (m_bRegions.size() == 0)
162 {
163 return;
164 }
165
166 std::vector<Array<OneD, const NekDouble>> points;
167
168 Array<OneD, Array<OneD, NekDouble>> x(3);
169 for (int i = 0; i < 3; i++)
170 {
171 x[i] = Array<OneD, NekDouble>(npts, 0.0);
172 }
173 pFields[0]->GetCoords(x[0], x[1], x[2]);
174 for (int i = 0; i < 3; i++)
175 {
176 points.push_back(x[i]);
177 }
178
179 Array<OneD, NekDouble> t(npts, 0.0);
180 points.push_back(t);
181
182 Array<OneD, NekDouble> r(npts, 0.0);
183 std::vector<unsigned int>::iterator it;
184 std::vector<BPointPair> inPoints;
185 Array<OneD, Array<OneD, NekDouble>> b(3);
186 for (it = m_bRegions.begin(); it != m_bRegions.end(); ++it)
187 {
188 int bpts = pFields[0]->GetBndCondExpansions()[*it]->GetNpoints();
189 for (int i = 0; i < 3; i++)
190 {
191 b[i] = Array<OneD, NekDouble>(bpts, 0.0);
192 }
193 pFields[0]->GetBndCondExpansions()[*it]->GetCoords(b[0], b[1],
194 b[2]);
195 for (int i = 0;
196 i < pFields[0]->GetBndCondExpansions()[*it]->GetNpoints(); ++i)
197 {
198 inPoints.push_back(
199 BPointPair(BPoint(b[0][i], b[1][i], b[2][i]), i));
200 }
201 }
203 m_rtree->insert(inPoints.begin(), inPoints.end());
204
205 for (int i = 0; i < npts; ++i)
206 {
207 std::vector<BPointPair> result;
208 BPoint sPoint(x[0][i], x[1][i], x[2][i]);
209 m_rtree->query(bgi::nearest(sPoint, 1), std::back_inserter(result));
210 r[i] = bg::distance(sPoint, result[0].first);
211 }
212 points.push_back(r);
213
214 std::string s_FieldStr;
215 for (int i = 0; i < m_NumVariable; ++i)
216 {
217 s_FieldStr = m_session->GetVariable(i);
218 ASSERTL0(m_session->DefinesFunction(funcName, s_FieldStr),
219 "Variable '" + s_FieldStr + "' not defined.");
220
222 m_session->GetFunction(funcName, s_FieldStr);
223 ASSERTL0(ffunc->GetVlist() == "x y z t r",
224 "EVARS of " + funcName + " must be 'r'");
225
226 ffunc->Evaluate(points, m_Absorption[i]);
227 }
228 }
229 else
230 {
231 for (int i = 0; i < m_NumVariable; ++i)
232 {
233 std::string s_FieldStr = m_session->GetVariable(i);
234 GetFunction(pFields, m_session, funcName)
235 ->Evaluate(s_FieldStr, m_Absorption[i]);
236 }
237 }
238
239 // If homogeneous expansion is used, transform the forcing term to
240 // be in the Fourier space
241 if (m_homogeneous)
242 {
243 for (int i = 0; i < m_NumVariable; ++i)
244 {
245 pFields[i]->HomogeneousFwdTrans(pFields[i]->GetTotPoints(),
247 }
248 }
249}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
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:130
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:121
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:192
LibUtilities::SessionReaderSharedPtr m_session
Session reader.
Definition: Forcing.h:115
std::shared_ptr< Equation > EquationSharedPtr
Definition: Equation.h:125
std::map< int, BoundaryRegionShPtr > BoundaryRegionCollection
Definition: Conditions.h:210
InputIterator find(InputIterator first, InputIterator last, InputIterator startingpoint, const EqualityComparable &value)
Definition: StdRegions.hpp:447

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::StdRegions::find(), Nektar::ParseUtils::GenerateVector(), Nektar::SpatialDomains::BoundaryConditions::GetBoundaryRegions(), Nektar::SolverUtils::Forcing::GetFunction(), m_Absorption, m_bRegions, m_homogeneous, 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 280 of file ForcingAbsorption.cpp.

283{
284 int nq = m_Forcing[0].size();
285
286 std::string s_FieldStr;
287 Array<OneD, NekDouble> TimeScale(1);
288 Array<OneD, Array<OneD, NekDouble>> RefflowScaled(m_NumVariable);
289
290 if (m_hasRefFlow)
291 {
292 for (int i = 0; i < m_NumVariable; i++)
293 {
294 RefflowScaled[i] = Array<OneD, NekDouble>(nq);
296 {
297 s_FieldStr = m_session->GetVariable(i);
298
299 std::string s_FieldStr = m_session->GetVariable(i);
301 ->Evaluate(s_FieldStr, m_Refflow[i], time);
302 Vmath::Vcopy(nq, m_Refflow[i], 1, RefflowScaled[i], 1);
303 }
304 else
305 {
306 Vmath::Vcopy(nq, m_Refflow[i], 1, RefflowScaled[i], 1);
307 }
308
309 Vmath::Vsub(nq, inarray[i], 1, RefflowScaled[i], 1, m_Forcing[i],
310 1);
311 Vmath::Vmul(nq, m_Absorption[i], 1, m_Forcing[i], 1, m_Forcing[i],
312 1);
313 }
314 }
315 else
316 {
317 for (int i = 0; i < m_NumVariable; i++)
318 {
319 Vmath::Vmul(nq, m_Absorption[i], 1, inarray[i], 1, m_Forcing[i], 1);
320 }
321 }
322}
Array< OneD, Array< OneD, NekDouble > > m_Refflow
Array< OneD, Array< OneD, NekDouble > > m_Forcing
Evaluated forcing function.
Definition: Forcing.h:119
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.hpp:72
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
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.hpp:220

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 61 of file ForcingAbsorption.h.

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

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 251 of file ForcingAbsorption.cpp.

255{
256 int nq = m_Forcing[0].size();
257 CalculateForcing(fields, inarray, time);
258 for (int i = 0; i < m_NumVariable; ++i)
259 {
260 Vmath::Vadd(nq, m_Forcing[i], 1, outarray[i], 1, outarray[i], 1);
261 }
262}
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.hpp:180

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 264 of file ForcingAbsorption.cpp.

268{
269 int ncoeff = outarray[m_NumVariable - 1].size();
270 Array<OneD, NekDouble> tmp(ncoeff, 0.0);
271 CalculateForcing(fields, inarray, time);
272
273 for (int i = 0; i < m_NumVariable; ++i)
274 {
275 fields[i]->FwdTrans(m_Forcing[i], tmp);
276 Vmath::Vadd(ncoeff, tmp, 1, outarray[i], 1, outarray[i], 1);
277 }
278}

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 56 of file ForcingAbsorption.cpp.

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

References ASSERTL0, CalcAbsorption(), Nektar::MultiRegions::e3DH1D, Nektar::MultiRegions::e3DH2D, Nektar::SolverUtils::Forcing::GetFunction(), Nektar::SolverUtils::Forcing::m_Forcing, m_funcNameTime, m_hasRefFlow, m_hasRefFlowTime, m_homogeneous, 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:197
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:42

Name of the class.

Definition at line 75 of file ForcingAbsorption.h.

◆ m_Absorption

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

Definition at line 84 of file ForcingAbsorption.h.

Referenced by CalcAbsorption(), and CalculateForcing().

◆ m_bRegions

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

Definition at line 87 of file ForcingAbsorption.h.

Referenced by CalcAbsorption().

◆ m_funcNameTime

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

Definition at line 86 of file ForcingAbsorption.h.

Referenced by CalculateForcing(), and v_InitObject().

◆ m_hasRefFlow

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

Definition at line 82 of file ForcingAbsorption.h.

Referenced by CalculateForcing(), and v_InitObject().

◆ m_hasRefFlowTime

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

Definition at line 83 of file ForcingAbsorption.h.

Referenced by CalculateForcing(), and v_InitObject().

◆ m_homogeneous

bool Nektar::SolverUtils::ForcingAbsorption::m_homogeneous
private

Definition at line 114 of file ForcingAbsorption.h.

Referenced by CalcAbsorption(), and v_InitObject().

◆ m_Refflow

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

Definition at line 85 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 88 of file ForcingAbsorption.h.

Referenced by CalcAbsorption().