Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Static Public Member Functions | Static Public Attributes | 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:
Inheritance graph
[legend]
Collaboration diagram for Nektar::SolverUtils::ForcingAbsorption:
Collaboration graph
[legend]

Static Public Member Functions

static SOLVER_UTILS_EXPORT
ForcingSharedPtr 
create (const LibUtilities::SessionReaderSharedPtr &pSession, 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< ForcingSharedPtr
Load (const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields=0)
 

Static Public Attributes

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

Protected Member Functions

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)
 
- Protected Member Functions inherited from Nektar::SolverUtils::Forcing
SOLVER_UTILS_EXPORT Forcing (const LibUtilities::SessionReaderSharedPtr &)
 Constructor. More...
 
SOLVER_UTILS_EXPORT void EvaluateFunction (Array< OneD, MultiRegions::ExpListSharedPtr > pFields, LibUtilities::SessionReaderSharedPtr pSession, std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, NekDouble pTime=NekDouble(0))
 
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
 
- Protected Attributes inherited from Nektar::SolverUtils::Forcing
LibUtilities::SessionReaderSharedPtr m_session
 Session reader. More...
 
Array< OneD, Array< OneD,
NekDouble > > 
m_Forcing
 Evaluated forcing function. More...
 
int m_NumVariable
 Number of variables. More...
 

Private Member Functions

 ForcingAbsorption (const LibUtilities::SessionReaderSharedPtr &pSession)
 
virtual ~ForcingAbsorption (void)
 

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...
 

Detailed Description

Definition at line 51 of file ForcingAbsorption.h.

Constructor & Destructor Documentation

Nektar::SolverUtils::ForcingAbsorption::ForcingAbsorption ( const LibUtilities::SessionReaderSharedPtr pSession)
private

Definition at line 50 of file ForcingAbsorption.cpp.

51  : Forcing(pSession),
52  m_hasRefFlow(false),
53  m_hasRefFlowTime(false)
54  {
55  }
SOLVER_UTILS_EXPORT Forcing(const LibUtilities::SessionReaderSharedPtr &)
Constructor.
Definition: Forcing.cpp:53
virtual Nektar::SolverUtils::ForcingAbsorption::~ForcingAbsorption ( void  )
inlineprivatevirtual

Definition at line 92 of file ForcingAbsorption.h.

92 {};

Member Function Documentation

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

Creates an instance of this class.

Definition at line 57 of file ForcingAbsorption.h.

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

62  {
64  AllocateSharedPtr(pSession);
65  p->InitObject(pFields, pNumForcingFields, pForce);
66  return p;
67  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
SOLVER_UTILS_EXPORT typedef boost::shared_ptr< Forcing > ForcingSharedPtr
A shared pointer to an EquationSystem object.
Definition: Forcing.h:51
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 117 of file ForcingAbsorption.cpp.

References Nektar::SolverUtils::Forcing::EvaluateTimeFunction(), 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::Smul(), Vmath::Vadd(), Vmath::Vcopy(), Vmath::Vmul(), and Vmath::Vsub().

122  {
123  int nq = m_Forcing[0].num_elements();
124 
125  std::string s_FieldStr;
126  Array<OneD, NekDouble> TimeScale(1);
127  Array<OneD, Array<OneD, NekDouble> > RefflowScaled(m_NumVariable);
128 
129  if (m_hasRefFlow)
130  {
131  for (int i = 0; i < m_NumVariable; i++)
132  {
133  RefflowScaled[i] = Array<OneD, NekDouble> (nq);
134  if (m_hasRefFlowTime)
135  {
136  s_FieldStr = m_session->GetVariable(i);
137  EvaluateTimeFunction(m_session, s_FieldStr, TimeScale, m_funcNameTime, time);
138  Vmath::Smul(nq, TimeScale[0], m_Refflow[i],1,RefflowScaled[i],1);
139  }
140  else
141  {
142  Vmath::Vcopy(nq, m_Refflow[i],1, RefflowScaled[i],1);
143  }
144 
145 
146  Vmath::Vsub(nq, inarray[i], 1,
147  RefflowScaled[i], 1, m_Forcing[i], 1);
148  Vmath::Vmul(nq, m_Absorption[i], 1,
149  m_Forcing[i], 1, m_Forcing[i], 1);
150  Vmath::Vadd(nq, m_Forcing[i], 1,
151  outarray[i], 1, outarray[i], 1);
152  }
153  }
154  else
155  {
156  for (int i = 0; i < m_NumVariable; i++)
157  {
158  Vmath::Vmul(nq, m_Absorption[i], 1,
159  inarray[i], 1, m_Forcing[i], 1);
160  Vmath::Vadd(nq, m_Forcing[i], 1,
161  outarray[i], 1, outarray[i], 1);
162  }
163  }
164  }
SOLVER_UTILS_EXPORT void EvaluateTimeFunction(LibUtilities::SessionReaderSharedPtr pSession, std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, NekDouble pTime=NekDouble(0))
Definition: Forcing.cpp:121
Array< OneD, Array< OneD, NekDouble > > m_Forcing
Evaluated forcing function.
Definition: Forcing.h:97
Array< OneD, Array< OneD, NekDouble > > m_Absorption
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
Array< OneD, Array< OneD, NekDouble > > m_Refflow
LibUtilities::SessionReaderSharedPtr m_session
Session reader.
Definition: Forcing.h:95
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:329
int m_NumVariable
Number of variables.
Definition: Forcing.h:99
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
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:285
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:169
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 57 of file ForcingAbsorption.cpp.

References ASSERTL0, Nektar::SolverUtils::Forcing::EvaluateFunction(), 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, and npts.

61  {
62  m_NumVariable = pNumForcingFields;
63  int npts = pFields[0]->GetTotPoints();
64 
65  const TiXmlElement* funcNameElmt;
66  funcNameElmt = pForce->FirstChildElement("COEFF");
67  ASSERTL0(funcNameElmt, "Requires COEFF tag, specifying function "
68  "name which prescribes absorption layer coefficient.");
69 
70  string funcName = funcNameElmt->GetText();
71  ASSERTL0(m_session->DefinesFunction(funcName),
72  "Function '" + funcName + "' not defined.");
73 
74  std::string s_FieldStr;
75  m_Absorption = Array<OneD, Array<OneD, NekDouble> > (m_NumVariable);
76  m_Forcing = Array<OneD, Array<OneD, NekDouble> > (m_NumVariable);
77  for (int i = 0; i < m_NumVariable; ++i)
78  {
79  s_FieldStr = m_session->GetVariable(i);
80  ASSERTL0(m_session->DefinesFunction(funcName, s_FieldStr),
81  "Variable '" + s_FieldStr + "' not defined.");
82  m_Absorption[i] = Array<OneD, NekDouble> (npts, 0.0);
83  m_Forcing[i] = Array<OneD, NekDouble> (npts, 0.0);
84  EvaluateFunction(pFields, m_session, s_FieldStr,
85  m_Absorption[i], funcName);
86  }
87 
88  funcNameElmt = pForce->FirstChildElement("REFFLOW");
89  if (funcNameElmt)
90  {
91  string funcName = funcNameElmt->GetText();
92  ASSERTL0(m_session->DefinesFunction(funcName),
93  "Function '" + funcName + "' not defined.");
94 
95  m_Refflow = Array<OneD, Array<OneD, NekDouble> > (m_NumVariable);
96  for (int i = 0; i < m_NumVariable; ++i)
97  {
98  s_FieldStr = m_session->GetVariable(i);
99  ASSERTL0(m_session->DefinesFunction(funcName, s_FieldStr),
100  "Variable '" + s_FieldStr + "' not defined.");
101  m_Refflow[i] = Array<OneD, NekDouble> (npts, 0.0);
102  EvaluateFunction(pFields, m_session, s_FieldStr,
103  m_Refflow[i], funcName);
104  }
105  m_hasRefFlow = true;
106  }
107 
108  funcNameElmt = pForce->FirstChildElement("REFFLOWTIME");
109  if (funcNameElmt)
110  {
111  m_funcNameTime = funcNameElmt->GetText();
112  m_hasRefFlowTime = true;
113  }
114 
115  }
Array< OneD, Array< OneD, NekDouble > > m_Forcing
Evaluated forcing function.
Definition: Forcing.h:97
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
SOLVER_UTILS_EXPORT void EvaluateFunction(Array< OneD, MultiRegions::ExpListSharedPtr > pFields, LibUtilities::SessionReaderSharedPtr pSession, std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, NekDouble pTime=NekDouble(0))
Definition: Forcing.cpp:151
Array< OneD, Array< OneD, NekDouble > > m_Absorption
Array< OneD, Array< OneD, NekDouble > > m_Refflow
static std::string npts
Definition: InputFld.cpp:43
LibUtilities::SessionReaderSharedPtr m_session
Session reader.
Definition: Forcing.h:95
int m_NumVariable
Number of variables.
Definition: Forcing.h:99

Friends And Related Function Documentation

friend class MemoryManager< ForcingAbsorption >
friend

Definition at line 54 of file ForcingAbsorption.h.

Member Data Documentation

std::string Nektar::SolverUtils::ForcingAbsorption::className
static
Initial value:
RegisterCreatorFunction("Absorption",
"Forcing Absorption")

Name of the class.

Definition at line 70 of file ForcingAbsorption.h.

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

Definition at line 75 of file ForcingAbsorption.h.

Referenced by v_Apply(), and v_InitObject().

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

Definition at line 77 of file ForcingAbsorption.h.

Referenced by v_Apply(), and v_InitObject().

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

Definition at line 73 of file ForcingAbsorption.h.

Referenced by v_Apply(), and v_InitObject().

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

Definition at line 74 of file ForcingAbsorption.h.

Referenced by v_Apply(), and v_InitObject().

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

Definition at line 76 of file ForcingAbsorption.h.

Referenced by v_Apply(), and v_InitObject().