Nektar++
WallBC.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: WallBC.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: Slip wall boundary condition
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include "WallBC.h"
36
37using namespace std;
38
39namespace Nektar
40{
41
43 "Wall", WallBC::create, "Slip wall boundary condition.");
44
47 const Array<OneD, Array<OneD, NekDouble>> &pTraceNormals,
48 const Array<OneD, Array<OneD, NekDouble>> &pGridVelocity,
49 const int pSpaceDim, const int bcRegion, const int cnt)
50 : CFSBndCond(pSession, pFields, pTraceNormals, pGridVelocity, pSpaceDim,
51 bcRegion, cnt)
52{
54}
55
58 [[maybe_unused]] const NekDouble &time)
59{
60 int i;
61 int nVariables = physarray.size();
62
63 const Array<OneD, const int> &traceBndMap = m_fields[0]->GetTraceBndMap();
64
65 // Adjust the physical values of the trace to take
66 // user defined boundaries into account
67 int e, id1, id2, nBCEdgePts, eMax;
68
69 eMax = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExpSize();
70
71 for (e = 0; e < eMax; ++e)
72 {
73 nBCEdgePts = m_fields[0]
74 ->GetBndCondExpansions()[m_bcRegion]
75 ->GetExp(e)
76 ->GetTotPoints();
77 id1 =
78 m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetPhys_Offset(e);
79 id2 =
80 m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[m_offset + e]);
81
82 // Boundary condition for epsilon term.
83 if (nVariables == m_spacedim + 3)
84 {
85 Vmath::Zero(nBCEdgePts, &Fwd[nVariables - 1][id2], 1);
86 }
87
88 // @TODO: Look at paper on this
89 // https://www.researchgate.net/publication/264044118_A_Guide_to_the_Implementation_of_Boundary_Conditions_in_Compact_High-Order_Methods_for_Compressible_Aerodynamics
90 // For 2D/3D, define: v* = v - 2(v.n)n
91 Array<OneD, NekDouble> tmp(nBCEdgePts, 0.0);
92
93 //@TODO: v - vg here... check nguyen paper, only issue is getting the vg
94 // for the trace in here
95 //@TODO: Update m_traceNormals, might be fine though.
96
97 if (m_fields[0]->GetGraph()->GetMovement()->GetMoveFlag())
98 {
99 for (i = 0; i < m_spacedim; ++i)
100 {
101 // This now does Vg * rho + Vin
102 for (int j = 0; j < nBCEdgePts; ++j)
103 {
104 Fwd[i + 1][id2 + j] +=
105 m_gridVelocityTrace[i][id2 + j] * Fwd[0][id2 + j];
106 }
107 }
108 }
109
110 // Calculate (v.n)
111 for (i = 0; i < m_spacedim; ++i)
112 {
113 Vmath::Vvtvp(nBCEdgePts, &Fwd[1 + i][id2], 1,
114 &m_traceNormals[i][id2], 1, &tmp[0], 1, &tmp[0], 1);
115 }
116
117 // Calculate 2.0(v.n)
118 Vmath::Smul(nBCEdgePts, -2.0, &tmp[0], 1, &tmp[0], 1);
119
120 // Calculate v* = v - 2.0(v.n)n
121 for (i = 0; i < m_spacedim; ++i)
122 {
123 Vmath::Vvtvp(nBCEdgePts, &tmp[0], 1, &m_traceNormals[i][id2], 1,
124 &Fwd[1 + i][id2], 1, &Fwd[1 + i][id2], 1);
125 }
126
127 // Copy boundary adjusted values into the boundary expansion
128 for (i = 0; i < nVariables; ++i)
129 {
130 Vmath::Vcopy(nBCEdgePts, &Fwd[i][id2], 1,
131 &(m_fields[i]
132 ->GetBndCondExpansions()[m_bcRegion]
133 ->UpdatePhys())[id1],
134 1);
135 }
136 }
137}
138
139} // 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
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Trace normals.
Definition: CFSBndCond.h:94
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
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
Definition: WallBC.cpp:56
static std::string className
Name of the class.
Definition: WallBC.h:66
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.
Definition: WallBC.h:52
WallBC(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)
Definition: WallBC.cpp:45
std::shared_ptr< SessionReader > SessionReaderSharedPtr
CFSBndCondFactory & GetCFSBndCondFactory()
Declaration of the boundary condition factory singleton.
Definition: CFSBndCond.cpp:41
double NekDouble
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 Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100
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.