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
36
37namespace Nektar
38{
39
40std::string ForcingQuasi1D::className =
42 "Quasi1D", ForcingQuasi1D::create, "Quasi-1D nozzle Forcing");
43
46 const std::weak_ptr<SolverUtils::EquationSystem> &pEquation)
47 : Forcing(pSession, pEquation)
48{
49}
50
53 const unsigned int &pNumForcingFields, const TiXmlElement *pForce)
54{
55 m_NumVariable = pNumForcingFields;
56 m_varConv =
58
59 ASSERTL0(pFields[0]->GetGraph()->GetSpaceDimension() == 1,
60 "ForcingQuasi1D requires a 1D problem.");
61
62 const TiXmlElement *funcNameElmt = pForce->FirstChildElement("AREAFCN");
63 if (!funcNameElmt)
64 {
65 ASSERTL0(funcNameElmt,
66 "Requires AREAFCN tag "
67 "specifying function name which prescribes nozzle area.");
68 }
69
70 std::string funcName = funcNameElmt->GetText();
71 ASSERTL0(m_session->DefinesFunction(funcName),
72 "Function '" + funcName + "' not defined.");
73
74 // Evaluate geometrical term -Ax/A for forcing
75 m_geomFactor = Array<OneD, NekDouble>(pFields[0]->GetTotPoints(), 0.0);
76 Array<OneD, NekDouble> tmp(pFields[0]->GetTotPoints(), 0.0);
77
78 std::string sFieldStr = m_session->GetVariable(0);
79 ASSERTL0(m_session->DefinesFunction(funcName, sFieldStr),
80 "Variable '" + sFieldStr + "' not defined.");
81 GetFunction(pFields, m_session, funcName, true)
82 ->Evaluate(sFieldStr, m_geomFactor, 0.0);
83
84 // Check if DADXFCN is defined
85 const TiXmlElement *dAFuncNameElmt = pForce->FirstChildElement("DADXFCN");
86 if (dAFuncNameElmt)
87 {
88 funcName = dAFuncNameElmt->GetText();
89 ASSERTL0(m_session->DefinesFunction(funcName),
90 "Function '" + funcName + "' not defined.");
91 ASSERTL0(m_session->DefinesFunction(funcName, sFieldStr),
92 "Variable '" + sFieldStr + "' not defined.");
93 GetFunction(pFields, m_session, funcName, true)
94 ->Evaluate(sFieldStr, tmp, 0.0);
95 }
96 else
97 {
98 // Numerically evaluate dA/dX
99 pFields[0]->PhysDeriv(MultiRegions::DirCartesianMap[0], m_geomFactor,
100 tmp);
101 }
102
103 Vmath::Vdiv(pFields[0]->GetTotPoints(), tmp, 1, m_geomFactor, 1,
104 m_geomFactor, 1);
105 Vmath::Neg(pFields[0]->GetTotPoints(), m_geomFactor, 1);
106
107 // Project m_geomFactor to solution space
108 Array<OneD, NekDouble> tmpCoeff(pFields[0]->GetNcoeffs(), 0.0);
109 pFields[0]->FwdTransLocalElmt(m_geomFactor, tmpCoeff);
110 pFields[0]->BwdTrans(tmpCoeff, m_geomFactor);
111
113 for (int i = 0; i < m_NumVariable; ++i)
114 {
115 m_Forcing[i] = Array<OneD, NekDouble>(pFields[0]->GetTotPoints(), 0.0);
116 }
117}
118
121 const Array<OneD, Array<OneD, NekDouble>> &inarray,
123 [[maybe_unused]] const NekDouble &time)
124{
125 int nPoints = pFields[0]->GetTotPoints();
126
127 // Get (E+p)
128 Array<OneD, NekDouble> tmp(nPoints, 0.0);
129 m_varConv->GetPressure(inarray, tmp);
130 Vmath::Vadd(nPoints, tmp, 1, inarray[2], 1, tmp, 1);
131
132 // F-rho = -Ax/A *rhou
133 Vmath::Vmul(nPoints, m_geomFactor, 1, inarray[1], 1, m_Forcing[0], 1);
134
135 // F-rhou = -Ax/A *rhou*u
136 Vmath::Vmul(nPoints, m_geomFactor, 1, inarray[1], 1, m_Forcing[1], 1);
137 Vmath::Vmul(nPoints, m_Forcing[1], 1, inarray[1], 1, m_Forcing[1], 1);
138 Vmath::Vdiv(nPoints, m_Forcing[1], 1, inarray[0], 1, m_Forcing[1], 1);
139
140 // F-E = -Ax/A *(E+p)*u
141 Vmath::Vmul(nPoints, m_geomFactor, 1, inarray[1], 1, m_Forcing[2], 1);
142 Vmath::Vmul(nPoints, m_Forcing[2], 1, tmp, 1, m_Forcing[2], 1);
143 Vmath::Vdiv(nPoints, m_Forcing[2], 1, inarray[0], 1, m_Forcing[2], 1);
144
145 // Apply forcing
146 for (int i = 0; i < m_NumVariable; i++)
147 {
148 Vmath::Vadd(nPoints, outarray[i], 1, m_Forcing[i], 1, outarray[i], 1);
149 }
150}
151
152} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
void v_InitObject(const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields, const TiXmlElement *pForce) override
VariableConverterSharedPtr m_varConv
Array< OneD, NekDouble > m_geomFactor
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.
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
ForcingQuasi1D(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< SolverUtils::EquationSystem > &pEquation)
static std::string className
Name of the class.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
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:71
int m_NumVariable
Number of variables.
Definition: Forcing.h:129
Array< OneD, Array< OneD, NekDouble > > m_Forcing
Evaluated forcing function.
Definition: Forcing.h:127
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:150
LibUtilities::SessionReaderSharedPtr m_session
Session reader.
Definition: Forcing.h:123
std::shared_ptr< SessionReader > SessionReaderSharedPtr
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:87
ForcingFactory & GetForcingFactory()
Declaration of the forcing factory singleton.
Definition: Forcing.cpp:41
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.hpp:72
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.hpp:292
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