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

virtual void v_InitObject (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields, const TiXmlElement *pForce)
 
virtual 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 &pSession, const std::weak_ptr< EquationSystem > &pEquation)
 Constructor. More...
 
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 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 48 of file ForcingQuasi1D.cpp.

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

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 51 of file ForcingQuasi1D.h.

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

57  {
60  AllocateSharedPtr(pSession, pEquation);
61  p->InitObject(pFields, pNumForcingFields, pForce);
62  return p;
63  }
SOLVER_UTILS_EXPORT typedef std::shared_ptr< Forcing > ForcingSharedPtr
A shared pointer to an EquationSystem object.
Definition: Forcing.h:52
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

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

Implements Nektar::SolverUtils::Forcing.

Definition at line 124 of file ForcingQuasi1D.cpp.

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

129 {
130  boost::ignore_unused(time);
131 
132  int nPoints = pFields[0]->GetTotPoints();
133 
134  // Get (E+p)
135  Array<OneD, NekDouble> tmp (nPoints, 0.0);
136  m_varConv->GetPressure(inarray, tmp);
137  Vmath::Vadd(nPoints, tmp, 1,
138  inarray[2], 1, tmp, 1);
139 
140  // F-rho = -Ax/A *rhou
141  Vmath::Vmul(nPoints, m_geomFactor, 1,
142  inarray[1], 1, m_Forcing[0], 1);
143 
144  // F-rhou = -Ax/A *rhou*u
145  Vmath::Vmul(nPoints, m_geomFactor, 1,
146  inarray[1], 1, m_Forcing[1], 1);
147  Vmath::Vmul(nPoints, m_Forcing[1], 1,
148  inarray[1], 1, m_Forcing[1], 1);
149  Vmath::Vdiv(nPoints, m_Forcing[1], 1,
150  inarray[0], 1, m_Forcing[1], 1);
151 
152  // F-E = -Ax/A *(E+p)*u
153  Vmath::Vmul(nPoints, m_geomFactor, 1,
154  inarray[1], 1, m_Forcing[2], 1);
155  Vmath::Vmul(nPoints, m_Forcing[2], 1,
156  tmp, 1, m_Forcing[2], 1);
157  Vmath::Vdiv(nPoints, m_Forcing[2], 1,
158  inarray[0], 1, m_Forcing[2], 1);
159 
160  // Apply forcing
161  for (int i = 0; i < m_NumVariable; i++)
162  {
163  Vmath::Vadd(nPoints, outarray[i], 1,
164  m_Forcing[i], 1, outarray[i], 1);
165  }
166 }
Array< OneD, Array< OneD, NekDouble > > m_Forcing
Evaluated forcing function.
Definition: Forcing.h:107
VariableConverterSharedPtr m_varConv
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.cpp:244
Array< OneD, NekDouble > m_geomFactor
int m_NumVariable
Number of variables.
Definition: Forcing.h:109
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:302
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:186

◆ v_InitObject()

void Nektar::ForcingQuasi1D::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 ForcingQuasi1D.cpp.

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

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

Friends And Related Function Documentation

◆ MemoryManager< ForcingQuasi1D >

friend class MemoryManager< ForcingQuasi1D >
friend

Definition at line 48 of file ForcingQuasi1D.h.

Member Data Documentation

◆ className

std::string Nektar::ForcingQuasi1D::className
static
Initial value:
RegisterCreatorFunction("Quasi1D",
"Quasi-1D nozzle Forcing")

Name of the class.

Definition at line 66 of file ForcingQuasi1D.h.

◆ m_geomFactor

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

Definition at line 86 of file ForcingQuasi1D.h.

Referenced by v_Apply(), and v_InitObject().

◆ m_varConv

VariableConverterSharedPtr Nektar::ForcingQuasi1D::m_varConv
private

Definition at line 87 of file ForcingQuasi1D.h.

Referenced by v_Apply(), and v_InitObject().