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

#include <ForcingAxiSymmetric.h>

Inheritance diagram for Nektar::ForcingAxiSymmetric:
[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

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

Private Attributes

Array< OneD, NekDoublem_geomFactor
 
VariableConverterSharedPtr m_varConv
 

Friends

class MemoryManager< ForcingAxiSymmetric >
 

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

Constructor & Destructor Documentation

◆ ForcingAxiSymmetric()

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

Definition at line 48 of file ForcingAxiSymmetric.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::ForcingAxiSymmetric::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 ForcingAxiSymmetric.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::ForcingAxiSymmetric::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 103 of file ForcingAxiSymmetric.cpp.

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

108 {
109  boost::ignore_unused(time);
110 
111  int nPoints = pFields[0]->GetTotPoints();
112 
113  // Get (E+p)
114  Array<OneD, NekDouble> tmp (nPoints, 0.0);
115  m_varConv->GetPressure(inarray, tmp);
116  Vmath::Vadd(nPoints, tmp, 1,
117  inarray[m_NumVariable-1], 1, tmp, 1);
118 
119  // F-rho = -1/r *rhou
120  Vmath::Vmul(nPoints, m_geomFactor, 1,
121  inarray[1], 1, m_Forcing[0], 1);
122 
123  // F-rhou_r = -1/r *rhou_r * u_r and F-rhou_y = -1/r *rhou_y * u_r
124  for (int i = 1; i < 3; ++i)
125  {
126  Vmath::Vmul(nPoints, inarray[1], 1,
127  inarray[i], 1, m_Forcing[i], 1);
128  Vmath::Vdiv(nPoints, m_Forcing[i], 1,
129  inarray[0], 1, m_Forcing[i], 1);
130  Vmath::Vmul(nPoints, m_Forcing[i], 1,
131  m_geomFactor, 1, m_Forcing[i], 1);
132  }
133 
134  // F-E = -1/r *(E+p)*u
135  Vmath::Vmul(nPoints, inarray[1], 1,
136  tmp, 1, m_Forcing[m_NumVariable-1], 1);
137  Vmath::Vdiv(nPoints, m_Forcing[m_NumVariable-1], 1,
138  inarray[0], 1, m_Forcing[m_NumVariable-1], 1);
139  Vmath::Vmul(nPoints, m_Forcing[m_NumVariable-1], 1,
141 
142  // Swirl
143  if (m_NumVariable == 5)
144  {
145  // F-rhou_r -= (-1/r) * rho * u_theta * u_theta
146  Vmath::Vmul(nPoints, inarray[3], 1,
147  inarray[3], 1, tmp, 1);
148  Vmath::Vdiv(nPoints, tmp, 1,
149  inarray[0], 1, tmp, 1);
150  Vmath::Vmul(nPoints, tmp, 1,
151  m_geomFactor, 1, tmp, 1);
152  Vmath::Vsub(nPoints, m_Forcing[1], 1,
153  tmp, 1, m_Forcing[1], 1);
154 
155  // F-rhou_theta = 2 * (-1/r *rhou_theta * u_r)
156  Vmath::Vmul(nPoints, inarray[1], 1,
157  inarray[3], 1, m_Forcing[3], 1);
158  Vmath::Vdiv(nPoints, m_Forcing[3], 1,
159  inarray[0], 1, m_Forcing[3], 1);
160  Vmath::Vmul(nPoints, m_Forcing[3], 1,
161  m_geomFactor, 1, m_Forcing[3], 1);
162  Vmath::Smul(nPoints, 2.0,
163  m_Forcing[3], 1, m_Forcing[3], 1);
164  }
165 
166  // Apply forcing
167  for (int i = 0; i < m_NumVariable; i++)
168  {
169  Vmath::Vadd(nPoints, outarray[i], 1,
170  m_Forcing[i], 1, outarray[i], 1);
171  }
172 }
Array< OneD, Array< OneD, NekDouble > > m_Forcing
Evaluated forcing function.
Definition: Forcing.h:107
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
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:216
Array< OneD, NekDouble > m_geomFactor
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:346
int m_NumVariable
Number of variables.
Definition: Forcing.h:109
VariableConverterSharedPtr m_varConv
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::ForcingAxiSymmetric::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 ForcingAxiSymmetric.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::NekConstants::kNekZeroTol, Nektar::SolverUtils::Forcing::m_Forcing, m_geomFactor, Nektar::SolverUtils::Forcing::m_NumVariable, Nektar::SolverUtils::Forcing::m_session, and m_varConv.

59 {
60  boost::ignore_unused(pForce);
61 
62  int spacedim = pFields[0]->GetGraph()->GetSpaceDimension();
63  int nPoints = pFields[0]->GetTotPoints();
64 
65  m_NumVariable = pNumForcingFields;
67  m_session, spacedim);
68 
69  // Get coordinates
70  Array<OneD, Array<OneD, NekDouble> > coords(3);
71  for (int i = 0; i < 3; i++)
72  {
73  coords[i] = Array<OneD, NekDouble> (nPoints);
74  }
75  pFields[0]->GetCoords(coords[0], coords[1], coords[2]);
76 
77  // Calculate fac = -1/r if r!=0, fac = 0 if r == 0
78  m_geomFactor = Array<OneD, NekDouble> (nPoints);
79  for (int i = 0; i < nPoints; ++i)
80  {
81  if (coords[0][i] < NekConstants::kNekZeroTol)
82  {
83  m_geomFactor[i] = 0;
84  }
85  else
86  {
87  m_geomFactor[i] = -1.0/coords[0][i];
88  }
89  }
90 
91  // Project m_geomFactor to solution space
92  Array<OneD, NekDouble> tmpCoeff (pFields[0]->GetNcoeffs(), 0.0);
93  pFields[0]->FwdTrans_IterPerExp(m_geomFactor, tmpCoeff);
94  pFields[0]->BwdTrans(tmpCoeff, m_geomFactor);
95 
96  m_Forcing = Array<OneD, Array<OneD, NekDouble> > (m_NumVariable);
97  for (int i = 0; i < m_NumVariable; ++i)
98  {
99  m_Forcing[i] = Array<OneD, NekDouble> (pFields[0]->GetTotPoints(), 0.0);
100  }
101 }
Array< OneD, Array< OneD, NekDouble > > m_Forcing
Evaluated forcing function.
Definition: Forcing.h:107
static const NekDouble kNekZeroTol
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
LibUtilities::SessionReaderSharedPtr m_session
Session reader.
Definition: Forcing.h:103
Array< OneD, NekDouble > m_geomFactor
int m_NumVariable
Number of variables.
Definition: Forcing.h:109
VariableConverterSharedPtr m_varConv

Friends And Related Function Documentation

◆ MemoryManager< ForcingAxiSymmetric >

friend class MemoryManager< ForcingAxiSymmetric >
friend

Definition at line 48 of file ForcingAxiSymmetric.h.

Member Data Documentation

◆ className

std::string Nektar::ForcingAxiSymmetric::className
static
Initial value:
RegisterCreatorFunction("AxiSymmetric",
"Forcing for axi-symmetric flow (around x=0)")

Name of the class.

Definition at line 66 of file ForcingAxiSymmetric.h.

◆ m_geomFactor

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

Definition at line 86 of file ForcingAxiSymmetric.h.

Referenced by v_Apply(), and v_InitObject().

◆ m_varConv

VariableConverterSharedPtr Nektar::ForcingAxiSymmetric::m_varConv
private

Definition at line 87 of file ForcingAxiSymmetric.h.

Referenced by v_Apply(), and v_InitObject().