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

virtual 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 int pSpaceDim, const int bcRegion, const int cnt)
 Constructor. More...
 
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 int pSpaceDim, const int bcRegion, const int cnt)
 
virtual ~SymmetryBC (void)
 

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...
 
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
 
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 int  pSpaceDim,
const int  bcRegion,
const int  cnt 
)
private

Definition at line 48 of file SymmetryBC.cpp.

53  : CFSBndCond(pSession, pFields, pTraceNormals, pSpaceDim, bcRegion, cnt)
54 {
56 }
NekDouble m_diffusionAveWeight
Weight for average calculation of diffusion term.
Definition: CFSBndCond.h:99
CFSBndCond(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble >> &pTraceNormals, const int pSpaceDim, const int bcRegion, const int cnt)
Constructor.
Definition: CFSBndCond.cpp:47

References Nektar::CFSBndCond::m_diffusionAveWeight.

◆ ~SymmetryBC()

virtual Nektar::SymmetryBC::~SymmetryBC ( void  )
inlineprivatevirtual

Definition at line 77 of file SymmetryBC.h.

77 {};

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 int  pSpaceDim,
const int  bcRegion,
const int  cnt 
)
inlinestatic

Creates an instance of this class.

Definition at line 52 of file SymmetryBC.h.

57  {
59  pSession, pFields, pTraceNormals, pSpaceDim, bcRegion, cnt);
60  return p;
61  }
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:48

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  boost::ignore_unused(time);
63 
64  int i;
65  int nVariables = physarray.size();
66 
67  const Array<OneD, const int> &traceBndMap = m_fields[0]->GetTraceBndMap();
68 
69  // Take into account that for PDE based shock capturing, eps = 0 at the
70  // wall.
71  int e, id1, id2, nBCEdgePts, eMax;
72 
73  eMax = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExpSize();
74 
75  for (e = 0; e < eMax; ++e)
76  {
77  nBCEdgePts = m_fields[0]
78  ->GetBndCondExpansions()[m_bcRegion]
79  ->GetExp(e)
80  ->GetTotPoints();
81  id1 =
82  m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetPhys_Offset(e);
83  id2 =
84  m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[m_offset + e]);
85 
86  // For 2D/3D, define: v* = v - 2(v.n)n
87  Array<OneD, NekDouble> tmp(nBCEdgePts, 0.0);
88 
89  // Calculate (v.n)
90  for (i = 0; i < m_spacedim; ++i)
91  {
92  Vmath::Vvtvp(nBCEdgePts, &Fwd[1 + i][id2], 1,
93  &m_traceNormals[i][id2], 1, &tmp[0], 1, &tmp[0], 1);
94  }
95 
96  // Calculate 2.0(v.n)
97  Vmath::Smul(nBCEdgePts, -2.0, &tmp[0], 1, &tmp[0], 1);
98 
99  // Calculate v* = v - 2.0(v.n)n
100  for (i = 0; i < m_spacedim; ++i)
101  {
102  Vmath::Vvtvp(nBCEdgePts, &tmp[0], 1, &m_traceNormals[i][id2], 1,
103  &Fwd[1 + i][id2], 1, &Fwd[1 + i][id2], 1);
104  }
105 
106  // Copy boundary adjusted values into the boundary expansion
107  for (i = 0; i < nVariables; ++i)
108  {
109  Vmath::Vcopy(nBCEdgePts, &Fwd[i][id2], 1,
110  &(m_fields[i]
111  ->GetBndCondExpansions()[m_bcRegion]
112  ->UpdatePhys())[id1],
113  1);
114  }
115  }
116 }
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
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 Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255

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.
Definition: NekFactory.hpp:198
static CFSBndCondSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble >> &pTraceNormals, 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 64 of file SymmetryBC.h.