Nektar++
Loading...
Searching...
No Matches
Static Public Member Functions | Static Public Attributes | Protected Member Functions | Private Member Functions | Friends | List of all members
Nektar::WallRotationalBC Class Reference

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

#include <WallRotationalBC.h>

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

Static Public Attributes

static std::string classNameRotational
 Name of the class.
 

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.
 
virtual ~CFSBndCond ()=default
 
virtual void v_ApplyBwdWeight ()
 

Private Member Functions

 WallRotationalBC (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)
 
 ~WallRotationalBC (void) override=default
 

Friends

class MemoryManager< WallRotationalBC >
 

Additional Inherited Members

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

Detailed Description

Wall boundary conditions for Rotational compressible flow problems.

Definition at line 46 of file WallRotationalBC.h.

Constructor & Destructor Documentation

◆ WallRotationalBC()

Nektar::WallRotationalBC::WallRotationalBC ( 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 WallRotationalBC.cpp.

51 : CFSBndCond(pSession, pFields, pTraceNormals, pGridVelocity, pSpaceDim,
52 bcRegion, cnt)
53{
55
56 // Set up rotational boundary edge velocities
57 const Array<OneD, const int> &traceBndMap = m_fields[0]->GetTraceBndMap();
58
59 int eMax = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExpSize();
60 for (int e = 0; e < eMax; ++e)
61 {
62 int nBCEdgePts = m_fields[0]
63 ->GetBndCondExpansions()[m_bcRegion]
64 ->GetExp(e)
65 ->GetTotPoints();
66 int id2 =
67 m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[m_offset + e]);
68
69 Array<OneD, NekDouble> x(nBCEdgePts, 0.0);
70 Array<OneD, NekDouble> y(nBCEdgePts, 0.0);
71 Array<OneD, NekDouble> z(nBCEdgePts, 0.0);
72 m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExp(e)->GetCoords(
73 x, y, z);
74 for (int j = 0; j < nBCEdgePts; ++j)
75 {
76 m_gridVelocityTrace[0][id2 + j] = m_angVel * -y[j];
77 m_gridVelocityTrace[1][id2 + j] = m_angVel * x[j];
78 }
79 }
80}
NekDouble m_diffusionAveWeight
Weight for average calculation of diffusion term.
Definition CFSBndCond.h:102
Array< OneD, Array< OneD, NekDouble > > m_gridVelocityTrace
Grid Velocity.
Definition CFSBndCond.h:96
int m_bcRegion
Id of the boundary region.
Definition CFSBndCond.h:113
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.
int m_offset
Offset.
Definition CFSBndCond.h:115
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array of fields.
Definition CFSBndCond.h:92
std::vector< double > z(NPUPPER)

References Nektar::CFSBndCond::m_angVel, Nektar::CFSBndCond::m_bcRegion, Nektar::CFSBndCond::m_diffusionAveWeight, Nektar::CFSBndCond::m_fields, Nektar::CFSBndCond::m_gridVelocityTrace, and Nektar::CFSBndCond::m_offset.

◆ ~WallRotationalBC()

Nektar::WallRotationalBC::~WallRotationalBC ( void  )
overrideprivatedefault

Member Function Documentation

◆ create()

static CFSBndCondSharedPtr Nektar::WallRotationalBC::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 WallRotationalBC.h.

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

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

◆ v_Apply()

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

Implements Nektar::CFSBndCond.

Definition at line 82 of file WallRotationalBC.cpp.

85{
86 int nVariables = physarray.size();
87
88 // @TODO: ALE we subtract the grid velocity ? "Set u = to ug for this one" -
89 // Dave
90
91 int i;
92 const Array<OneD, const int> &traceBndMap = m_fields[0]->GetTraceBndMap();
93
94 // Take into account that for PDE based shock capturing, eps = 0 at the
95 // wall. Adjust the physical values of the trace to take user defined
96 // boundaries into account
97 int e, id1, id2, nBCEdgePts, eMax;
98
99 eMax = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExpSize();
100
101 for (e = 0; e < eMax; ++e)
102 {
103 nBCEdgePts = m_fields[0]
104 ->GetBndCondExpansions()[m_bcRegion]
105 ->GetExp(e)
106 ->GetTotPoints();
107 id1 =
108 m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetPhys_Offset(e);
109 id2 =
110 m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[m_offset + e]);
111
112 // Boundary condition for epsilon term. @TODO: Is this correct, or
113 // should I do E = p/(gamma -1) + 1/2*rho(u^2 +v^2 + w^2)... or
114 if (nVariables == m_spacedim + 3)
115 {
116 Vmath::Zero(nBCEdgePts, &Fwd[nVariables - 1][id2], 1);
117 }
118
119 for (i = 0; i < m_spacedim; i++)
120 {
121 // V = -Vin
122 // Vmath::Neg(nBCEdgePts, &Fwd[i+1][id2], 1);
123
124 // This now does Vg * rho + Vin
125 // Vmath::Vvtvp(nBCEdgePts, &m_gridVelocityTrace[i][id2], 1,
126 // &Fwd[0][id2], 1, &Fwd[i+1][id2], 1, &Fwd[i+1][id2], 1);
127
128 for (int j = 0; j < nBCEdgePts; ++j)
129 {
130 // Fwd[i+1][id2 + j] = 2 * m_gridVelocityTrace[i][id2 + j] *
131 // Fwd[0][id2 + j] - Fwd[i+1][id2 + j];
132 Fwd[i + 1][id2 + j] =
133 2 * m_gridVelocityTrace[i][id2 + j] * Fwd[0][id2 + j] -
134 Fwd[i + 1][id2 + j];
135 }
136 }
137
138 // Copy boundary adjusted values into the boundary expansion
139 for (i = 0; i < nVariables; ++i)
140 {
141 Vmath::Vcopy(nBCEdgePts, &Fwd[i][id2], 1,
142 &(m_fields[i]
143 ->GetBndCondExpansions()[m_bcRegion]
144 ->UpdatePhys())[id1],
145 1);
146 }
147 }
148}
int m_spacedim
Space dimension.
Definition CFSBndCond.h:98
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, Vmath::Vcopy(), and Vmath::Zero().

Friends And Related Symbol Documentation

◆ MemoryManager< WallRotationalBC >

friend class MemoryManager< WallRotationalBC >
friend

Definition at line 1 of file WallRotationalBC.h.

Member Data Documentation

◆ classNameRotational

std::string Nektar::WallRotationalBC::classNameRotational
static
Initial value:
=
"WallRotational", WallRotationalBC::create,
"Adiabatic rotational 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.
CFSBndCondFactory & GetCFSBndCondFactory()
Declaration of the boundary condition factory singleton.

Name of the class.

Definition at line 67 of file WallRotationalBC.h.