Nektar++
Static Public Member Functions | Static Public Attributes | Protected Member Functions | Private Member Functions | Private Attributes | Friends | List of all members
Nektar::ForcingQuasi1D Class Reference

#include <ForcingQuasi1D.h>

Inheritance diagram for Nektar::ForcingQuasi1D:
[legend]

Static Public Member Functions

static SolverUtils::ForcingSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< SolverUtils::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 Member Functions

void v_InitObject (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields, const TiXmlElement *pForce) override
 
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
 
- 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)
 

Private Member Functions

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

Private Attributes

Array< OneD, NekDoublem_geomFactor
 
VariableConverterSharedPtr m_varConv
 

Friends

class MemoryManager< ForcingQuasi1D >
 

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

Detailed Description

Definition at line 44 of file ForcingQuasi1D.h.

Constructor & Destructor Documentation

◆ ForcingQuasi1D()

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

Definition at line 45 of file ForcingQuasi1D.cpp.

48 : Forcing(pSession, pEquation)
49{
50}
SOLVER_UTILS_EXPORT Forcing(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation)
Constructor.
Definition: Forcing.cpp:48

Member Function Documentation

◆ create()

static SolverUtils::ForcingSharedPtr Nektar::ForcingQuasi1D::create ( const LibUtilities::SessionReaderSharedPtr pSession,
const std::weak_ptr< SolverUtils::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 50 of file ForcingQuasi1D.h.

55 {
58 pEquation);
59 p->InitObject(pFields, pNumForcingFields, pForce);
60 return p;
61 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
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::ForcingQuasi1D::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 120 of file ForcingQuasi1D.cpp.

125{
126 int nPoints = pFields[0]->GetTotPoints();
127
128 // Get (E+p)
129 Array<OneD, NekDouble> tmp(nPoints, 0.0);
130 m_varConv->GetPressure(inarray, tmp);
131 Vmath::Vadd(nPoints, tmp, 1, inarray[2], 1, tmp, 1);
132
133 // F-rho = -Ax/A *rhou
134 Vmath::Vmul(nPoints, m_geomFactor, 1, inarray[1], 1, m_Forcing[0], 1);
135
136 // F-rhou = -Ax/A *rhou*u
137 Vmath::Vmul(nPoints, m_geomFactor, 1, inarray[1], 1, m_Forcing[1], 1);
138 Vmath::Vmul(nPoints, m_Forcing[1], 1, inarray[1], 1, m_Forcing[1], 1);
139 Vmath::Vdiv(nPoints, m_Forcing[1], 1, inarray[0], 1, m_Forcing[1], 1);
140
141 // F-E = -Ax/A *(E+p)*u
142 Vmath::Vmul(nPoints, m_geomFactor, 1, inarray[1], 1, m_Forcing[2], 1);
143 Vmath::Vmul(nPoints, m_Forcing[2], 1, tmp, 1, m_Forcing[2], 1);
144 Vmath::Vdiv(nPoints, m_Forcing[2], 1, inarray[0], 1, m_Forcing[2], 1);
145
146 // Apply forcing
147 for (int i = 0; i < m_NumVariable; i++)
148 {
149 Vmath::Vadd(nPoints, outarray[i], 1, m_Forcing[i], 1, outarray[i], 1);
150 }
151}
VariableConverterSharedPtr m_varConv
Array< OneD, NekDouble > m_geomFactor
int m_NumVariable
Number of variables.
Definition: Forcing.h:121
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 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
void Vdiv(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:126

References Nektar::SolverUtils::Forcing::m_Forcing, m_geomFactor, Nektar::SolverUtils::Forcing::m_NumVariable, m_varConv, Vmath::Vadd(), Vmath::Vdiv(), and Vmath::Vmul().

◆ v_InitObject()

void Nektar::ForcingQuasi1D::v_InitObject ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  pFields,
const unsigned int &  pNumForcingFields,
const TiXmlElement *  pForce 
)
overrideprotectedvirtual

Implements Nektar::SolverUtils::Forcing.

Definition at line 52 of file ForcingQuasi1D.cpp.

55{
56 m_NumVariable = pNumForcingFields;
57 m_varConv =
59
60 ASSERTL0(pFields[0]->GetGraph()->GetSpaceDimension() == 1,
61 "ForcingQuasi1D requires a 1D problem.");
62
63 const TiXmlElement *funcNameElmt = pForce->FirstChildElement("AREAFCN");
64 if (!funcNameElmt)
65 {
66 ASSERTL0(funcNameElmt,
67 "Requires AREAFCN tag "
68 "specifying function name which prescribes nozzle area.");
69 }
70
71 string funcName = funcNameElmt->GetText();
72 ASSERTL0(m_session->DefinesFunction(funcName),
73 "Function '" + funcName + "' not defined.");
74
75 // Evaluate geometrical term -Ax/A for forcing
76 m_geomFactor = Array<OneD, NekDouble>(pFields[0]->GetTotPoints(), 0.0);
77 Array<OneD, NekDouble> tmp(pFields[0]->GetTotPoints(), 0.0);
78
79 std::string sFieldStr = m_session->GetVariable(0);
80 ASSERTL0(m_session->DefinesFunction(funcName, sFieldStr),
81 "Variable '" + sFieldStr + "' not defined.");
82 GetFunction(pFields, m_session, funcName, true)
83 ->Evaluate(sFieldStr, m_geomFactor, 0.0);
84
85 // Check if DADXFCN is defined
86 const TiXmlElement *dAFuncNameElmt = pForce->FirstChildElement("DADXFCN");
87 if (dAFuncNameElmt)
88 {
89 funcName = dAFuncNameElmt->GetText();
90 ASSERTL0(m_session->DefinesFunction(funcName),
91 "Function '" + funcName + "' not defined.");
92 ASSERTL0(m_session->DefinesFunction(funcName, sFieldStr),
93 "Variable '" + sFieldStr + "' not defined.");
94 GetFunction(pFields, m_session, funcName, true)
95 ->Evaluate(sFieldStr, tmp, 0.0);
96 }
97 else
98 {
99 // Numerically evaluate dA/dX
100 pFields[0]->PhysDeriv(MultiRegions::DirCartesianMap[0], m_geomFactor,
101 tmp);
102 }
103
104 Vmath::Vdiv(pFields[0]->GetTotPoints(), tmp, 1, m_geomFactor, 1,
105 m_geomFactor, 1);
106 Vmath::Neg(pFields[0]->GetTotPoints(), m_geomFactor, 1);
107
108 // Project m_geomFactor to solution space
109 Array<OneD, NekDouble> tmpCoeff(pFields[0]->GetNcoeffs(), 0.0);
110 pFields[0]->FwdTransLocalElmt(m_geomFactor, tmpCoeff);
111 pFields[0]->BwdTrans(tmpCoeff, m_geomFactor);
112
113 m_Forcing = Array<OneD, Array<OneD, NekDouble>>(m_NumVariable);
114 for (int i = 0; i < m_NumVariable; ++i)
115 {
116 m_Forcing[i] = Array<OneD, NekDouble>(pFields[0]->GetTotPoints(), 0.0);
117 }
118}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
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
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:86
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.hpp:292

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::MultiRegions::DirCartesianMap, Nektar::SolverUtils::Forcing::GetFunction(), Nektar::SolverUtils::Forcing::m_Forcing, m_geomFactor, Nektar::SolverUtils::Forcing::m_NumVariable, Nektar::SolverUtils::Forcing::m_session, m_varConv, Vmath::Neg(), and Vmath::Vdiv().

Friends And Related Function Documentation

◆ MemoryManager< ForcingQuasi1D >

friend class MemoryManager< ForcingQuasi1D >
friend

Definition at line 1 of file ForcingQuasi1D.h.

Member Data Documentation

◆ className

std::string Nektar::ForcingQuasi1D::className
static
Initial value:
=
"Quasi1D", ForcingQuasi1D::create, "Quasi-1D nozzle Forcing")
static SolverUtils::ForcingSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< SolverUtils::EquationSystem > &pEquation, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields, const TiXmlElement *pForce)
Creates an instance of this class.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
ForcingFactory & GetForcingFactory()
Declaration of the forcing factory singleton.
Definition: Forcing.cpp:42

Name of the class.

Definition at line 64 of file ForcingQuasi1D.h.

◆ m_geomFactor

Array<OneD, NekDouble> Nektar::ForcingQuasi1D::m_geomFactor
private

Definition at line 81 of file ForcingQuasi1D.h.

Referenced by v_Apply(), and v_InitObject().

◆ m_varConv

VariableConverterSharedPtr Nektar::ForcingQuasi1D::m_varConv
private

Definition at line 82 of file ForcingQuasi1D.h.

Referenced by v_Apply(), and v_InitObject().