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

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.
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))
void EvaluateTimeFunction (LibUtilities::SessionReaderSharedPtr pSession, std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, NekDouble pTime=NekDouble(0))

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.
Array< OneD, Array< OneD,
NekDouble > > 
m_Forcing
 Evaluated forcing function.
int m_NumVariable
 Number of variables.

Private Member Functions

 ForcingAbsorption (const LibUtilities::SessionReaderSharedPtr &pSession)

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

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

: Forcing(pSession),
m_hasRefFlow(false),
{
}

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::SolverUtils::ForcingSharedPtr.

{
p->InitObject(pFields, pNumForcingFields, pForce);
return p;
}
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().

{
int nq = m_Forcing[0].num_elements();
std::string s_FieldStr;
Array<OneD, NekDouble> TimeScale(1);
Array<OneD, Array<OneD, NekDouble> > RefflowScaled(m_NumVariable);
{
for (int i = 0; i < m_NumVariable; i++)
{
RefflowScaled[i] = Array<OneD, NekDouble> (nq);
{
s_FieldStr = m_session->GetVariable(i);
EvaluateTimeFunction(m_session, s_FieldStr, TimeScale, m_funcNameTime, time);
Vmath::Smul(nq, TimeScale[0], m_Refflow[i],1,RefflowScaled[i],1);
}
else
{
Vmath::Vcopy(nq, m_Refflow[i],1, RefflowScaled[i],1);
}
Vmath::Vsub(nq, inarray[i], 1,
RefflowScaled[i], 1, m_Forcing[i], 1);
m_Forcing[i], 1, m_Forcing[i], 1);
outarray[i], 1, outarray[i], 1);
}
}
else
{
for (int i = 0; i < m_NumVariable; i++)
{
inarray[i], 1, m_Forcing[i], 1);
outarray[i], 1, outarray[i], 1);
}
}
}
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 55 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.

{
m_NumVariable = pNumForcingFields;
int npts = pFields[0]->GetTotPoints();
const TiXmlElement* funcNameElmt;
funcNameElmt = pForce->FirstChildElement("COEFF");
ASSERTL0(funcNameElmt, "Requires COEFF tag, specifying function "
"name which prescribes absorption layer coefficient.");
string funcName = funcNameElmt->GetText();
ASSERTL0(m_session->DefinesFunction(funcName),
"Function '" + funcName + "' not defined.");
std::string s_FieldStr;
m_Absorption = Array<OneD, Array<OneD, NekDouble> > (m_NumVariable);
m_Forcing = Array<OneD, Array<OneD, NekDouble> > (m_NumVariable);
for (int i = 0; i < m_NumVariable; ++i)
{
s_FieldStr = m_session->GetVariable(i);
ASSERTL0(m_session->DefinesFunction(funcName, s_FieldStr),
"Variable '" + s_FieldStr + "' not defined.");
m_Absorption[i] = Array<OneD, NekDouble> (npts, 0.0);
m_Forcing[i] = Array<OneD, NekDouble> (npts, 0.0);
EvaluateFunction(pFields, m_session, s_FieldStr,
m_Absorption[i], funcName);
}
funcNameElmt = pForce->FirstChildElement("REFFLOW");
if (funcNameElmt)
{
string funcName = funcNameElmt->GetText();
ASSERTL0(m_session->DefinesFunction(funcName),
"Function '" + funcName + "' not defined.");
m_Refflow = Array<OneD, Array<OneD, NekDouble> > (m_NumVariable);
for (int i = 0; i < m_NumVariable; ++i)
{
s_FieldStr = m_session->GetVariable(i);
ASSERTL0(m_session->DefinesFunction(funcName, s_FieldStr),
"Variable '" + s_FieldStr + "' not defined.");
m_Refflow[i] = Array<OneD, NekDouble> (npts, 0.0);
EvaluateFunction(pFields, m_session, s_FieldStr,
m_Refflow[i], funcName);
}
m_hasRefFlow = true;
}
funcNameElmt = pForce->FirstChildElement("REFFLOWTIME");
if (funcNameElmt)
{
m_funcNameTime = funcNameElmt->GetText();
ASSERTL0(m_session->DefinesFunction(funcName),
"Function '" + funcName + "' not defined.");
}
}

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