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().
45  RegisterCreatorFunction("Wall",
46  WallBC::create,
47  "Slip wall boundary condition.");
48 
49 WallBC::WallBC(const LibUtilities::SessionReaderSharedPtr& pSession,
51  const Array<OneD, Array<OneD, NekDouble> >& pTraceNormals,
52  const int pSpaceDim,
53  const int bcRegion,
54  const int cnt)
55  : CFSBndCond(pSession, pFields, pTraceNormals, pSpaceDim, bcRegion, cnt)
56 {
57 }
58 
61  Array<OneD, Array<OneD, NekDouble> > &physarray,
62  const NekDouble &time)
63 {
64  boost::ignore_unused(time);
65 
66  int i;
67  int nVariables = physarray.num_elements();
68 
69  const Array<OneD, const int> &traceBndMap
70  = m_fields[0]->GetTraceBndMap();
71 
72  // Adjust the physical values of the trace to take
73  // user defined boundaries into account
74  int e, id1, id2, nBCEdgePts, eMax;
75 
76  eMax = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExpSize();
77 
78  for (e = 0; e < eMax; ++e)
79  {
80  nBCEdgePts = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->
81  GetExp(e)->GetTotPoints();
82  id1 = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->
83  GetPhys_Offset(e);
84  id2 = m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[m_offset+e]);
85 
86  // Boundary condition for epsilon term.
87  if (nVariables == m_spacedim+3)
88  {
89  Vmath::Zero(nBCEdgePts, &Fwd[nVariables-1][id2], 1);
90  }
91 
92  // For 2D/3D, define: v* = v - 2(v.n)n
93  Array<OneD, NekDouble> tmp(nBCEdgePts, 0.0);
94 
95  // Calculate (v.n)
96  for (i = 0; i < m_spacedim; ++i)
97  {
98  Vmath::Vvtvp(nBCEdgePts,
99  &Fwd[1+i][id2], 1,
100  &m_traceNormals[i][id2], 1,
101  &tmp[0], 1,
102  &tmp[0], 1);
103  }
104 
105  // Calculate 2.0(v.n)
106  Vmath::Smul(nBCEdgePts, -2.0, &tmp[0], 1, &tmp[0], 1);
107 
108  // Calculate v* = v - 2.0(v.n)n
109  for (i = 0; i < m_spacedim; ++i)
110  {
111  Vmath::Vvtvp(nBCEdgePts,
112  &tmp[0], 1,
113  &m_traceNormals[i][id2], 1,
114  &Fwd[1+i][id2], 1,
115  &Fwd[1+i][id2], 1);
116  }
117 
118  // Copy boundary adjusted values into the boundary expansion
119  for (i = 0; i < nVariables; ++i)
120  {
121  Vmath::Vcopy(nBCEdgePts, &Fwd[i][id2], 1,
122  &(m_fields[i]->GetBndCondExpansions()[m_bcRegion]->
123  UpdatePhys())[id1], 1);
124  }
125  }
126 }
127 
128 }
int m_spacedim
Space dimension.
Definition: CFSBndCond.h:89
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:445
STL namespace.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array of fields.
Definition: CFSBndCond.h:85
CFSBndCondFactory & GetCFSBndCondFactory()
Declaration of the boundary condition factory singleton.
Definition: CFSBndCond.cpp:41
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:216
Encapsulates the user-defined boundary conditions for compressible flow solver.
Definition: CFSBndCond.h:70
double NekDouble
virtual void v_Apply(Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray, const NekDouble &time)
Definition: WallBC.cpp:59
int m_offset
Offset.
Definition: CFSBndCond.h:102
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:376
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Trace normals.
Definition: CFSBndCond.h:87
std::shared_ptr< SessionReader > SessionReaderSharedPtr
int m_bcRegion
Id of the boundary region.
Definition: CFSBndCond.h:100