Nektar++
ForcingQuasi1D.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: ForcingQuasi1D.cpp
4 //
5 // For more information, please see: http://www.nektar.info
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: Forcing for quasi-1D nozzle flow.
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #include <boost/core/ignore_unused.hpp>
36 
38 
39 using namespace std;
40 
41 namespace Nektar
42 {
43 std::string ForcingQuasi1D::className = SolverUtils::GetForcingFactory().
44  RegisterCreatorFunction("Quasi1D",
45  ForcingQuasi1D::create,
46  "Quasi-1D nozzle Forcing");
47 
48 ForcingQuasi1D::ForcingQuasi1D(
50  const std::weak_ptr<SolverUtils::EquationSystem> &pEquation)
51  : Forcing(pSession, pEquation)
52 {
53 }
54 
57  const unsigned int& pNumForcingFields,
58  const TiXmlElement* pForce)
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 
118  for (int i = 0; i < m_NumVariable; ++i)
119  {
120  m_Forcing[i] = Array<OneD, NekDouble> (pFields[0]->GetTotPoints(), 0.0);
121  }
122 }
123 
126  const Array<OneD, Array<OneD, NekDouble> >& inarray,
127  Array<OneD, Array<OneD, NekDouble> >& outarray,
128  const NekDouble& time)
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 }
167 
168 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
VariableConverterSharedPtr m_varConv
Array< OneD, NekDouble > m_geomFactor
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)
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
Defines a forcing term to be explicitly applied.
Definition: Forcing.h:73
int m_NumVariable
Number of variables.
Definition: Forcing.h:124
Array< OneD, Array< OneD, NekDouble > > m_Forcing
Evaluated forcing function.
Definition: Forcing.h:122
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:202
LibUtilities::SessionReaderSharedPtr m_session
Session reader.
Definition: Forcing.h:118
std::shared_ptr< SessionReader > SessionReaderSharedPtr
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:90
ForcingFactory & GetForcingFactory()
Declaration of the forcing factory singleton.
Definition: Forcing.cpp:44
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
double NekDouble
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:192
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:461
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:322
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:257