Nektar++
APE.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: APE.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2018 Kilian Lackhove
10// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
11// Department of Aeronautics, Imperial College London (UK), and Scientific
12// Computing and Imaging Institute, University of Utah (USA).
13//
14// Permission is hereby granted, free of charge, to any person obtaining a
15// copy of this software and associated documentation files (the "Software"),
16// to deal in the Software without restriction, including without limitation
17// the rights to use, copy, modify, merge, publish, distribute, sublicense,
18// and/or sell copies of the Software, and to permit persons to whom the
19// Software is furnished to do so, subject to the following conditions:
20//
21// The above copyright notice and this permission notice shall be included
22// in all copies or substantial portions of the Software.
23//
24// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30// DEALINGS IN THE SOFTWARE.
31//
32// Description: APE1/APE4 (Acoustic Perturbation Equations)
33//
34///////////////////////////////////////////////////////////////////////////////
35
36#include <iostream>
37
39
40using namespace std;
41
42namespace Nektar
43{
45 "APE", APE::create, "APE1/APE4 (Acoustic Perturbation Equations)");
46
49 : UnsteadySystem(pSession, pGraph), AcousticSystem(pSession, pGraph)
50{
51 m_ip = 0;
52 m_irho = -1;
53 m_iu = 1;
54
55 m_conservative = false;
56}
57
58/**
59 * @brief Initialization object for the APE class.
60 */
61void APE::v_InitObject(bool DeclareFields)
62{
63 AcousticSystem::v_InitObject(DeclareFields);
64
65 // Initialize basefield again
67 for (int i = 0; i < m_bf.size(); ++i)
68 {
70 }
71 GetFunction("Baseflow", m_fields[0], true)
72 ->Evaluate(m_bfNames, m_bf, m_time);
73
74 // Define the normal velocity fields
76 for (int i = 0; i < m_bfFwdBwd.size(); i++)
77 {
79 }
80
81 string riemName;
82 m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
83 if (boost::to_lower_copy(riemName) == "characteristics" ||
84 boost::to_lower_copy(riemName) == "apeupwind" ||
85 boost::to_lower_copy(riemName) == "upwind")
86 {
87 riemName = "APEUpwind";
88 }
89 if (boost::to_lower_copy(riemName) == "laxfriedrichs")
90 {
91 riemName = "APELaxFriedrichs";
92 }
94 riemName, m_session);
95 m_riemannSolver->SetVector("N", &APE::GetNormals, this);
96 m_riemannSolver->SetVector("basefieldFwdBwd", &APE::GetBasefieldFwdBwd,
97 this);
98 m_riemannSolver->SetAuxVec("vecLocs", &APE::GetVecLocs, this);
99
100 // Set up advection operator
101 string advName;
102 m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
105 m_advection->SetFluxVector(&APE::v_GetFluxVector, this);
106 m_advection->SetRiemannSolver(m_riemannSolver);
107 m_advection->InitObject(m_session, m_fields);
108
110 {
113 }
114 else
115 {
116 ASSERTL0(false, "Implicit APE not set up.");
117 }
118}
119
120/**
121 * @brief Destructor for APE class.
122 */
124{
125}
126
127/**
128 * @brief Return the flux vector for the APE equations.
129 *
130 * @param physfield Fields.
131 * @param flux Resulting flux. flux[eq][dir][pt]
132 */
134 const Array<OneD, Array<OneD, NekDouble>> &physfield,
136{
137 int nq = physfield[0].size();
138 Array<OneD, NekDouble> tmp1(nq);
139 Array<OneD, NekDouble> tmp2(nq);
140
141 ASSERTL1(flux[0].size() == m_spacedim,
142 "Dimension of flux array and velocity array do not match");
143
144 // F_{adv,p',j} = \bar{rho} \bar{c^2} u'_j + p' \bar{u}_j
145 for (int j = 0; j < m_spacedim; ++j)
146 {
147 Vmath::Zero(nq, flux[0][j], 1);
148
149 // construct \bar{rho} \bar{c^2} u'_j
150 Vmath::Vmul(nq, m_bf[0], 1, m_bf[1], 1, tmp1, 1);
151 Vmath::Vmul(nq, tmp1, 1, physfield[j + 1], 1, tmp1, 1);
152
153 // construct p' \bar{u}_j term
154 Vmath::Vmul(nq, physfield[0], 1, m_bf[j + 2], 1, tmp2, 1);
155
156 // add both terms
157 Vmath::Vadd(nq, tmp1, 1, tmp2, 1, flux[0][j], 1);
158 }
159
160 for (int i = 1; i < flux.size(); ++i)
161 {
162 ASSERTL1(flux[i].size() == m_spacedim,
163 "Dimension of flux array and velocity array do not match");
164
165 // F_{adv,u'_i,j} = (p'/ \bar{rho} + \bar{u}_k u'_k) \delta_{ij}
166 for (int j = 0; j < m_spacedim; ++j)
167 {
168 Vmath::Zero(nq, flux[i][j], 1);
169
170 if (i - 1 == j)
171 {
172 // contruct p'/ \bar{rho} term
173 Vmath::Vdiv(nq, physfield[0], 1, m_bf[1], 1, flux[i][j], 1);
174
175 // construct \bar{u}_k u'_k term
176 Vmath::Zero(nq, tmp1, 1);
177 for (int k = 0; k < m_spacedim; ++k)
178 {
179 Vmath::Vvtvp(nq, physfield[k + 1], 1, m_bf[k + 2], 1, tmp1,
180 1, tmp1, 1);
181 }
182
183 // add terms
184 Vmath::Vadd(nq, flux[i][j], 1, tmp1, 1, flux[i][j], 1);
185 }
186 }
187 }
188}
189
190/**
191 * @brief Outflow characteristic boundary conditions for compressible
192 * flow problems.
193 */
194void APE::v_RiemannInvariantBC(int bcRegion, int cnt,
197 Array<OneD, Array<OneD, NekDouble>> &physarray)
198{
199 int id1, id2, nBCEdgePts;
200 int nVariables = physarray.size();
201
202 const Array<OneD, const int> &traceBndMap = m_fields[0]->GetTraceBndMap();
203
204 int eMax = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
205
206 for (int e = 0; e < eMax; ++e)
207 {
208 nBCEdgePts = m_fields[0]
209 ->GetBndCondExpansions()[bcRegion]
210 ->GetExp(e)
211 ->GetTotPoints();
212 id1 = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetPhys_Offset(e);
213 id2 = m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[cnt + e]);
214
215 // Calculate (v.n)
216 Array<OneD, NekDouble> Vn(nBCEdgePts, 0.0);
217 for (int i = 0; i < m_spacedim; ++i)
218 {
219 Vmath::Vvtvp(nBCEdgePts, &Fwd[m_iu + i][id2], 1,
220 &m_traceNormals[i][id2], 1, &Vn[0], 1, &Vn[0], 1);
221 }
222
223 // Calculate (v0.n)
224 Array<OneD, NekDouble> Vn0(nBCEdgePts, 0.0);
225 for (int i = 0; i < m_spacedim; ++i)
226 {
227 Vmath::Vvtvp(nBCEdgePts, &BfFwd[2 + i][id2], 1,
228 &m_traceNormals[i][id2], 1, &Vn0[0], 1, &Vn0[0], 1);
229 }
230
231 for (int i = 0; i < nBCEdgePts; ++i)
232 {
233 NekDouble c = sqrt(BfFwd[0][id2 + i]);
234
235 // LODI
236 NekDouble h1, h2;
237
238 // outgoing
239 if (Vn0[i] - c > 0)
240 {
241 // u/2 - p/(2*rho0*sqr(c0sq))
242 h1 = Vn[i] / 2 -
243 Fwd[m_ip][id2 + i] / (2 * BfFwd[1][id2 + i] * c);
244 }
245 // incoming
246 else
247 {
248 h1 = 0.0;
249 }
250
251 // outgoing
252 if (Vn0[i] + c > 0)
253 {
254 // u/2 + p/(2*rho0*sqr(c0sq))
255 h2 = Vn[i] / 2 +
256 Fwd[m_ip][id2 + i] / (2 * BfFwd[1][id2 + i] * c);
257 }
258 // incoming
259 else
260 {
261 h2 = 0.0;
262 }
263
264 // compute primitive variables
265 // p = rho0*sqr(c0sq) * (h2 - h1)
266 Fwd[m_ip][id2 + i] = BfFwd[1][id2 + i] * c * (h2 - h1);
267 // u = h1 + h2
268 NekDouble VnNew = h1 + h2;
269
270 // adjust velocity pert. according to new value
271 for (int j = 0; j < m_spacedim; ++j)
272 {
273 Fwd[m_iu + j][id2 + i] =
274 Fwd[m_iu + j][id2 + i] +
275 (VnNew - Vn[i]) * m_traceNormals[j][id2 + i];
276 }
277 }
278
279 // Copy boundary adjusted values into the boundary expansion
280 for (int i = 0; i < nVariables; ++i)
281 {
282 Vmath::Vcopy(nBCEdgePts, &Fwd[i][id2], 1,
283 &(m_fields[i]
284 ->GetBndCondExpansions()[bcRegion]
285 ->UpdatePhys())[id1],
286 1);
287 }
288 }
289}
290
291} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
static std::string className
Name of class.
Definition: APE.h:62
APE(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises UnsteadySystem class members.
Definition: APE.cpp:47
void v_GetFluxVector(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux) override
Return the flux vector for the APE equations.
Definition: APE.cpp:133
static EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
Definition: APE.h:52
~APE() override
Destructor.
Definition: APE.cpp:123
void v_InitObject(bool DeclareFields=true) override
Initialization object for the APE class.
Definition: APE.cpp:61
void v_RiemannInvariantBC(int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &BfFwd, Array< OneD, Array< OneD, NekDouble > > &physarray) override
Outflow characteristic boundary conditions for compressible flow problems.
Definition: APE.cpp:194
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the right-hand side.
const Array< OneD, const Array< OneD, NekDouble > > & GetVecLocs()
Get the locations of the components of the directed fields within the fields array.
const Array< OneD, const Array< OneD, NekDouble > > & GetNormals()
Get the normal vectors.
std::vector< std::string > m_bfNames
bool m_conservative
we are dealing with a conservative formualtion
SolverUtils::AdvectionSharedPtr m_advection
int m_ip
indices of the fields
Array< OneD, Array< OneD, NekDouble > > m_bfFwdBwd
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the projection and call the method for imposing the boundary conditions in case of discontinu...
Array< OneD, Array< OneD, NekDouble > > m_bf
void v_InitObject(bool DeclareFields=true) override
Initialization object for the AcousticSystem class.
const Array< OneD, const Array< OneD, NekDouble > > & GetBasefieldFwdBwd()
Get the baseflow field.
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
int m_spacedim
Spatial dimension (>= expansion dim).
SOLVER_UTILS_EXPORT int GetTraceNpoints()
NekDouble m_time
Current time of simulation.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
Get a SessionFunction by name.
SOLVER_UTILS_EXPORT int GetTotPoints()
Base class for unsteady solvers.
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
bool m_explicitAdvection
Indicates if explicit or implicit treatment of advection is used.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:43
EquationSystemFactory & GetEquationSystemFactory()
RiemannSolverFactory & GetRiemannSolverFactory()
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
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 Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.hpp:366
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
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
STL namespace.
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294