Nektar++
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
Nektar::SolverUtils::Forcing Class Referenceabstract

Defines a forcing term to be explicitly applied. More...

#include <Forcing.h>

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

Public Member Functions

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

Static Public Member Functions

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)
 

Protected Member Functions

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

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

Detailed Description

Defines a forcing term to be explicitly applied.

Definition at line 70 of file Forcing.h.

Constructor & Destructor Documentation

◆ ~Forcing()

virtual SOLVER_UTILS_EXPORT Nektar::SolverUtils::Forcing::~Forcing ( )
inlinevirtual

Definition at line 73 of file Forcing.h.

74 {
75 }

◆ Forcing()

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

Constructor.

Definition at line 48 of file Forcing.cpp.

50 : m_session(pSession), m_equ(pEquation)
51{
52}
const std::weak_ptr< EquationSystem > m_equ
Weak pointer to equation system using this forcing.
Definition: Forcing.h:117
LibUtilities::SessionReaderSharedPtr m_session
Session reader.
Definition: Forcing.h:115

Member Function Documentation

◆ Apply()

void Nektar::SolverUtils::Forcing::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.

Parameters
fieldsExpansion lists corresponding to input arrays
inarrayu^n from previous timestep
outarrayoutput array to append forcing to

Definition at line 66 of file Forcing.cpp.

70{
71 v_Apply(fields, inarray, outarray, time);
72}
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

References v_Apply().

◆ ApplyCoeff()

void Nektar::SolverUtils::Forcing::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.

Parameters
fieldsExpansion lists corresponding to input arrays
inarrayu^n from previous timestep
outarrayoutput array to append forcing to

Definition at line 107 of file Forcing.cpp.

111{
112 v_ApplyCoeff(fields, inarray, outarray, time);
113}
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)
Definition: Forcing.cpp:217

References v_ApplyCoeff().

◆ EvaluateTimeFunction() [1/2]

void Nektar::SolverUtils::Forcing::EvaluateTimeFunction ( const NekDouble  pTime,
const LibUtilities::EquationSharedPtr pEqn,
Array< OneD, NekDouble > &  pArray 
)
protected

Definition at line 182 of file Forcing.cpp.

185{
186 // dummy array of zero pts.
187 Array<OneD, NekDouble> x0(pArray.size(), 0.0);
188
189 pEqn->Evaluate(x0, x0, x0, pTime, pArray);
190}

◆ EvaluateTimeFunction() [2/2]

void Nektar::SolverUtils::Forcing::EvaluateTimeFunction ( LibUtilities::SessionReaderSharedPtr  pSession,
std::string  pFieldName,
Array< OneD, NekDouble > &  pArray,
const std::string &  pFunctionName,
NekDouble  pTime = NekDouble(0) 
)
protected

Definition at line 168 of file Forcing.cpp.

172{
173 ASSERTL0(pSession->DefinesFunction(pFunctionName),
174 "Function '" + pFunctionName + "' does not exist.");
175
177 pSession->GetFunction(pFunctionName, pFieldName);
178
179 EvaluateTimeFunction(pTime, ffunc, pArray);
180}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
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:168
std::shared_ptr< Equation > EquationSharedPtr
Definition: Equation.h:125

References ASSERTL0, and EvaluateTimeFunction().

Referenced by EvaluateTimeFunction(), Nektar::SolverUtils::ForcingBody::v_Apply(), and Nektar::SolverUtils::ForcingBody::v_ApplyCoeff().

◆ GetForces()

const Nektar::Array< OneD, const Array< OneD, NekDouble > > & Nektar::SolverUtils::Forcing::GetForces ( )

Definition at line 163 of file Forcing.cpp.

164{
165 return m_Forcing;
166}
Array< OneD, Array< OneD, NekDouble > > m_Forcing
Evaluated forcing function.
Definition: Forcing.h:119

References m_Forcing.

◆ GetFunction()

SessionFunctionSharedPtr Nektar::SolverUtils::Forcing::GetFunction ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields,
const LibUtilities::SessionReaderSharedPtr pSession,
std::string  pName,
bool  pCache = false 
)
protected

Get a SessionFunction by name.

Definition at line 192 of file Forcing.cpp.

196{
197 if (pCache)
198 {
199 if ((m_sessionFunctions.find(pName) == m_sessionFunctions.end()) ||
200 (m_sessionFunctions[pName]->GetSession() != pSession) ||
201 (m_sessionFunctions[pName]->GetExpansion() != pFields[0]))
202 {
203 m_sessionFunctions[pName] =
205 pSession, pFields[0], pName, pCache);
206 }
207
208 return m_sessionFunctions[pName];
209 }
210 else
211 {
213 new SessionFunction(pSession, pFields[0], pName, pCache));
214 }
215}
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::map< std::string, SolverUtils::SessionFunctionSharedPtr > m_sessionFunctions
Map of known SessionFunctions.
Definition: Forcing.h:124
std::shared_ptr< SessionFunction > SessionFunctionSharedPtr

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

Referenced by Nektar::SolverUtils::ForcingAbsorption::CalcAbsorption(), Nektar::SolverUtils::ForcingAbsorption::CalculateForcing(), Nektar::SolverUtils::ForcingBody::Update(), Nektar::ForcingMovingBody::v_Apply(), Nektar::SolverUtils::ForcingAbsorption::v_InitObject(), and Nektar::ForcingQuasi1D::v_InitObject().

◆ InitObject()

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

Initialise the forcing object.

Definition at line 54 of file Forcing.cpp.

57{
58 v_InitObject(pFields, pNumForcingFields, pForce);
59}
virtual SOLVER_UTILS_EXPORT void v_InitObject(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields, const TiXmlElement *pForce)=0

References v_InitObject().

◆ Load()

vector< ForcingSharedPtr > Nektar::SolverUtils::Forcing::Load ( const LibUtilities::SessionReaderSharedPtr pSession,
const std::weak_ptr< EquationSystem > &  pEquation,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields,
const unsigned int &  pNumForcingFields = 0 
)
static

Definition at line 118 of file Forcing.cpp.

123{
124 vector<ForcingSharedPtr> vForceList;
125
126 if (!pSession->DefinesElement("Nektar/Forcing"))
127 {
128 return vForceList;
129 }
130
131 TiXmlElement *vForcing = pSession->GetElement("Nektar/Forcing");
132 if (vForcing)
133 {
134 unsigned int vNumForcingFields = pNumForcingFields;
135 if (!pNumForcingFields)
136 {
137 vNumForcingFields = pFields.size();
138 }
139
140 TiXmlElement *vForce = vForcing->FirstChildElement("FORCE");
141 while (vForce)
142 {
143 string vType = vForce->Attribute("TYPE");
144
145 TiXmlElement *vForceParam = vForce;
147 vForceParam, pSession->GetTimeLevel(), false);
148 vForceList.push_back(GetForcingFactory().CreateInstance(
149 vType, pSession, pEquation, pFields, vNumForcingFields,
150 vForceParam));
151
152 vForce = vForce->NextSiblingElement("FORCE");
153 }
154 }
155 return vForceList;
156}
static void GetXMLElementTimeLevel(TiXmlElement *&element, const size_t timeLevel, const bool enableCheck=true)
Get XML elment time level (Parallel-in-Time)
ForcingFactory & GetForcingFactory()
Declaration of the forcing factory singleton.
Definition: Forcing.cpp:42

References Nektar::SolverUtils::GetForcingFactory(), and Nektar::LibUtilities::SessionReader::GetXMLElementTimeLevel().

Referenced by Nektar::IncNavierStokes::v_InitObject(), Nektar::AcousticSystem::v_InitObject(), Nektar::UnsteadyAdvection::v_InitObject(), Nektar::UnsteadyInviscidBurgers::v_InitObject(), Nektar::UnsteadyReactionDiffusion::v_InitObject(), Nektar::CompressibleFlowSystem::v_InitObject(), and Nektar::Dummy::v_InitObject().

◆ PreApply()

void Nektar::SolverUtils::Forcing::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.

Definition at line 74 of file Forcing.cpp.

78{
79 v_PreApply(fields, inarray, outarray, time);
80}
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)
Definition: Forcing.cpp:82

References v_PreApply().

◆ UpdateForces()

Nektar::Array< OneD, Array< OneD, NekDouble > > & Nektar::SolverUtils::Forcing::UpdateForces ( )

Definition at line 158 of file Forcing.cpp.

159{
160 return m_Forcing;
161}

References m_Forcing.

◆ v_Apply()

virtual SOLVER_UTILS_EXPORT void Nektar::SolverUtils::Forcing::v_Apply ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble time 
)
protectedpure virtual

◆ v_ApplyCoeff()

void Nektar::SolverUtils::Forcing::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 in Nektar::SolverUtils::ForcingAbsorption, and Nektar::SolverUtils::ForcingBody.

Definition at line 217 of file Forcing.cpp.

222{
223 ASSERTL0(false, "v_ApplyCoeff not defined");
224}

References ASSERTL0.

Referenced by ApplyCoeff().

◆ v_InitObject()

virtual SOLVER_UTILS_EXPORT void Nektar::SolverUtils::Forcing::v_InitObject ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields,
const unsigned int &  pNumForcingFields,
const TiXmlElement *  pForce 
)
protectedpure virtual

◆ v_PreApply()

void Nektar::SolverUtils::Forcing::v_PreApply ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble time 
)
protectedvirtual

Reimplemented in Nektar::SolverUtils::ForcingMovingReferenceFrame.

Definition at line 82 of file Forcing.cpp.

87{
88 if (&inarray != &outarray)
89 {
90 int nvar = std::min(inarray.size(), outarray.size());
91 for (int i = 0; i < nvar; ++i)
92 {
93 if (&inarray[i] != &outarray[i])
94 {
95 int np = std::min(inarray[i].size(), outarray[i].size());
96 Vmath::Vcopy(np, inarray[i], 1, outarray[i], 1);
97 }
98 }
99 }
100}
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

References Vmath::Vcopy().

Referenced by PreApply().

Member Data Documentation

◆ m_equ

const std::weak_ptr<EquationSystem> Nektar::SolverUtils::Forcing::m_equ
protected

◆ m_Forcing

Array<OneD, Array<OneD, NekDouble> > Nektar::SolverUtils::Forcing::m_Forcing
protected

◆ m_NumVariable

int Nektar::SolverUtils::Forcing::m_NumVariable
protected

◆ m_session

LibUtilities::SessionReaderSharedPtr Nektar::SolverUtils::Forcing::m_session
protected

◆ m_sessionFunctions

std::map<std::string, SolverUtils::SessionFunctionSharedPtr> Nektar::SolverUtils::Forcing::m_sessionFunctions
protected

Map of known SessionFunctions.

Definition at line 124 of file Forcing.h.

Referenced by GetFunction().