Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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:
Inheritance graph
[legend]
Collaboration diagram for Nektar::SymmetryBC:
Collaboration graph
[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)
 
- 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...
 

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...
 
- Protected Attributes inherited from Nektar::CFSBndCond
LibUtilities::SessionReaderSharedPtr m_session
 Session reader. More...
 
Array< OneD,
MultiRegions::ExpListSharedPtr
m_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_gamma
 Parameters of the flow. More...
 
NekDouble m_rhoInf
 
NekDouble m_pInf
 
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 48 of file SymmetryBC.h.

Constructor & Destructor Documentation

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.

54  : CFSBndCond(pSession, pFields, pTraceNormals, pSpaceDim, bcRegion, cnt)
55 {
56 }
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:51
virtual Nektar::SymmetryBC::~SymmetryBC ( void  )
inlineprivatevirtual

Definition at line 85 of file SymmetryBC.h.

85 {};

Member Function Documentation

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 55 of file SymmetryBC.h.

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

60  {
62  AllocateSharedPtr(pSession, pFields,
63  pTraceNormals, pSpaceDim, bcRegion, cnt);
64  return p;
65  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
boost::shared_ptr< CFSBndCond > CFSBndCondSharedPtr
A shared pointer to a boundary condition object.
Definition: CFSBndCond.h:49
void Nektar::SymmetryBC::v_Apply ( Array< OneD, Array< OneD, NekDouble > > &  Fwd,
Array< OneD, Array< OneD, NekDouble > > &  physarray,
const NekDouble time 
)
protectedvirtual

Implements Nektar::CFSBndCond.

Definition at line 58 of file SymmetryBC.cpp.

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().

62 {
63  int i;
64  int nVariables = physarray.num_elements();
65 
66  const Array<OneD, const int> &traceBndMap
67  = 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]->GetBndCondExpansions()[m_bcRegion]->
78  GetExp(e)->GetTotPoints();
79  id1 = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->
80  GetPhys_Offset(e);
81  id2 = m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[m_offset+e]);
82 
83  // For 2D/3D, define: v* = v - 2(v.n)n
84  Array<OneD, NekDouble> tmp(nBCEdgePts, 0.0);
85 
86  // Calculate (v.n)
87  for (i = 0; i < m_spacedim; ++i)
88  {
89  Vmath::Vvtvp(nBCEdgePts,
90  &Fwd[1+i][id2], 1,
91  &m_traceNormals[i][id2], 1,
92  &tmp[0], 1,
93  &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,
103  &tmp[0], 1,
104  &m_traceNormals[i][id2], 1,
105  &Fwd[1+i][id2], 1,
106  &Fwd[1+i][id2], 1);
107  }
108 
109  // Copy boundary adjusted values into the boundary expansion
110  for (i = 0; i < nVariables; ++i)
111  {
112  Vmath::Vcopy(nBCEdgePts, &Fwd[i][id2], 1,
113  &(m_fields[i]->GetBndCondExpansions()[m_bcRegion]->
114  UpdatePhys())[id1], 1);
115  }
116  }
117 }
int m_spacedim
Space dimension.
Definition: CFSBndCond.h:90
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:442
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array of fields.
Definition: CFSBndCond.h:86
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:213
int m_offset
Offset.
Definition: CFSBndCond.h:103
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Trace normals.
Definition: CFSBndCond.h:88
int m_bcRegion
Id of the boundary region.
Definition: CFSBndCond.h:101

Friends And Related Function Documentation

friend class MemoryManager< SymmetryBC >
friend

Definition at line 52 of file SymmetryBC.h.

Member Data Documentation

std::string Nektar::SymmetryBC::className
static
Initial value:
RegisterCreatorFunction("Symmetry",
"Symmetry boundary condition.")

Name of the class.

Definition at line 68 of file SymmetryBC.h.