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

Static Public Attributes

static std::string classNameRotational
 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

 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
 

Friends

class MemoryManager< WallRotationalBC >
 

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 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 49 of file WallRotationalBC.cpp.

55 : CFSBndCond(pSession, pFields, pTraceNormals, pGridVelocity, pSpaceDim,
56 bcRegion, cnt)
57{
59
60 // Set up rotational boundary edge velocities
61 const Array<OneD, const int> &traceBndMap = m_fields[0]->GetTraceBndMap();
62
63 int eMax = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExpSize();
64 for (int e = 0; e < eMax; ++e)
65 {
66 int nBCEdgePts = m_fields[0]
67 ->GetBndCondExpansions()[m_bcRegion]
68 ->GetExp(e)
69 ->GetTotPoints();
70 int id2 =
71 m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[m_offset + e]);
72
73 Array<OneD, NekDouble> x(nBCEdgePts, 0.0);
74 Array<OneD, NekDouble> y(nBCEdgePts, 0.0);
75 m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExp(e)->GetCoords(
76 x, y);
77 for (int j = 0; j < nBCEdgePts; ++j)
78 {
79 m_gridVelocityTrace[0][id2 + j] = m_angVel * -y[j];
80 m_gridVelocityTrace[1][id2 + j] = m_angVel * x[j];
81 }
82 }
83}
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.
Definition: CFSBndCond.cpp:47
int m_offset
Offset.
Definition: CFSBndCond.h:115
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array of fields.
Definition: CFSBndCond.h:92
NekDouble m_angVel
Definition: CFSBndCond.h:110

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

Definition at line 81 of file WallRotationalBC.h.

81{};

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::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::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 85 of file WallRotationalBC.cpp.

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

Name of the class.

Definition at line 67 of file WallRotationalBC.h.