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

#include <NavierStokesAdvection.h>

Inheritance diagram for Nektar::NavierStokesAdvection:
[legend]

Public Member Functions

void SetSpecHPDealiasing (bool value)
 
- Public Member Functions inherited from Nektar::SolverUtils::Advection
virtual SOLVER_UTILS_EXPORT ~Advection ()
 
SOLVER_UTILS_EXPORT void InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
 Interface function to initialise the advection object. More...
 
SOLVER_UTILS_EXPORT void Advect (const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &advVel, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble &time, const Array< OneD, Array< OneD, NekDouble > > &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayOfArray)
 Interface function to advect the vector field. More...
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetFluxVector (FuncPointerT func, ObjectPointerT obj)
 Set the flux vector callback function. More...
 
void SetRiemannSolver (RiemannSolverSharedPtr riemann)
 Set a Riemann solver object for this advection object. More...
 
void SetFluxVector (AdvectionFluxVecCB fluxVector)
 Set the flux vector callback function. More...
 
void SetBaseFlow (const Array< OneD, Array< OneD, NekDouble > > &inarray, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields)
 Set the base flow used for linearised advection objects. More...
 
template<typename DataType , typename TypeNekBlkMatSharedPtr >
SOLVER_UTILS_EXPORT void AddTraceJacToMat (const int nConvectiveFields, const int nSpaceDim, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, TypeNekBlkMatSharedPtr > &TracePntJacCons, Array< OneD, Array< OneD, TypeNekBlkMatSharedPtr > > &gmtxarray, const Array< OneD, TypeNekBlkMatSharedPtr > &TracePntJacGrad, const Array< OneD, Array< OneD, DataType > > &TracePntJacGradSign)
 
template<typename DataType , typename TypeNekBlkMatSharedPtr >
void CalcJacobTraceInteg (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int m, const int n, const Array< OneD, const TypeNekBlkMatSharedPtr > &PntJac, const Array< OneD, const Array< OneD, DataType > > &PntJacSign, Array< OneD, DNekMatSharedPtr > &TraceJacFwd, Array< OneD, DNekMatSharedPtr > &TraceJacBwd)
 
template<typename DataType , typename TypeNekBlkMatSharedPtr >
void AddTraceJacToMat (const int nConvectiveFields, const int nSpaceDim, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, TypeNekBlkMatSharedPtr > &TracePntJacCons, Array< OneD, Array< OneD, TypeNekBlkMatSharedPtr > > &gmtxarray, const Array< OneD, TypeNekBlkMatSharedPtr > &TracePntJacGrad, const Array< OneD, Array< OneD, DataType > > &TracePntJacGradSign)
 

Static Public Member Functions

static SolverUtils::AdvectionSharedPtr create (std::string)
 Creates an instance of this class. More...
 

Static Public Attributes

static std::string className
 Name of class. More...
 
static std::string className2
 

Protected Member Functions

 NavierStokesAdvection ()
 
 ~NavierStokesAdvection () override
 
void v_InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields) override
 Initialises the advection object. More...
 
void v_Advect (const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &advVel, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble &time, const Array< OneD, Array< OneD, NekDouble > > &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayOfArray) override
 Advects a vector field. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::Advection
virtual SOLVER_UTILS_EXPORT void v_InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
 Initialises the advection object. More...
 
virtual SOLVER_UTILS_EXPORT void v_Advect (const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &advVel, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble &time, const Array< OneD, Array< OneD, NekDouble > > &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayOfArray)=0
 Advects a vector field. More...
 
virtual SOLVER_UTILS_EXPORT void v_SetBaseFlow (const Array< OneD, Array< OneD, NekDouble > > &inarray, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields)
 Overrides the base flow used during linearised advection. More...
 

Static Protected Attributes

static std::string navierStokesAdvectionTypeLookupIds []
 

Private Attributes

bool m_specHP_dealiasing
 
bool m_homogen_dealiasing
 
bool m_SingleMode
 
bool m_HalfMode
 

Friends

class MemoryManager< NavierStokesAdvection >
 

Additional Inherited Members

- Protected Attributes inherited from Nektar::SolverUtils::Advection
AdvectionFluxVecCB m_fluxVector
 Callback function to the flux vector (set when advection is in conservative form). More...
 
RiemannSolverSharedPtr m_riemann
 Riemann solver for DG-type schemes. More...
 
int m_spaceDim
 Storage for space dimension. Used for homogeneous extension. More...
 

Detailed Description

Definition at line 43 of file NavierStokesAdvection.h.

Constructor & Destructor Documentation

◆ NavierStokesAdvection()

Nektar::NavierStokesAdvection::NavierStokesAdvection ( )
protected

Constructor. Creates ...

Parameters

param

Definition at line 63 of file NavierStokesAdvection.cpp.

63 : Advection()
64
65{
66}

◆ ~NavierStokesAdvection()

Nektar::NavierStokesAdvection::~NavierStokesAdvection ( )
overrideprotected

Definition at line 68 of file NavierStokesAdvection.cpp.

69{
70}

Member Function Documentation

◆ create()

static SolverUtils::AdvectionSharedPtr Nektar::NavierStokesAdvection::create ( std::string  )
inlinestatic

Creates an instance of this class.

Definition at line 50 of file NavierStokesAdvection.h.

51 {
53 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

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

◆ SetSpecHPDealiasing()

void Nektar::NavierStokesAdvection::SetSpecHPDealiasing ( bool  value)
inline

Definition at line 59 of file NavierStokesAdvection.h.

References m_specHP_dealiasing.

◆ v_Advect()

void Nektar::NavierStokesAdvection::v_Advect ( const int  nConvectiveFields,
const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble > > &  advVel,
const Array< OneD, Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble time,
const Array< OneD, Array< OneD, NekDouble > > &  pFwd = NullNekDoubleArrayOfArray,
const Array< OneD, Array< OneD, NekDouble > > &  pBwd = NullNekDoubleArrayOfArray 
)
overrideprotectedvirtual

Advects a vector field.

Implements Nektar::SolverUtils::Advection.

Definition at line 86 of file NavierStokesAdvection.cpp.

95{
96 int nqtot = fields[0]->GetTotPoints();
97 ASSERTL1(nConvectiveFields == inarray.size(),
98 "Number of convective fields and Inarray are not compatible");
99
100 // use dimension of Velocity vector to dictate dimension of operation
101 int ndim = advVel.size();
102 Array<OneD, Array<OneD, NekDouble>> AdvVel(advVel.size());
103
104 Array<OneD, Array<OneD, NekDouble>> velocity(ndim);
105
106 LibUtilities::Timer timer;
107
108 for (int i = 0; i < ndim; ++i)
109 {
110 if (fields[i]->GetWaveSpace() && !m_SingleMode && !m_HalfMode &&
112 {
113 velocity[i] = Array<OneD, NekDouble>(nqtot, 0.0);
114 fields[i]->HomogeneousBwdTrans(nqtot, advVel[i], velocity[i]);
115 }
116 else
117 {
118 velocity[i] = advVel[i];
119 }
120 }
121
122 int nPointsTot = fields[0]->GetNpoints();
123 Array<OneD, NekDouble> grad0, grad1, grad2, wkSp;
124
125 NekDouble OneDptscale = 1.5; // factor to rescale 1d points in dealiasing
126
128 {
129 // Get number of points to dealias a quadratic non-linearity
130 nPointsTot = fields[0]->Get1DScaledTotPoints(OneDptscale);
131 }
132
133 // interpolate Advection velocity
134 if (m_specHP_dealiasing) // interpolate advection field to higher space.
135 {
136 for (int i = 0; i < ndim; ++i)
137 {
138 AdvVel[i] = Array<OneD, NekDouble>(nPointsTot);
139 // interpolate infield to 3/2 dimension
140 timer.Start();
141 fields[0]->PhysInterp1DScaled(OneDptscale, velocity[i], AdvVel[i]);
142 timer.Stop();
143 timer.AccumulateRegion("Interp1DScaled");
144 }
145 }
146 else
147 {
148 for (int i = 0; i < ndim; ++i)
149 {
150 AdvVel[i] = velocity[i];
151 }
152 }
153
154 wkSp = Array<OneD, NekDouble>(nPointsTot);
155
156 // Evaluate V\cdot Grad(u)
157 switch (ndim)
158 {
159 case 1:
160 grad0 = Array<OneD, NekDouble>(fields[0]->GetNpoints());
161 for (int n = 0; n < nConvectiveFields; ++n)
162 {
163 fields[0]->PhysDeriv(inarray[n], grad0);
164 if (m_specHP_dealiasing) // interpolate gradient field
165 {
166 Array<OneD, NekDouble> Outarray(nPointsTot);
167 fields[0]->PhysInterp1DScaled(OneDptscale, grad0, wkSp);
168 Vmath::Vmul(nPointsTot, wkSp, 1, AdvVel[0], 1, Outarray, 1);
169 // Galerkin project solution back to origianl spac
170 timer.Start();
171 fields[0]->PhysGalerkinProjection1DScaled(
172 OneDptscale, Outarray, outarray[n]);
173 timer.Stop();
174 timer.AccumulateRegion("GalerinProject");
175 }
176 else
177 {
178 Vmath::Vmul(nPointsTot, grad0, 1, AdvVel[0], 1, outarray[n],
179 1);
180 }
181 }
182 break;
183 case 2:
184 grad0 = Array<OneD, NekDouble>(nqtot);
185 grad1 = Array<OneD, NekDouble>(nqtot);
186 for (int n = 0; n < nConvectiveFields; ++n)
187 {
188 fields[0]->PhysDeriv(inarray[n], grad0, grad1);
189
190 if (m_specHP_dealiasing) // interpolate gradient field
191 {
192 Array<OneD, NekDouble> Outarray(nPointsTot);
193 fields[0]->PhysInterp1DScaled(OneDptscale, grad0, wkSp);
194 Vmath::Vmul(nPointsTot, wkSp, 1, AdvVel[0], 1, Outarray, 1);
195 timer.Start();
196 fields[0]->PhysInterp1DScaled(OneDptscale, grad1, wkSp);
197 timer.Stop();
198 timer.AccumulateRegion("Interp1DScaled");
199 Vmath::Vvtvp(nPointsTot, wkSp, 1, AdvVel[1], 1, Outarray, 1,
200 Outarray, 1);
201 // Galerkin project solution back to original space
202 timer.Start();
203 fields[0]->PhysGalerkinProjection1DScaled(
204 OneDptscale, Outarray, outarray[n]);
205 timer.Stop();
206 timer.AccumulateRegion("GalerinProject");
207 }
208 else
209 {
210 Vmath::Vmul(nPointsTot, grad0, 1, AdvVel[0], 1, outarray[n],
211 1);
212 Vmath::Vvtvp(nPointsTot, grad1, 1, AdvVel[1], 1,
213 outarray[n], 1, outarray[n], 1);
214 }
215 }
216 break;
217 case 3:
218 if (m_homogen_dealiasing == true && m_specHP_dealiasing == true)
219 {
220 Array<OneD, Array<OneD, NekDouble>> grad(ndim);
221 Array<OneD, Array<OneD, NekDouble>> gradScaled(
222 ndim * nConvectiveFields);
223 Array<OneD, Array<OneD, NekDouble>> Outarray(nConvectiveFields);
224 for (int i = 0; i < ndim; i++)
225 {
226 grad[i] = Array<OneD, NekDouble>(nqtot);
227 }
228 for (int i = 0; i < ndim * nConvectiveFields; i++)
229 {
230 gradScaled[i] = Array<OneD, NekDouble>(nPointsTot);
231 }
232 for (int i = 0; i < nConvectiveFields; i++)
233 {
234 Outarray[i] = Array<OneD, NekDouble>(nPointsTot);
235 }
236
237 for (int n = 0; n < nConvectiveFields; n++)
238 {
239 fields[0]->PhysDeriv(inarray[n], grad[0], grad[1], grad[2]);
240 for (int i = 0; i < ndim; i++)
241 {
242 timer.Start();
243 fields[0]->PhysInterp1DScaled(OneDptscale, grad[i],
244 gradScaled[n * ndim + i]);
245 timer.Stop();
246 timer.AccumulateRegion("Interp1DScaled");
247 }
248 }
249
250 fields[0]->DealiasedDotProd(nPointsTot, AdvVel, gradScaled,
251 Outarray);
252
253 timer.Start();
254 for (int n = 0; n < nConvectiveFields; n++)
255 {
256 fields[0]->PhysGalerkinProjection1DScaled(
257 OneDptscale, Outarray[n], outarray[n]);
258 }
259 timer.Stop();
260 timer.AccumulateRegion("GalerinProject");
261 }
262 else if (m_homogen_dealiasing == true &&
263 m_specHP_dealiasing == false)
264 {
265 Array<OneD, Array<OneD, NekDouble>> grad(ndim *
266 nConvectiveFields);
267 Array<OneD, Array<OneD, NekDouble>> Outarray(nConvectiveFields);
268 for (int i = 0; i < ndim * nConvectiveFields; i++)
269 {
270 grad[i] = Array<OneD, NekDouble>(nPointsTot);
271 }
272 for (int i = 0; i < nConvectiveFields; i++)
273 {
274 Outarray[i] = Array<OneD, NekDouble>(nPointsTot);
275 }
276
277 for (int n = 0; n < nConvectiveFields; n++)
278 {
279 fields[0]->PhysDeriv(inarray[n], grad[n * ndim + 0],
280 grad[n * ndim + 1],
281 grad[n * ndim + 2]);
282 }
283
284 fields[0]->DealiasedDotProd(nPointsTot, AdvVel, grad, outarray);
285 }
286 else
287 {
288 grad0 = Array<OneD, NekDouble>(nqtot);
289 grad1 = Array<OneD, NekDouble>(nqtot);
290 grad2 = Array<OneD, NekDouble>(nqtot);
291 Array<OneD, NekDouble> tmp = grad2;
292 for (int n = 0; n < nConvectiveFields; ++n)
293 {
294 if (fields[0]->GetWaveSpace() == true &&
295 fields[0]->GetExpType() == MultiRegions::e3DH1D)
296 {
297 fields[0]->HomogeneousBwdTrans(nqtot, inarray[n], tmp);
298 fields[0]->PhysDeriv(tmp, grad0, grad1);
299 // Take d/dz derivative using wave space field
300 fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[2],
301 inarray[n], outarray[n]);
302 fields[0]->HomogeneousBwdTrans(nqtot, outarray[n],
303 grad2);
304 }
305 else if (fields[0]->GetWaveSpace() == true &&
306 fields[0]->GetExpType() == MultiRegions::e3DH2D)
307 {
308 fields[0]->HomogeneousBwdTrans(nqtot, inarray[n], tmp);
309 fields[0]->PhysDeriv(tmp, grad0);
310 // Take d/dy derivative using wave space field
311 fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[1],
312 inarray[n], outarray[n]);
313 fields[0]->HomogeneousBwdTrans(nqtot, outarray[n],
314 grad1);
315 // Take d/dz derivative using wave space field
316 fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[2],
317 inarray[n], outarray[n]);
318 fields[0]->HomogeneousBwdTrans(nqtot, outarray[n],
319 grad2);
320 }
321 else
322 {
323 fields[0]->PhysDeriv(inarray[n], grad0, grad1, grad2);
324 }
325 if (m_specHP_dealiasing) // interpolate spectral/hp gradient
326 // field
327 {
328 Array<OneD, NekDouble> Outarray(nPointsTot);
329 timer.Start();
330 fields[0]->PhysInterp1DScaled(OneDptscale, grad0, wkSp);
331 timer.Stop();
332 timer.AccumulateRegion("Interp1DScaled");
333 Vmath::Vmul(nPointsTot, wkSp, 1, AdvVel[0], 1, Outarray,
334 1);
335
336 timer.Start();
337 fields[0]->PhysInterp1DScaled(OneDptscale, grad1, wkSp);
338 timer.Stop();
339 timer.AccumulateRegion("Interp1DScaled");
340 Vmath::Vvtvp(nPointsTot, wkSp, 1, AdvVel[1], 1,
341 Outarray, 1, Outarray, 1);
342
343 timer.Start();
344 fields[0]->PhysInterp1DScaled(OneDptscale, grad2, wkSp);
345 timer.Stop();
346 timer.AccumulateRegion("Interp1DScaled");
347 Vmath::Vvtvp(nPointsTot, wkSp, 1, AdvVel[2], 1,
348 Outarray, 1, Outarray, 1);
349 timer.Start();
350 fields[0]->PhysGalerkinProjection1DScaled(
351 OneDptscale, Outarray, outarray[n]);
352 timer.Stop();
353 timer.AccumulateRegion("GalerinProject");
354 }
355 else
356 {
357 Vmath::Vmul(nPointsTot, grad0, 1, AdvVel[0], 1,
358 outarray[n], 1);
359 Vmath::Vvtvp(nPointsTot, grad1, 1, AdvVel[1], 1,
360 outarray[n], 1, outarray[n], 1);
361 Vmath::Vvtvp(nPointsTot, grad2, 1, AdvVel[2], 1,
362 outarray[n], 1, outarray[n], 1);
363 }
364
365 if (fields[0]->GetWaveSpace() == true)
366 {
367 fields[0]->HomogeneousFwdTrans(nqtot, outarray[n],
368 outarray[n]);
369 }
370 }
371 }
372 break;
373 default:
374 ASSERTL0(false, "dimension unknown");
375 }
376
377 for (int n = 0; n < nConvectiveFields; ++n)
378 {
379 Vmath::Neg(nqtot, outarray[n], 1);
380 }
381}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:87
double NekDouble
void Vmul(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:72
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.hpp:292
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

References Nektar::LibUtilities::Timer::AccumulateRegion(), ASSERTL0, ASSERTL1, Nektar::MultiRegions::DirCartesianMap, Nektar::MultiRegions::e3DH1D, Nektar::MultiRegions::e3DH2D, m_HalfMode, m_homogen_dealiasing, m_SingleMode, m_specHP_dealiasing, Vmath::Neg(), Nektar::LibUtilities::Timer::Start(), Nektar::LibUtilities::Timer::Stop(), Vmath::Vmul(), and Vmath::Vvtvp().

◆ v_InitObject()

void Nektar::NavierStokesAdvection::v_InitObject ( LibUtilities::SessionReaderSharedPtr  pSession,
Array< OneD, MultiRegions::ExpListSharedPtr pFields 
)
overrideprotectedvirtual

Initialises the advection object.

This function should be overridden in derived classes to initialise the specific advection data members. However, this base class function should be called as the first statement of the overridden function to ensure the base class is correctly initialised in order.

Parameters
pSessionSession information.
pFieldsExpansion lists for scalar fields.

Reimplemented from Nektar::SolverUtils::Advection.

Definition at line 72 of file NavierStokesAdvection.cpp.

75{
76 m_homogen_dealiasing = pSession->DefinesSolverInfo("dealiasing");
77
78 pSession->MatchSolverInfo("SPECTRALHPDEALIASING", "True",
79 m_specHP_dealiasing, false);
80 pSession->MatchSolverInfo("ModeType", "SingleMode", m_SingleMode, false);
81 pSession->MatchSolverInfo("ModeType", "HalfMode", m_HalfMode, false);
82
83 Advection::v_InitObject(pSession, pFields);
84}
virtual SOLVER_UTILS_EXPORT void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Initialises the advection object.
Definition: Advection.cpp:295

References m_HalfMode, m_homogen_dealiasing, m_SingleMode, m_specHP_dealiasing, and Nektar::SolverUtils::Advection::v_InitObject().

Friends And Related Function Documentation

◆ MemoryManager< NavierStokesAdvection >

friend class MemoryManager< NavierStokesAdvection >
friend

Definition at line 1 of file NavierStokesAdvection.h.

Member Data Documentation

◆ className

string Nektar::NavierStokesAdvection::className
static
Initial value:
=
"Convective", NavierStokesAdvection::create, "Convective")
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
static SolverUtils::AdvectionSharedPtr create(std::string)
Creates an instance of this class.
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:43

Name of class.

Definition at line 56 of file NavierStokesAdvection.h.

◆ className2

string Nektar::NavierStokesAdvection::className2
static
Initial value:

Definition at line 57 of file NavierStokesAdvection.h.

◆ m_HalfMode

bool Nektar::NavierStokesAdvection::m_HalfMode
private

Definition at line 90 of file NavierStokesAdvection.h.

Referenced by v_Advect(), and v_InitObject().

◆ m_homogen_dealiasing

bool Nektar::NavierStokesAdvection::m_homogen_dealiasing
private

Definition at line 88 of file NavierStokesAdvection.h.

Referenced by v_Advect(), and v_InitObject().

◆ m_SingleMode

bool Nektar::NavierStokesAdvection::m_SingleMode
private

Definition at line 89 of file NavierStokesAdvection.h.

Referenced by v_Advect(), and v_InitObject().

◆ m_specHP_dealiasing

bool Nektar::NavierStokesAdvection::m_specHP_dealiasing
private

Definition at line 87 of file NavierStokesAdvection.h.

Referenced by SetSpecHPDealiasing(), v_Advect(), and v_InitObject().

◆ navierStokesAdvectionTypeLookupIds

std::string Nektar::NavierStokesAdvection::navierStokesAdvectionTypeLookupIds
staticprotected
Initial value:
= {
"True", 0),
"False", 1)}
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.

Definition at line 84 of file NavierStokesAdvection.h.