Nektar++
WallRotationalBC.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: WallRotationalBC.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: No-slip wall boundary condition
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include <boost/core/ignore_unused.hpp>
36
37#include "WallRotationalBC.h"
38
39using namespace std;
40
41namespace Nektar
42{
43
46 "WallRotational", WallRotationalBC::create,
47 "Adiabatic rotational wall boundary condition.");
48
52 const Array<OneD, Array<OneD, NekDouble>> &pTraceNormals,
53 const Array<OneD, Array<OneD, NekDouble>> &pGridVelocity,
54 const int pSpaceDim, const int bcRegion, const int cnt)
55 : CFSBndCond(pSession, pFields, pTraceNormals, pGridVelocity, pSpaceDim,
56 bcRegion, cnt)
57{
59
60 // Set up rotational boundary edge velocities
61 const Array<OneD, const int> &traceBndMap = m_fields[0]->GetTraceBndMap();
62
63 int eMax = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExpSize();
64 for (int e = 0; e < eMax; ++e)
65 {
66 int nBCEdgePts = m_fields[0]
67 ->GetBndCondExpansions()[m_bcRegion]
68 ->GetExp(e)
69 ->GetTotPoints();
70 int id2 =
71 m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[m_offset + e]);
72
73 Array<OneD, NekDouble> x(nBCEdgePts, 0.0);
74 Array<OneD, NekDouble> y(nBCEdgePts, 0.0);
75 m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExp(e)->GetCoords(
76 x, y);
77 for (int j = 0; j < nBCEdgePts; ++j)
78 {
79 m_gridVelocityTrace[0][id2 + j] = m_angVel * -y[j];
80 m_gridVelocityTrace[1][id2 + j] = m_angVel * x[j];
81 }
82 }
83}
84
87 const NekDouble &time)
88{
89 boost::ignore_unused(time);
90 int nVariables = physarray.size();
91
92 // @TODO: ALE we subtract the grid velocity ? "Set u = to ug for this one" -
93 // Dave
94
95 int i;
96 const Array<OneD, const int> &traceBndMap = m_fields[0]->GetTraceBndMap();
97
98 // Take into account that for PDE based shock capturing, eps = 0 at the
99 // wall. Adjust the physical values of the trace to take user defined
100 // boundaries into account
101 int e, id1, id2, nBCEdgePts, eMax;
102
103 eMax = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExpSize();
104
105 for (e = 0; e < eMax; ++e)
106 {
107 nBCEdgePts = m_fields[0]
108 ->GetBndCondExpansions()[m_bcRegion]
109 ->GetExp(e)
110 ->GetTotPoints();
111 id1 =
112 m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetPhys_Offset(e);
113 id2 =
114 m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[m_offset + e]);
115
116 // Boundary condition for epsilon term. @TODO: Is this correct, or
117 // should I do E = p/(gamma -1) + 1/2*rho(u^2 +v^2 + w^2)... or
118 if (nVariables == m_spacedim + 3)
119 {
120 Vmath::Zero(nBCEdgePts, &Fwd[nVariables - 1][id2], 1);
121 }
122
123 for (i = 0; i < m_spacedim; i++)
124 {
125 // V = -Vin
126 // Vmath::Neg(nBCEdgePts, &Fwd[i+1][id2], 1);
127
128 // This now does Vg * rho + Vin
129 // Vmath::Vvtvp(nBCEdgePts, &m_gridVelocityTrace[i][id2], 1,
130 // &Fwd[0][id2], 1, &Fwd[i+1][id2], 1, &Fwd[i+1][id2], 1);
131
132 for (int j = 0; j < nBCEdgePts; ++j)
133 {
134 // Fwd[i+1][id2 + j] = 2 * m_gridVelocityTrace[i][id2 + j] *
135 // Fwd[0][id2 + j] - Fwd[i+1][id2 + j];
136 Fwd[i + 1][id2 + j] =
137 2 * m_gridVelocityTrace[i][id2 + j] * Fwd[0][id2 + j] -
138 Fwd[i + 1][id2 + j];
139 }
140 }
141
142 // Copy boundary adjusted values into the boundary expansion
143 for (i = 0; i < nVariables; ++i)
144 {
145 Vmath::Vcopy(nBCEdgePts, &Fwd[i][id2], 1,
146 &(m_fields[i]
147 ->GetBndCondExpansions()[m_bcRegion]
148 ->UpdatePhys())[id1],
149 1);
150 }
151 }
152}
153
154} // namespace Nektar
Encapsulates the user-defined boundary conditions for compressible flow solver.
Definition: CFSBndCond.h:71
NekDouble m_diffusionAveWeight
Weight for average calculation of diffusion term.
Definition: CFSBndCond.h:102
int m_spacedim
Space dimension.
Definition: CFSBndCond.h:98
Array< OneD, Array< OneD, NekDouble > > m_gridVelocityTrace
Grid Velocity.
Definition: CFSBndCond.h:96
int m_bcRegion
Id of the boundary region.
Definition: CFSBndCond.h:113
int m_offset
Offset.
Definition: CFSBndCond.h:115
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array of fields.
Definition: CFSBndCond.h:92
NekDouble m_angVel
Definition: CFSBndCond.h:110
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
void v_Apply(Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray, const NekDouble &time) override
static CFSBndCondSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble > > &pTraceNormals, const Array< OneD, Array< OneD, NekDouble > > &pGridVelocity, const int pSpaceDim, const int bcRegion, const int cnt)
Creates an instance of this class.
static std::string classNameRotational
Name of the class.
WallRotationalBC(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble > > &pTraceNormals, const Array< OneD, Array< OneD, NekDouble > > &pGridVelocity, const int pSpaceDim, const int bcRegion, const int cnt)
std::shared_ptr< SessionReader > SessionReaderSharedPtr
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
STL namespace.