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...
 
SOLVER_UTILS_EXPORT void AdvectVolumeFlux (const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble >> &pAdvVel, const Array< OneD, Array< OneD, NekDouble >> &pInarray, TensorOfArray3D< NekDouble > &pVolumeFlux, const NekDouble &pTime)
 Interface function to advect the Volume field. More...
 
SOLVER_UTILS_EXPORT void AdvectTraceFlux (const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble >> &pAdvVel, const Array< OneD, Array< OneD, NekDouble >> &pInarray, Array< OneD, Array< OneD, NekDouble >> &pTraceFlux, const NekDouble &pTime, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray)
 Interface function to advect the Trace field. More...
 
SOLVER_UTILS_EXPORT void AdvectCoeffs (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)
 Similar with Advection::Advect(): calculate the advection flux The difference is in the outarray: it is the coefficients of basis for AdvectCoeffs() it is the physics on quadrature points for Advect() 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)
 
SOLVER_UTILS_EXPORT void AddVolumJacToMat (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int &nConvectiveFields, const TensorOfArray5D< NekDouble > &ElmtJacArray, Array< OneD, Array< OneD, SNekBlkMatSharedPtr >> &gmtxarray)
 
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 = SolverUtils::GetAdvectionFactory().RegisterCreatorFunction("Convective", NavierStokesAdvection::create)
 Name of class. More...
 
static std::string className2 = SolverUtils::GetAdvectionFactory().RegisterCreatorFunction("NonConservative", NavierStokesAdvection::create)
 

Protected Member Functions

 NavierStokesAdvection ()
 
virtual ~NavierStokesAdvection ()
 
virtual void v_InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
 Initialises the advection object. More...
 
virtual 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)
 Advects a vector field. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::Advection
virtual SOLVER_UTILS_EXPORT void v_AdvectVolumeFlux (const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &advVel, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &pVolumeFlux, const NekDouble &time)
 Advects Volume Flux. More...
 
virtual SOLVER_UTILS_EXPORT void v_AdvectTraceFlux (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 >> &pTraceFlux, const NekDouble &time, const Array< OneD, Array< OneD, NekDouble > > &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayOfArray)
 Advects Trace Flux. More...
 
virtual SOLVER_UTILS_EXPORT void v_AdvectCoeffs (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)
 
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...
 
virtual SOLVER_UTILS_EXPORT void v_AddVolumJacToMat (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int &nConvectiveFields, const TensorOfArray5D< NekDouble > &ElmtJacArray, Array< OneD, Array< OneD, SNekBlkMatSharedPtr >> &gmtxarray)
 

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 44 of file NavierStokesAdvection.h.

Constructor & Destructor Documentation

◆ NavierStokesAdvection()

Nektar::NavierStokesAdvection::NavierStokesAdvection ( )
protected

Constructor. Creates ...

Parameters

Definition at line 59 of file NavierStokesAdvection.cpp.

59  :
60  Advection()
61 
62  {
63 
64  }

◆ ~NavierStokesAdvection()

Nektar::NavierStokesAdvection::~NavierStokesAdvection ( )
protectedvirtual

Definition at line 66 of file NavierStokesAdvection.cpp.

67  {
68  }

Member Function Documentation

◆ create()

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

Creates an instance of this class.

Definition at line 51 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.

60  {
61  m_specHP_dealiasing = value;
62  }

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

Advects a vector field.

Implements Nektar::SolverUtils::Advection.

Definition at line 85 of file NavierStokesAdvection.cpp.

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

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

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 71 of file NavierStokesAdvection.cpp.

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

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 = SolverUtils::GetAdvectionFactory().RegisterCreatorFunction("Convective", NavierStokesAdvection::create)
static

Name of class.

Definition at line 56 of file NavierStokesAdvection.h.

◆ className2

string Nektar::NavierStokesAdvection::className2 = SolverUtils::GetAdvectionFactory().RegisterCreatorFunction("NonConservative", NavierStokesAdvection::create)
static

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.