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

Outflow characteristic boundary conditions for compressible flow problems. More...

#include <EnforceEntropyVelocity.h>

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

void v_Apply (Array< OneD, Array< OneD, NekDouble > > &Fwd, Array< OneD, Array< OneD, NekDouble > > &physarray, const NekDouble &time) override
 
void GenerateRotationMatrices (const Array< OneD, const Array< OneD, NekDouble > > &normals)
 
void FromToRotation (Array< OneD, const NekDouble > &from, Array< OneD, const NekDouble > &to, NekDouble *mat)
 
SOLVER_UTILS_EXPORT void rotateToNormal (const Array< OneD, const Array< OneD, NekDouble > > &inarray, const Array< OneD, const Array< OneD, NekDouble > > &normals, const Array< OneD, const Array< OneD, NekDouble > > &vecLocs, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
SOLVER_UTILS_EXPORT void rotateFromNormal (const Array< OneD, const Array< OneD, NekDouble > > &inarray, const Array< OneD, const Array< OneD, NekDouble > > &normals, const Array< OneD, const Array< OneD, NekDouble > > &vecLocs, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
- 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...
 
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

int m_npts
 
Array< OneD, NekDoublem_rhoBC
 
Array< OneD, Array< OneD, NekDouble > > m_velBC
 
Array< OneD, NekDoublem_pBC
 
Array< OneD, NekDoublem_VnInf
 Reference normal velocity. More...
 
Array< OneD, int > m_bndToTraceMap
 
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...
 
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
 
int m_bcRegion
 Id of the boundary region. More...
 
int m_offset
 Offset. More...
 

Private Member Functions

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

Friends

class MemoryManager< EnforceEntropyVelocity >
 

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

Outflow characteristic boundary conditions for compressible flow problems.

Definition at line 51 of file EnforceEntropyVelocity.h.

Constructor & Destructor Documentation

◆ EnforceEntropyVelocity()

Nektar::EnforceEntropyVelocity::EnforceEntropyVelocity ( 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 50 of file EnforceEntropyVelocity.cpp.

55 : CFSBndCond(pSession, pFields, pTraceNormals, pSpaceDim, bcRegion, cnt)
56{
57
59 m_fields[0]->GetBndCondExpansions()[m_bcRegion];
60
61 //-> Gather a list of index from trace to this boundary
62 m_npts = bndexp->GetTotPoints();
63
64 m_bndToTraceMap = Array<OneD, int>(m_npts, -1);
65
66 const Array<OneD, const int> &traceBndMap = m_fields[0]->GetTraceBndMap();
67
68 // Construct a map for the boundary to trace map for easy acess to
69 // phys space points
70 int cnt1 = 0;
71 for (int e = 0; e < bndexp->GetNumElmts(); ++e)
72 {
73 int nTracePts = bndexp->GetExp(e)->GetTotPoints();
74
75 int id =
76 m_fields[0]->GetTrace()->GetPhys_Offset(traceBndMap[m_offset + e]);
77
78 // Loop on the points of the m_bcRegion
79 for (int i = 0; i < nTracePts; i++)
80 {
81 // the ith point in region e
82 m_bndToTraceMap[cnt1++] = id + i;
83 }
84 }
85
86 Array<OneD, Array<OneD, NekDouble>> BCvals(m_fields.size());
87 m_bndPhys = Array<OneD, Array<OneD, NekDouble>>(m_fields.size());
88
89 for (int i = 0; i < m_fields.size(); ++i)
90 {
91 m_bndPhys[i] =
92 m_fields[i]->GetBndCondExpansions()[m_bcRegion]->UpdatePhys();
93
94 BCvals[i] = Array<OneD, NekDouble>(m_npts);
95 Vmath::Vcopy(m_npts, m_bndPhys[i], 1, BCvals[i], 1);
96 }
97
98 // Set up boudnary required BCs
99 m_rhoBC = Array<OneD, NekDouble>(m_npts);
100 Vmath::Vcopy(m_npts, BCvals[0], 1, m_rhoBC, 1);
101
102 m_velBC = Array<OneD, Array<OneD, NekDouble>>(m_spacedim);
103 // Evaluate velocity on boundary
104 for (int i = 0; i < m_spacedim; ++i)
105 {
106 m_velBC[i] = Array<OneD, NekDouble>(m_npts);
107 Vmath::Vcopy(m_npts, BCvals[i + 1], 1, m_velBC[i], 1);
108 Vmath::Vdiv(m_npts, m_velBC[i], 1, m_rhoBC, 1, m_velBC[i], 1);
109 }
110 m_pBC = Array<OneD, NekDouble>(m_npts);
111 m_varConv->GetPressure(BCvals, m_pBC);
112
113 // Computing the normal velocity for characteristics coming
114 // from outside the computational domain
115 m_VnInf = Array<OneD, NekDouble>(m_npts, 0.0);
116 for (int i = 0; i < m_spacedim; i++)
117 {
118 for (int j = 0; j < m_npts; ++j)
119 {
120 m_VnInf[j] += m_traceNormals[i][m_bndToTraceMap[j]] * m_velBC[i][j];
121 }
122 }
123}
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:47
int m_spacedim
Space dimension.
Definition: CFSBndCond.h:95
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Trace normals.
Definition: CFSBndCond.h:93
int m_bcRegion
Id of the boundary region.
Definition: CFSBndCond.h:109
VariableConverterSharedPtr m_varConv
Auxiliary object to convert variables.
Definition: CFSBndCond.h:97
int m_offset
Offset.
Definition: CFSBndCond.h:111
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array of fields.
Definition: CFSBndCond.h:91
Array< OneD, NekDouble > m_VnInf
Reference normal velocity.
Array< OneD, NekDouble > m_rhoBC
Array< OneD, Array< OneD, NekDouble > > m_bndPhys
Array< OneD, Array< OneD, NekDouble > > m_velBC
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
void Vdiv(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x/y.
Definition: Vmath.hpp:126
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, m_bndToTraceMap, Nektar::CFSBndCond::m_fields, m_npts, Nektar::CFSBndCond::m_offset, m_pBC, m_rhoBC, Nektar::CFSBndCond::m_spacedim, Nektar::CFSBndCond::m_traceNormals, Nektar::CFSBndCond::m_varConv, m_velBC, m_VnInf, Vmath::Vcopy(), and Vmath::Vdiv().

◆ ~EnforceEntropyVelocity()

Nektar::EnforceEntropyVelocity::~EnforceEntropyVelocity ( void  )
inlineoverrideprivate

Definition at line 118 of file EnforceEntropyVelocity.h.

118{};

Member Function Documentation

◆ create()

static CFSBndCondSharedPtr Nektar::EnforceEntropyVelocity::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 57 of file EnforceEntropyVelocity.h.

62 {
65 pSession, pFields, pTraceNormals, pSpaceDim, bcRegion, cnt);
66 return p;
67 }
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.

◆ FromToRotation()

void Nektar::EnforceEntropyVelocity::FromToRotation ( Array< OneD, const NekDouble > &  from,
Array< OneD, const NekDouble > &  to,
NekDouble mat 
)
protected

◆ GenerateRotationMatrices()

void Nektar::EnforceEntropyVelocity::GenerateRotationMatrices ( const Array< OneD, const Array< OneD, NekDouble > > &  normals)
protected

◆ rotateFromNormal()

SOLVER_UTILS_EXPORT void Nektar::EnforceEntropyVelocity::rotateFromNormal ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
const Array< OneD, const Array< OneD, NekDouble > > &  normals,
const Array< OneD, const Array< OneD, NekDouble > > &  vecLocs,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
protected

◆ rotateToNormal()

SOLVER_UTILS_EXPORT void Nektar::EnforceEntropyVelocity::rotateToNormal ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
const Array< OneD, const Array< OneD, NekDouble > > &  normals,
const Array< OneD, const Array< OneD, NekDouble > > &  vecLocs,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
protected

◆ v_Apply()

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

Implements Nektar::CFSBndCond.

Definition at line 125 of file EnforceEntropyVelocity.cpp.

129{
130 int i, j;
131 int nDimensions = m_spacedim;
132
133 Array<OneD, Array<OneD, NekDouble>> FwdBnd(Fwd.size());
134 Array<OneD, Array<OneD, NekDouble>> bndPhys(Fwd.size());
135
136 // make a local copy of Fwd along boundary of interest
137 for (i = 0; i < Fwd.size(); ++i)
138 {
139 FwdBnd[i] = Array<OneD, NekDouble>(m_npts);
140 for (j = 0; j < m_npts; ++j)
141 {
142 FwdBnd[i][j] = Fwd[i][m_bndToTraceMap[j]];
143 }
144 }
145
146 // Computing the normal velocity for characteristics coming
147 // from inside the computational domain
148 Array<OneD, NekDouble> Vn(m_npts, 0.0);
149
150 for (i = 0; i < nDimensions; ++i)
151 {
152 for (j = 0; j < m_npts; ++j)
153 {
154 Vn[j] += m_traceNormals[i][m_bndToTraceMap[j]] * FwdBnd[i + 1][j];
155 }
156 }
157 // divide by density.
158 Vmath::Vdiv(m_npts, Vn, 1, FwdBnd[0], 1, Vn, 1);
159
160 // Get speed of sound
161 Array<OneD, NekDouble> pressure(m_npts);
162 Array<OneD, NekDouble> soundSpeed(m_npts);
163
164 m_varConv->GetPressure(FwdBnd, pressure);
165 m_varConv->GetSoundSpeed(FwdBnd, soundSpeed);
166
167 // Get Mach. Note: it is computed by Vn/c
168 Array<OneD, NekDouble> Mach(m_npts, 0.0);
169 Vmath::Vdiv(m_npts, Vn, 1, soundSpeed, 1, Mach, 1);
170 Vmath::Vabs(m_npts, Mach, 1, Mach, 1);
171
172 // Auxiliary variables
173 Array<OneD, NekDouble> velBC(nDimensions, 0.0);
174
175 // L represents properties outside boundary (ghost state)
176 // R represents properties inside boundary (numerical state)
177 NekDouble rhoL, uL, pL, cL, sL;
178 NekDouble cR, uR; // rR, r = RiemannInvariant
179 NekDouble EBC; // cstar, pstar, rhostar, ustar, BC = star
180
181 NekDouble gamMinOne = m_gamma - 1.0;
182 NekDouble oneOverGamMinOne = 1.0 / gamMinOne;
183 NekDouble twoOverGamMinOne = 2.0 * oneOverGamMinOne;
184 NekDouble gamInv = 1.0 / m_gamma;
185
186 // Loop on m_bcRegions
187 for (int pnt = 0; pnt < m_npts; ++pnt)
188 {
189 // Impose inflow Riemann invariant
190 if (Vn[pnt] <= 0.0)
191 {
192 // Subsonic flows
193 if (Mach[pnt] < 1.00)
194 {
195 // right state variables
196 cR = sqrt(m_gamma * pressure[pnt] / FwdBnd[0][pnt]);
197 uR = -Vn[pnt];
198
199 uL = 2 * (-m_VnInf[pnt]) - uR;
200 cL = cR;
201 sL = m_pBC[pnt] / (pow(m_rhoBC[pnt], m_gamma));
202 rhoL = pow((cL * cL) / (m_gamma * sL), oneOverGamMinOne);
203 pL = cL * cL * rhoL * gamInv;
204 }
205 else // Supersonic inflow
206 {
207 // all characteristics are from left so just impose
208 // star state to left values
209 // Note: m_vnInf is the negative of the normal velocity
210 // across boundary
211 rhoL = m_rhoBC[pnt];
212 uL = -m_VnInf[pnt];
213 pL = m_pBC[pnt];
214 }
215
216 // Boundary energy
217 EBC = pL * twoOverGamMinOne * 0.5;
218
219 // evaluate the different between the left state normal
220 // velocity and that from the desired condition (note
221 // m_VnInf is using an outwards normal definition.
222 NekDouble VnDiff = uL + m_VnInf[pnt];
223
224 // Note: Can just use the BC values directly!!
225 for (j = 0; j < nDimensions; ++j)
226 {
227 // Set velocity to the desired conditions modified to
228 // take account of the normal state for Riemann
229 // problem. (Negative accounts for outwards normal definition)
230 // velBC[j] = m_velBC[j][pnt];
231 velBC[j] = m_velBC[j][pnt] -
232 VnDiff * m_traceNormals[j][m_bndToTraceMap[pnt]];
233
234 EBC += 0.5 * rhoL * velBC[j] * velBC[j];
235 }
236
237 //-------------------------------------------------------------------------
238 // Impose Left hand Riemann Invariant boundary conditions
239 m_bndPhys[0][pnt] = rhoL;
240 for (j = 0; j < nDimensions; ++j)
241 {
242 m_bndPhys[j + 1][pnt] = rhoL * velBC[j];
243 }
244 m_bndPhys[nDimensions + 1][pnt] = EBC;
245 }
246 else // Outflow
247 {
248
249 // Note: Allowing the switch can cause worse convergence in this
250 // type BC.
251 // So improve it later.
252 if (Mach[pnt] < 1.00)
253 {
254 rhoL = m_rhoBC[pnt];
255 uL = -m_VnInf[pnt];
256 pL = m_pBC[pnt];
257 }
258 else
259 {
260 // supersonic outflow
261 // Just set to imposed state and let Riemann BC dictate values
262 rhoL = m_rhoBC[pnt];
263 uL = -m_VnInf[pnt];
264 pL = m_pBC[pnt];
265 }
266
267 // Boundary energy
268 EBC = pL * twoOverGamMinOne * 0.5;
269
270 // Boundary velocities & Kinite energy
271 // Note: normals are negated since they point outwards in
272 // the domain
273 for (j = 0; j < nDimensions; ++j)
274 {
275 velBC[j] = -1.0 * uL * m_traceNormals[j][m_bndToTraceMap[pnt]];
276 EBC += 0.5 * rhoL * velBC[j] * velBC[j];
277 }
278
279 // Impose Left hand Riemann Invariant boundary conditions
280 m_bndPhys[0][pnt] = rhoL;
281 for (j = 0; j < nDimensions; ++j)
282 {
283 m_bndPhys[j + 1][pnt] = rhoL * velBC[j];
284 }
285 m_bndPhys[nDimensions + 1][pnt] = EBC;
286 }
287 }
288}
NekDouble m_gamma
Parameters of the flow.
Definition: CFSBndCond.h:102
double NekDouble
void Vabs(int n, const T *x, const int incx, T *y, const int incy)
vabs: y = |x|
Definition: Vmath.hpp:352
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References m_bndPhys, m_bndToTraceMap, Nektar::CFSBndCond::m_gamma, m_npts, m_pBC, m_rhoBC, Nektar::CFSBndCond::m_spacedim, Nektar::CFSBndCond::m_traceNormals, Nektar::CFSBndCond::m_varConv, m_velBC, m_VnInf, CG_Iterations::pressure, tinysimd::sqrt(), Vmath::Vabs(), and Vmath::Vdiv().

Friends And Related Function Documentation

◆ MemoryManager< EnforceEntropyVelocity >

friend class MemoryManager< EnforceEntropyVelocity >
friend

Definition at line 1 of file EnforceEntropyVelocity.h.

Member Data Documentation

◆ className

std::string Nektar::EnforceEntropyVelocity::className
static
Initial value:
=
"EnforceEntropyVelocity", EnforceEntropyVelocity::create,
"Riemann boundary condition enforcing entropy and velocity.")
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.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
CFSBndCondFactory & GetCFSBndCondFactory()
Declaration of the boundary condition factory singleton.
Definition: CFSBndCond.cpp:41

Name of the class.

Definition at line 70 of file EnforceEntropyVelocity.h.

◆ m_bndPhys

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

Definition at line 87 of file EnforceEntropyVelocity.h.

Referenced by EnforceEntropyVelocity(), and v_Apply().

◆ m_bndToTraceMap

Array<OneD, int> Nektar::EnforceEntropyVelocity::m_bndToTraceMap
protected

Definition at line 83 of file EnforceEntropyVelocity.h.

Referenced by EnforceEntropyVelocity(), and v_Apply().

◆ m_npts

int Nektar::EnforceEntropyVelocity::m_npts
protected

Definition at line 73 of file EnforceEntropyVelocity.h.

Referenced by EnforceEntropyVelocity(), and v_Apply().

◆ m_pBC

Array<OneD, NekDouble> Nektar::EnforceEntropyVelocity::m_pBC
protected

Definition at line 79 of file EnforceEntropyVelocity.h.

Referenced by EnforceEntropyVelocity(), and v_Apply().

◆ m_rhoBC

Array<OneD, NekDouble> Nektar::EnforceEntropyVelocity::m_rhoBC
protected

Definition at line 75 of file EnforceEntropyVelocity.h.

Referenced by EnforceEntropyVelocity(), and v_Apply().

◆ m_velBC

Array<OneD, Array<OneD, NekDouble> > Nektar::EnforceEntropyVelocity::m_velBC
protected

Definition at line 77 of file EnforceEntropyVelocity.h.

Referenced by EnforceEntropyVelocity(), and v_Apply().

◆ m_VnInf

Array<OneD, NekDouble> Nektar::EnforceEntropyVelocity::m_VnInf
protected

Reference normal velocity.

Definition at line 81 of file EnforceEntropyVelocity.h.

Referenced by EnforceEntropyVelocity(), and v_Apply().