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 <boost/core/ignore_unused.hpp>
36 
37 #include "WallBC.h"
38 
39 using namespace std;
40 
41 namespace Nektar
42 {
43 
44 std::string WallBC::className = GetCFSBndCondFactory().RegisterCreatorFunction(
45  "Wall", WallBC::create, "Slip wall boundary condition.");
46 
47 WallBC::WallBC(const LibUtilities::SessionReaderSharedPtr &pSession,
49  const Array<OneD, Array<OneD, NekDouble>> &pTraceNormals,
50  const int pSpaceDim, const int bcRegion, const int cnt)
51  : CFSBndCond(pSession, pFields, pTraceNormals, pSpaceDim, bcRegion, cnt)
52 {
54 }
55 
57  Array<OneD, Array<OneD, NekDouble>> &physarray,
58  const NekDouble &time)
59 {
60  boost::ignore_unused(time);
61 
62  int i;
63  int nVariables = physarray.size();
64 
65  const Array<OneD, const int> &traceBndMap = m_fields[0]->GetTraceBndMap();
66 
67  // Adjust the physical values of the trace to take
68  // user defined boundaries into account
69  int e, id1, id2, nBCEdgePts, eMax;
70 
71  eMax = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExpSize();
72 
73  for (e = 0; e < eMax; ++e)
74  {
75  nBCEdgePts = m_fields[0]
76  ->GetBndCondExpansions()[m_bcRegion]
77  ->GetExp(e)
78  ->GetTotPoints();
79  id1 =
80  m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetPhys_Offset(e);
81  id2 =
82  m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[m_offset + e]);
83 
84  // Boundary condition for epsilon term.
85  if (nVariables == m_spacedim + 3)
86  {
87  Vmath::Zero(nBCEdgePts, &Fwd[nVariables - 1][id2], 1);
88  }
89 
90  // For 2D/3D, define: v* = v - 2(v.n)n
91  Array<OneD, NekDouble> tmp(nBCEdgePts, 0.0);
92 
93  // Calculate (v.n)
94  for (i = 0; i < m_spacedim; ++i)
95  {
96  Vmath::Vvtvp(nBCEdgePts, &Fwd[1 + i][id2], 1,
97  &m_traceNormals[i][id2], 1, &tmp[0], 1, &tmp[0], 1);
98  }
99 
100  // Calculate 2.0(v.n)
101  Vmath::Smul(nBCEdgePts, -2.0, &tmp[0], 1, &tmp[0], 1);
102 
103  // Calculate v* = v - 2.0(v.n)n
104  for (i = 0; i < m_spacedim; ++i)
105  {
106  Vmath::Vvtvp(nBCEdgePts, &tmp[0], 1, &m_traceNormals[i][id2], 1,
107  &Fwd[1 + i][id2], 1, &Fwd[1 + i][id2], 1);
108  }
109 
110  // Copy boundary adjusted values into the boundary expansion
111  for (i = 0; i < nVariables; ++i)
112  {
113  Vmath::Vcopy(nBCEdgePts, &Fwd[i][id2], 1,
114  &(m_fields[i]
115  ->GetBndCondExpansions()[m_bcRegion]
116  ->UpdatePhys())[id1],
117  1);
118  }
119  }
120 }
121 
122 } // namespace Nektar
Encapsulates the user-defined boundary conditions for compressible flow solver.
Definition: CFSBndCond.h:70
NekDouble m_diffusionAveWeight
Weight for average calculation of diffusion term.
Definition: CFSBndCond.h:99
int m_spacedim
Space dimension.
Definition: CFSBndCond.h:95
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Trace normals.
Definition: CFSBndCond.h:93
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
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
virtual void v_Apply(Array< OneD, Array< OneD, NekDouble >> &Fwd, Array< OneD, Array< OneD, NekDouble >> &physarray, const NekDouble &time) override
Definition: WallBC.cpp:56
std::shared_ptr< SessionReader > SessionReaderSharedPtr
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
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.cpp:574
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.cpp:248
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255