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

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

#include <WallViscousBC.h>

Inheritance diagram for Nektar::WallViscousBC:
[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 classNameViscous
 Name of the class. More...
 
static std::string classNameAdiabatic
 

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

Protected Attributes

Array< OneD, Array< OneD, NekDouble > > m_bndPhys
 
- 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...
 

Private Member Functions

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

Friends

class MemoryManager< WallViscousBC >
 

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

Detailed Description

Wall boundary conditions for viscous compressible flow problems.

Definition at line 46 of file WallViscousBC.h.

Constructor & Destructor Documentation

◆ WallViscousBC()

Nektar::WallViscousBC::WallViscousBC ( 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 52 of file WallViscousBC.cpp.

58 : CFSBndCond(pSession, pFields, pTraceNormals, pGridVelocity, pSpaceDim,
59 bcRegion, cnt)
60{
62
63 m_bndPhys = Array<OneD, Array<OneD, NekDouble>>(m_fields.size());
64}
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
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array of fields.
Definition: CFSBndCond.h:92
Array< OneD, Array< OneD, NekDouble > > m_bndPhys
Definition: WallViscousBC.h:72

References m_bndPhys, Nektar::CFSBndCond::m_diffusionAveWeight, and Nektar::CFSBndCond::m_fields.

◆ ~WallViscousBC()

Nektar::WallViscousBC::~WallViscousBC ( void  )
inlineoverrideprivate

Definition at line 85 of file WallViscousBC.h.

85{};

Member Function Documentation

◆ create()

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

Implements Nektar::CFSBndCond.

Definition at line 66 of file WallViscousBC.cpp.

69{
70 int i;
71 int nVariables = physarray.size();
72
73 // Find the fields whose WallViscous/Adiabatic-BC is time-dependent
74 // Update variables on the boundaries of these fields
75 // Get the updated variables on the WallViscous/Adiabatic boundary
76 //
77 // Maybe the EvaluateBoundaryConditions() should be put upstream to
78 // CompressibleFlowSystem::NumCalRiemFluxJac(), So that the BCs will not
79 // be repeatedly updated when there are more than one time-dependent BC.
80 std::string varName;
81 for (i = 0; i < nVariables; ++i)
82 {
83 if (m_fields[i]->GetBndConditions()[m_bcRegion]->IsTimeDependent())
84 {
85 varName = m_session->GetVariable(i);
86 m_fields[i]->EvaluateBoundaryConditions(time, varName);
87
88 m_bndPhys[i] =
89 m_fields[i]->GetBndCondExpansions()[m_bcRegion]->UpdatePhys();
90 }
91 }
92
93 const Array<OneD, const int> &traceBndMap = m_fields[0]->GetTraceBndMap();
94
95 // Take into account that for PDE based shock capturing, eps = 0 at the
96 // wall. Adjust the physical values of the trace to take user defined
97 // boundaries into account
98 int e, id1, id2, nBCEdgePts, eMax;
99
100 eMax = m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetExpSize();
101
102 for (e = 0; e < eMax; ++e)
103 {
104 nBCEdgePts = m_fields[0]
105 ->GetBndCondExpansions()[m_bcRegion]
106 ->GetExp(e)
107 ->GetTotPoints();
108 id1 =
109 m_fields[0]->GetBndCondExpansions()[m_bcRegion]->GetPhys_Offset(e);
110 id2 =
111 m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[m_offset + e]);
112
113 // Boundary condition for epsilon term.
114 if (nVariables == m_spacedim + 3)
115 {
116 Vmath::Zero(nBCEdgePts, &Fwd[nVariables - 1][id2], 1);
117 }
118
119 // V = - Vin
120 for (i = 0; i < m_spacedim; i++)
121 {
122 // V = -Vin
123 // Vmath::Neg(nBCEdgePts, &Fwd[i+1][id2], 1);
124
125 // This now does Vg * rho + Vin
126 // Vmath::Vvtvp(nBCEdgePts, &m_gridVelocityTrace[i][id2], 1,
127 // &Fwd[0][id2], 1, &Fwd[i+1][id2], 1, &Fwd[i+1][id2], 1);
128
129 for (int j = 0; j < nBCEdgePts; ++j)
130 {
131 Fwd[i + 1][id2 + j] =
132 2 * m_gridVelocityTrace[i][id2 + j] * Fwd[0][id2 + j] -
133 Fwd[i + 1][id2 + j];
134 }
135 }
136
137 // Superimpose the perturbation
138 for (i = 0; i < nVariables; ++i)
139 {
140 if (m_fields[i]->GetBndConditions()[m_bcRegion]->IsTimeDependent())
141 {
142 Vmath::Vadd(nBCEdgePts, &m_bndPhys[i][id1], 1, &Fwd[i][id2], 1,
143 &Fwd[i][id2], 1);
144 }
145 }
146
147 // Copy boundary adjusted values into the boundary expansion
148 for (i = 0; i < nVariables; ++i)
149 {
150 Vmath::Vcopy(nBCEdgePts, &Fwd[i][id2], 1,
151 &(m_fields[i]
152 ->GetBndCondExpansions()[m_bcRegion]
153 ->UpdatePhys())[id1],
154 1);
155 }
156 }
157}
LibUtilities::SessionReaderSharedPtr m_session
Session reader.
Definition: CFSBndCond.h:90
int m_spacedim
Space dimension.
Definition: CFSBndCond.h:98
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
int m_offset
Offset.
Definition: CFSBndCond.h:115
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.hpp:180
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, m_bndPhys, Nektar::CFSBndCond::m_fields, Nektar::CFSBndCond::m_gridVelocityTrace, Nektar::CFSBndCond::m_offset, Nektar::CFSBndCond::m_session, Nektar::CFSBndCond::m_spacedim, Vmath::Vadd(), Vmath::Vcopy(), and Vmath::Zero().

Friends And Related Function Documentation

◆ MemoryManager< WallViscousBC >

friend class MemoryManager< WallViscousBC >
friend

Definition at line 1 of file WallViscousBC.h.

Member Data Documentation

◆ classNameAdiabatic

std::string Nektar::WallViscousBC::classNameAdiabatic
static
Initial value:
=
"WallAdiabatic", WallViscousBC::create,
"Adiabatic 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: WallViscousBC.h:52
CFSBndCondFactory & GetCFSBndCondFactory()
Declaration of the boundary condition factory singleton.
Definition: CFSBndCond.cpp:41

Definition at line 67 of file WallViscousBC.h.

◆ classNameViscous

std::string Nektar::WallViscousBC::classNameViscous
static
Initial value:
=
"WallViscous", WallViscousBC::create,
"No-slip (viscous) wall boundary condition.")

Name of the class.

Definition at line 66 of file WallViscousBC.h.

◆ m_bndPhys

Array<OneD, Array<OneD, NekDouble> > Nektar::WallViscousBC::m_bndPhys
protected

Definition at line 72 of file WallViscousBC.h.

Referenced by v_Apply(), and WallViscousBC().