Nektar++
SymmetryBC.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: SymmetryBC.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: Symmetry boundary condition
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include "SymmetryBC.h"
36
37using namespace std;
38
39namespace Nektar
40{
41
42std::string SymmetryBC::className =
44 "Symmetry", SymmetryBC::create, "Symmetry boundary condition.");
45
49 const Array<OneD, Array<OneD, NekDouble>> &pTraceNormals,
50 const Array<OneD, Array<OneD, NekDouble>> &pGridVelocity,
51 const int pSpaceDim, const int bcRegion, const int cnt)
52 : CFSBndCond(pSession, pFields, pTraceNormals, pGridVelocity, pSpaceDim,
53 bcRegion, cnt)
54{
56}
57
60 [[maybe_unused]] const NekDouble &time)
61{
62 int i;
63 int nVariables = physarray.size();
64
65 const Array<OneD, const int> &traceBndMap = m_fields[0]->GetTraceBndMap();
66
67 // Take into account that for PDE based shock capturing, eps = 0 at the
68 // wall.
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 // For 2D/3D, define: v* = v - 2(v.n)n
85 Array<OneD, NekDouble> tmp(nBCEdgePts, 0.0);
86
87 // Calculate (v.n)
88 for (i = 0; i < m_spacedim; ++i)
89 {
90 Vmath::Vvtvp(nBCEdgePts, &Fwd[1 + i][id2], 1,
91 &m_traceNormals[i][id2], 1, &tmp[0], 1, &tmp[0], 1);
92 }
93
94 // Calculate 2.0(v.n)
95 Vmath::Smul(nBCEdgePts, -2.0, &tmp[0], 1, &tmp[0], 1);
96
97 // Calculate v* = v - 2.0(v.n)n
98 for (i = 0; i < m_spacedim; ++i)
99 {
100 Vmath::Vvtvp(nBCEdgePts, &tmp[0], 1, &m_traceNormals[i][id2], 1,
101 &Fwd[1 + i][id2], 1, &Fwd[1 + i][id2], 1);
102 }
103
104 // Copy boundary adjusted values into the boundary expansion
105 for (i = 0; i < nVariables; ++i)
106 {
107 Vmath::Vcopy(nBCEdgePts, &Fwd[i][id2], 1,
108 &(m_fields[i]
109 ->GetBndCondExpansions()[m_bcRegion]
110 ->UpdatePhys())[id1],
111 1);
112 }
113 }
114}
115
116} // 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_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: SymmetryBC.cpp:58
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: SymmetryBC.h:52
SymmetryBC(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: SymmetryBC.cpp:46
static std::string className
Name of the class.
Definition: SymmetryBC.h:66
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 Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
STL namespace.