Nektar++
Static Public Member Functions | Static Public Attributes | Protected Member Functions | Private Member Functions | Friends | List of all members
Nektar::SymmetryBC Class Reference

Symmetry boundary conditions for compressible flow problems. More...

#include <SymmetryBC.h>

Inheritance diagram for Nektar::SymmetryBC:
[legend]

Static Public Member Functions

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. More...
 

Static Public Attributes

static std::string className
 Name of the class. More...
 

Protected Member Functions

void v_Apply (Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray, const NekDouble &time) override
 
- Protected Member Functions inherited from Nektar::CFSBndCond
 CFSBndCond (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)
 Constructor. More...
 
virtual void v_Apply (Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray, const NekDouble &time)=0
 
virtual void v_ApplyBwdWeight ()
 

Private Member Functions

 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)
 
 ~SymmetryBC (void) override
 

Friends

class MemoryManager< SymmetryBC >
 

Additional Inherited Members

- Public Member Functions inherited from Nektar::CFSBndCond
virtual ~CFSBndCond ()
 
void Apply (Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray, const NekDouble &time=0)
 Apply the boundary condition. More...
 
void ApplyBwdWeight ()
 Apply the Weight of boundary condition. More...
 
- Protected Attributes inherited from Nektar::CFSBndCond
LibUtilities::SessionReaderSharedPtr m_session
 Session reader. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 Array of fields. More...
 
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
 Trace normals. More...
 
Array< OneD, Array< OneD, NekDouble > > m_gridVelocityTrace
 Grid Velocity. More...
 
int m_spacedim
 Space dimension. More...
 
VariableConverterSharedPtr m_varConv
 Auxiliary object to convert variables. More...
 
NekDouble m_diffusionAveWeight
 Weight for average calculation of diffusion term. More...
 
NekDouble m_gamma
 Parameters of the flow. More...
 
NekDouble m_rhoInf
 
NekDouble m_pInf
 
NekDouble m_pOut
 
Array< OneD, NekDoublem_velInf
 
NekDouble m_angVel
 
int m_bcRegion
 Id of the boundary region. More...
 
int m_offset
 Offset. More...
 

Detailed Description

Symmetry boundary conditions for compressible flow problems.

Definition at line 46 of file SymmetryBC.h.

Constructor & Destructor Documentation

◆ SymmetryBC()

Nektar::SymmetryBC::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 
)
private

Definition at line 46 of file SymmetryBC.cpp.

52 : CFSBndCond(pSession, pFields, pTraceNormals, pGridVelocity, pSpaceDim,
53 bcRegion, cnt)
54{
56}
NekDouble m_diffusionAveWeight
Weight for average calculation of diffusion term.
Definition: CFSBndCond.h:102
CFSBndCond(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)
Constructor.
Definition: CFSBndCond.cpp:47

References Nektar::CFSBndCond::m_diffusionAveWeight.

◆ ~SymmetryBC()

Nektar::SymmetryBC::~SymmetryBC ( void  )
inlineoverrideprivate

Definition at line 80 of file SymmetryBC.h.

80{};

Member Function Documentation

◆ create()

static CFSBndCondSharedPtr Nektar::SymmetryBC::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 
)
inlinestatic

Creates an instance of this class.

Definition at line 52 of file SymmetryBC.h.

58 {
60 pSession, pFields, pTraceNormals, pGridVelocity, pSpaceDim,
61 bcRegion, cnt);
62 return p;
63 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< CFSBndCond > CFSBndCondSharedPtr
A shared pointer to a boundary condition object.
Definition: CFSBndCond.h:51

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and CellMLToNektar.cellml_metadata::p.

◆ v_Apply()

void Nektar::SymmetryBC::v_Apply ( Array< OneD, Array< OneD, NekDouble > > &  Fwd,
Array< OneD, Array< OneD, NekDouble > > &  physarray,
const NekDouble time 
)
overrideprotectedvirtual

Implements Nektar::CFSBndCond.

Definition at line 58 of file SymmetryBC.cpp.

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}
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
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

References Nektar::CFSBndCond::m_bcRegion, Nektar::CFSBndCond::m_fields, Nektar::CFSBndCond::m_offset, Nektar::CFSBndCond::m_spacedim, Nektar::CFSBndCond::m_traceNormals, Vmath::Smul(), Vmath::Vcopy(), and Vmath::Vvtvp().

Friends And Related Function Documentation

◆ MemoryManager< SymmetryBC >

friend class MemoryManager< SymmetryBC >
friend

Definition at line 1 of file SymmetryBC.h.

Member Data Documentation

◆ className

std::string Nektar::SymmetryBC::className
static
Initial value:
=
"Symmetry", SymmetryBC::create, "Symmetry boundary condition.")
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
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
CFSBndCondFactory & GetCFSBndCondFactory()
Declaration of the boundary condition factory singleton.
Definition: CFSBndCond.cpp:41

Name of the class.

Definition at line 66 of file SymmetryBC.h.