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

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

#include <WallBC.h>

Inheritance diagram for Nektar::WallBC:
[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

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

Friends

class MemoryManager< WallBC >
 

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

Wall boundary conditions for compressible flow problems.

Definition at line 46 of file WallBC.h.

Constructor & Destructor Documentation

◆ WallBC()

Nektar::WallBC::WallBC ( 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 45 of file WallBC.cpp.

50 : CFSBndCond(pSession, pFields, pTraceNormals, pGridVelocity, pSpaceDim,
51 bcRegion, cnt)
52{
54}
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.

◆ ~WallBC()

Nektar::WallBC::~WallBC ( void  )
inlineoverrideprivate

Definition at line 80 of file WallBC.h.

80{};

Member Function Documentation

◆ create()

static CFSBndCondSharedPtr Nektar::WallBC::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 WallBC.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::WallBC::v_Apply ( Array< OneD, Array< OneD, NekDouble > > &  Fwd,
Array< OneD, Array< OneD, NekDouble > > &  physarray,
const NekDouble time 
)
overrideprotectedvirtual

Implements Nektar::CFSBndCond.

Definition at line 56 of file WallBC.cpp.

59{
60 int i;
61 int nVariables = physarray.size();
62
63 const Array<OneD, const int> &traceBndMap = m_fields[0]->GetTraceBndMap();
64
65 // Adjust the physical values of the trace to take
66 // user defined boundaries into account
67 int e, id1, id2, nBCEdgePts, eMax;
68
69 eMax = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExpSize();
70
71 for (e = 0; e < eMax; ++e)
72 {
73 nBCEdgePts = m_fields[0]
74 ->GetBndCondExpansions()[m_bcRegion]
75 ->GetExp(e)
76 ->GetTotPoints();
77 id1 =
78 m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetPhys_Offset(e);
79 id2 =
80 m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[m_offset + e]);
81
82 // Boundary condition for epsilon term.
83 if (nVariables == m_spacedim + 3)
84 {
85 Vmath::Zero(nBCEdgePts, &Fwd[nVariables - 1][id2], 1);
86 }
87
88 // @TODO: Look at paper on this
89 // https://www.researchgate.net/publication/264044118_A_Guide_to_the_Implementation_of_Boundary_Conditions_in_Compact_High-Order_Methods_for_Compressible_Aerodynamics
90 // For 2D/3D, define: v* = v - 2(v.n)n
91 Array<OneD, NekDouble> tmp(nBCEdgePts, 0.0);
92
93 //@TODO: v - vg here... check nguyen paper, only issue is getting the vg
94 // for the trace in here
95 //@TODO: Update m_traceNormals, might be fine though.
96
97 if (m_fields[0]->GetGraph()->GetMovement()->GetMoveFlag())
98 {
99 for (i = 0; i < m_spacedim; ++i)
100 {
101 // This now does Vg * rho + Vin
102 for (int j = 0; j < nBCEdgePts; ++j)
103 {
104 Fwd[i + 1][id2 + j] +=
105 m_gridVelocityTrace[i][id2 + j] * Fwd[0][id2 + j];
106 }
107 }
108 }
109
110 // Calculate (v.n)
111 for (i = 0; i < m_spacedim; ++i)
112 {
113 Vmath::Vvtvp(nBCEdgePts, &Fwd[1 + i][id2], 1,
114 &m_traceNormals[i][id2], 1, &tmp[0], 1, &tmp[0], 1);
115 }
116
117 // Calculate 2.0(v.n)
118 Vmath::Smul(nBCEdgePts, -2.0, &tmp[0], 1, &tmp[0], 1);
119
120 // Calculate v* = v - 2.0(v.n)n
121 for (i = 0; i < m_spacedim; ++i)
122 {
123 Vmath::Vvtvp(nBCEdgePts, &tmp[0], 1, &m_traceNormals[i][id2], 1,
124 &Fwd[1 + i][id2], 1, &Fwd[1 + i][id2], 1);
125 }
126
127 // Copy boundary adjusted values into the boundary expansion
128 for (i = 0; i < nVariables; ++i)
129 {
130 Vmath::Vcopy(nBCEdgePts, &Fwd[i][id2], 1,
131 &(m_fields[i]
132 ->GetBndCondExpansions()[m_bcRegion]
133 ->UpdatePhys())[id1],
134 1);
135 }
136 }
137}
int m_spacedim
Space dimension.
Definition: CFSBndCond.h:98
Array< OneD, Array< OneD, NekDouble > > m_gridVelocityTrace
Grid Velocity.
Definition: CFSBndCond.h:96
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 Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
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_gridVelocityTrace, Nektar::CFSBndCond::m_offset, Nektar::CFSBndCond::m_spacedim, Nektar::CFSBndCond::m_traceNormals, Vmath::Smul(), Vmath::Vcopy(), Vmath::Vvtvp(), and Vmath::Zero().

Friends And Related Function Documentation

◆ MemoryManager< WallBC >

friend class MemoryManager< WallBC >
friend

Definition at line 1 of file WallBC.h.

Member Data Documentation

◆ className

std::string Nektar::WallBC::className
static
Initial value:
"Wall", WallBC::create, "Slip wall 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: WallBC.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 WallBC.h.