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)
 
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)
 
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)
 
- 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_PreApply (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 81 of file ForcingAbsorption.h.

◆ BPointPair

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

Definition at line 82 of file ForcingAbsorption.h.

◆ BRTree

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

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

58  : Forcing(pSession, pEquation),
59  m_hasRefFlow(false),
60  m_hasRefFlowTime(false)
61  {
62  }
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 115 of file ForcingAbsorption.cpp.

119  {
120  const TiXmlElement *funcNameElmt = pForce->FirstChildElement("COEFF");
121  ASSERTL0(funcNameElmt,
122  "Requires COEFF tag, specifying function "
123  "name which prescribes absorption layer coefficient.");
124  string funcName = funcNameElmt->GetText();
125  ASSERTL0(m_session->DefinesFunction(funcName),
126  "Function '" + funcName + "' not defined.");
127 
128  int npts = pFields[0]->GetTotPoints();
129 
130  m_Absorption = Array<OneD, Array<OneD, NekDouble> >(m_NumVariable);
131  for (int i = 0; i < m_NumVariable; ++i)
132  {
133  m_Absorption[i] = Array<OneD, NekDouble>(npts, 0.0);
134  }
135 
136  funcNameElmt = pForce->FirstChildElement("BOUNDARYREGIONS");
137  if (funcNameElmt)
138  {
139  ASSERTL0(ParseUtils::GenerateVector(funcNameElmt->GetText(),
140  m_bRegions),
141  "Unable to process list of BOUNDARYREGIONS in Absorption "
142  "Forcing: " +
143  std::string(funcNameElmt->GetText()));
144 
145  // alter m_bRegions so that it contains the boundaryRegions of this rank
146  std::vector<unsigned int> localBRegions;
147  SpatialDomains::BoundaryConditions bcs(m_session, pFields[0]->GetGraph());
148  SpatialDomains::BoundaryRegionCollection regions = bcs.GetBoundaryRegions();
149  SpatialDomains::BoundaryRegionCollection::iterator it1;
150  int n = 0;
151  for (it1 = regions.begin(); it1 != regions.end(); ++it1)
152  {
153  if (std::find(m_bRegions.begin(), m_bRegions.end(), it1->first) != 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(
194  b[0], b[1], b[2]);
195  for (int i = 0;
196  i < pFields[0]->GetBndCondExpansions()[*it]->GetNpoints();
197  ++i)
198  {
199  inPoints.push_back(
200  BPointPair(BPoint(b[0][i], b[1][i], b[2][i]), i));
201  }
202  }
204  m_rtree->insert(inPoints.begin(), inPoints.end());
205 
206  for (int i = 0; i < npts; ++i)
207  {
208  std::vector<BPointPair> result;
209  BPoint sPoint(x[0][i], x[1][i], x[2][i]);
210  m_rtree->query(bgi::nearest(sPoint, 1), std::back_inserter(result));
211  r[i] = bg::distance(sPoint, result[0].first);
212  }
213  points.push_back(r);
214 
215  std::string s_FieldStr;
216  for (int i = 0; i < m_NumVariable; ++i)
217  {
218  s_FieldStr = m_session->GetVariable(i);
219  ASSERTL0(m_session->DefinesFunction(funcName, s_FieldStr),
220  "Variable '" + s_FieldStr + "' not defined.");
221 
223  m_session->GetFunction(funcName, s_FieldStr);
224  ASSERTL0(ffunc->GetVlist() == "x y z t r",
225  "EVARS of " + funcName + " must be 'r'");
226 
227  ffunc->Evaluate(points, m_Absorption[i]);
228  }
229  }
230  else
231  {
232  for (int i = 0; i < m_NumVariable; ++i)
233  {
234  std::string s_FieldStr = m_session->GetVariable(i);
235  GetFunction(pFields, m_session, funcName)->Evaluate(s_FieldStr, m_Absorption[i]);
236  }
237  }
238  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
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:135
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:124
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:202
LibUtilities::SessionReaderSharedPtr m_session
Session reader.
Definition: Forcing.h:118
std::shared_ptr< Equation > EquationSharedPtr
Definition: Equation.h:131
std::map< int, BoundaryRegionShPtr > BoundaryRegionCollection
Definition: Conditions.h:217
InputIterator find(InputIterator first, InputIterator last, InputIterator startingpoint, const EqualityComparable &value)
Definition: StdRegions.hpp:362

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

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

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.

69  {
71  AllocateSharedPtr(pSession, 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:52

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 
)
protectedvirtual

Implements Nektar::SolverUtils::Forcing.

Definition at line 240 of file ForcingAbsorption.cpp.

245  {
246  int nq = m_Forcing[0].size();
247  CalculateForcing(fields, inarray, time);
248  for (int i = 0; i < m_NumVariable; ++i)
249  {
250  Vmath::Vadd(nq, m_Forcing[i], 1,
251  outarray[i], 1, outarray[i], 1);
252  }
253  }
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:322

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 
)
protectedvirtual

Reimplemented from Nektar::SolverUtils::Forcing.

Definition at line 255 of file ForcingAbsorption.cpp.

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

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 
)
protectedvirtual

Implements Nektar::SolverUtils::Forcing.

Definition at line 64 of file ForcingAbsorption.cpp.

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

Referenced by CalcAbsorption(), and CalculateForcing().

◆ m_bRegions

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

Definition at line 90 of file ForcingAbsorption.h.

Referenced by CalcAbsorption().

◆ m_funcNameTime

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

Definition at line 89 of file ForcingAbsorption.h.

Referenced by CalculateForcing(), and v_InitObject().

◆ m_hasRefFlow

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

Definition at line 85 of file ForcingAbsorption.h.

Referenced by CalculateForcing(), and v_InitObject().

◆ m_hasRefFlowTime

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

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

Referenced by CalcAbsorption().