Nektar++
IsentropicVortexBC.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: IsentropicVortexBC.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: Isentropic vortex boundary condition
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include "IsentropicVortexBC.h"
36
37using namespace std;
38
39namespace Nektar
40{
41
44 "IsentropicVortex", IsentropicVortexBC::create,
45 "Isentropic vortex boundary condition.");
46
50 const Array<OneD, Array<OneD, NekDouble>> &pTraceNormals,
51 const int pSpaceDim, const int bcRegion, const int cnt)
52 : CFSBndCond(pSession, pFields, pTraceNormals, pSpaceDim, bcRegion, cnt)
53{
54}
55
58 const NekDouble &time)
59{
60 int nvariables = physarray.size();
61
62 const Array<OneD, const int> &bndTraceMap = m_fields[0]->GetTraceBndMap();
63 // loop over Boundary Regions
64 int npoints, id1, id2, e_max;
65
66 e_max = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExpSize();
67
68 for (int e = 0; e < e_max; ++e)
69 {
70 npoints = m_fields[0]
71 ->GetBndCondExpansions()[m_bcRegion]
72 ->GetExp(e)
73 ->GetTotPoints();
74 id1 =
75 m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetPhys_Offset(e);
76 id2 =
77 m_fields[0]->GetTrace()->GetPhys_Offset(bndTraceMap[m_offset + e]);
78
79 Array<OneD, NekDouble> x(npoints, 0.0);
80 Array<OneD, NekDouble> y(npoints, 0.0);
81 Array<OneD, NekDouble> z(npoints, 0.0);
82
83 m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExp(e)->GetCoords(
84 x, y, z);
85
86 EvaluateIsentropicVortex(x, y, z, Fwd, time, id2);
87
88 for (int i = 0; i < nvariables; ++i)
89 {
90 Vmath::Vcopy(npoints, &Fwd[i][id2], 1,
91 &(m_fields[i]
92 ->GetBndCondExpansions()[m_bcRegion]
93 ->UpdatePhys())[id1],
94 1);
95 }
96 }
97}
98
101 [[maybe_unused]] const Array<OneD, NekDouble> &z,
102 Array<OneD, Array<OneD, NekDouble>> &u, NekDouble time, const int o)
103{
104 int nq = x.size();
105
106 // Flow parameters
107 const NekDouble x0 = 5.0;
108 const NekDouble y0 = 0.0;
109 const NekDouble beta = 5.0;
110 const NekDouble u0 = 1.0;
111 const NekDouble v0 = 0.5;
112 const NekDouble gamma = m_gamma;
113 NekDouble r, xbar, ybar, tmp;
114 NekDouble fac = 1.0 / (16.0 * gamma * M_PI * M_PI);
115
116 // In 3D zero rhow field.
117 if (m_spacedim == 3)
118 {
119 Vmath::Zero(nq, &u[3][o], 1);
120 }
121
122 // Fill storage
123 for (int i = 0; i < nq; ++i)
124 {
125 xbar = x[i] - u0 * time - x0;
126 ybar = y[i] - v0 * time - y0;
127 r = sqrt(xbar * xbar + ybar * ybar);
128 tmp = beta * exp(1 - r * r);
129 u[0][i + o] =
130 pow(1.0 - (gamma - 1.0) * tmp * tmp * fac, 1.0 / (gamma - 1.0));
131 u[1][i + o] = u[0][i + o] * (u0 - tmp * ybar / (2 * M_PI));
132 u[2][i + o] = u[0][i + o] * (v0 + tmp * xbar / (2 * M_PI));
133 u[m_spacedim + 1][i + o] =
134 pow(u[0][i + o], gamma) / (gamma - 1.0) +
135 0.5 * (u[1][i + o] * u[1][i + o] + u[2][i + o] * u[2][i + o]) /
136 u[0][i + o];
137 }
138}
139
140} // namespace Nektar
Encapsulates the user-defined boundary conditions for compressible flow solver.
Definition: CFSBndCond.h:70
int m_spacedim
Space dimension.
Definition: CFSBndCond.h:95
NekDouble m_gamma
Parameters of the flow.
Definition: CFSBndCond.h:102
int m_bcRegion
Id of the boundary region.
Definition: CFSBndCond.h:109
int m_offset
Offset.
Definition: CFSBndCond.h:111
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array of fields.
Definition: CFSBndCond.h:91
IsentropicVortexBC(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble > > &pTraceNormals, const int pSpaceDim, const int bcRegion, const int cnt)
static std::string className
Name of the class.
static CFSBndCondSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble > > &pTraceNormals, const int pSpaceDim, const int bcRegion, const int cnt)
Creates an instance of this class.
void EvaluateIsentropicVortex(const Array< OneD, NekDouble > &x, const Array< OneD, NekDouble > &y, const Array< OneD, NekDouble > &z, Array< OneD, Array< OneD, NekDouble > > &u, NekDouble time, const int o=0)
void v_Apply(Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray, const NekDouble &time) override
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
std::shared_ptr< SessionReader > SessionReaderSharedPtr
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:59
std::vector< double > z(NPUPPER)
CFSBndCondFactory & GetCFSBndCondFactory()
Declaration of the boundary condition factory singleton.
Definition: CFSBndCond.cpp:41
double NekDouble
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
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294