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

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

MultiRegions::CoeffState m_CoeffState
 
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 58 of file NavierStokesAdvection.cpp.

Referenced by SetSpecHPDealiasing().

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

◆ ~NavierStokesAdvection()

Nektar::NavierStokesAdvection::~NavierStokesAdvection ( )
protectedvirtual

Definition at line 65 of file NavierStokesAdvection.cpp.

Referenced by SetSpecHPDealiasing().

66  {
67  }

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.

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

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

◆ SetSpecHPDealiasing()

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

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

References ASSERTL0, ASSERTL1, Nektar::MultiRegions::DirCartesianMap, Nektar::MultiRegions::e3DH1D, Nektar::MultiRegions::e3DH2D, m_CoeffState, m_HalfMode, m_homogen_dealiasing, m_SingleMode, m_specHP_dealiasing, Vmath::Neg(), Vmath::Vmul(), and Vmath::Vvtvp().

Referenced by SetSpecHPDealiasing().

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

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

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

Referenced by SetSpecHPDealiasing().

73  {
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  }
Local coefficients.
MultiRegions::CoeffState m_CoeffState
virtual SOLVER_UTILS_EXPORT void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Initialises the advection object.
Definition: Advection.cpp:98

Friends And Related Function Documentation

◆ MemoryManager< NavierStokesAdvection >

friend class MemoryManager< NavierStokesAdvection >
friend

Definition at line 48 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_CoeffState

MultiRegions::CoeffState Nektar::NavierStokesAdvection::m_CoeffState
private

Definition at line 87 of file NavierStokesAdvection.h.

Referenced by v_Advect(), and v_InitObject().

◆ m_HalfMode

bool Nektar::NavierStokesAdvection::m_HalfMode
private

Definition at line 91 of file NavierStokesAdvection.h.

Referenced by v_Advect(), and v_InitObject().

◆ m_homogen_dealiasing

bool Nektar::NavierStokesAdvection::m_homogen_dealiasing
private

Definition at line 89 of file NavierStokesAdvection.h.

Referenced by v_Advect(), and v_InitObject().

◆ m_SingleMode

bool Nektar::NavierStokesAdvection::m_SingleMode
private

Definition at line 90 of file NavierStokesAdvection.h.

Referenced by v_Advect(), and v_InitObject().

◆ m_specHP_dealiasing

bool Nektar::NavierStokesAdvection::m_specHP_dealiasing
private

Definition at line 88 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)}

Definition at line 84 of file NavierStokesAdvection.h.