Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Private Attributes | Friends | List of all members
Nektar::NavierStokesAdvection Class Reference

#include <NavierStokesAdvection.h>

Inheritance diagram for Nektar::NavierStokesAdvection:
Inheritance graph
[legend]
Collaboration diagram for Nektar::NavierStokesAdvection:
Collaboration graph
[legend]

Public Member Functions

void SetSpecHPDealiasing (bool value)
 
- Public Member Functions inherited from Nektar::SolverUtils::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)
 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)
 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)
 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)
 Overrides the base flow used during linearised advection. More...
 

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

Constructor & Destructor Documentation

Nektar::NavierStokesAdvection::NavierStokesAdvection ( )
protected

Constructor. Creates ...

Parameters
param

Definition at line 52 of file NavierStokesAdvection.cpp.

52  :
53  Advection()
54 
55  {
56 
57  }
Nektar::NavierStokesAdvection::~NavierStokesAdvection ( )
protectedvirtual

Definition at line 59 of file NavierStokesAdvection.cpp.

60  {
61  }

Member Function Documentation

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

Creates an instance of this class.

Definition at line 52 of file NavierStokesAdvection.h.

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

52  {
54  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
void Nektar::NavierStokesAdvection::SetSpecHPDealiasing ( bool  value)
inline

Definition at line 60 of file NavierStokesAdvection.h.

References m_specHP_dealiasing.

61  {
62  m_specHP_dealiasing = value;
63  }
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 
)
protectedvirtual

Advects a vector field.

Implements Nektar::SolverUtils::Advection.

Definition at line 83 of file NavierStokesAdvection.cpp.

References ASSERTL0, ASSERTL1, Nektar::MultiRegions::DirCartesianMap, m_CoeffState, m_homogen_dealiasing, m_specHP_dealiasing, Vmath::Neg(), Vmath::Vadd(), Vmath::Vcopy(), Vmath::Vmul(), and Vmath::Vvtvp().

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

67  {
69  m_homogen_dealiasing = pSession->DefinesSolverInfo("dealiasing");
70 
71  pSession->MatchSolverInfo("SPECTRALHPDEALIASING","True",m_specHP_dealiasing,false);
72  if(m_specHP_dealiasing == false)
73  {
74  pSession->MatchSolverInfo("SPECTRALHPDEALIASING","On",m_specHP_dealiasing,false);
75  }
76  pSession->MatchSolverInfo("ModeType","SingleMode",m_SingleMode,false);
77  pSession->MatchSolverInfo("ModeType","HalfMode",m_HalfMode,false);
78 
79  Advection::v_InitObject(pSession, pFields);
80  }
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:97

Friends And Related Function Documentation

friend class MemoryManager< NavierStokesAdvection >
friend

Definition at line 49 of file NavierStokesAdvection.h.

Member Data Documentation

string Nektar::NavierStokesAdvection::className = SolverUtils::GetAdvectionFactory().RegisterCreatorFunction("Convective", NavierStokesAdvection::create)
static

Name of class.

Definition at line 57 of file NavierStokesAdvection.h.

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

Definition at line 58 of file NavierStokesAdvection.h.

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

Definition at line 84 of file NavierStokesAdvection.h.

Referenced by v_Advect(), and v_InitObject().

bool Nektar::NavierStokesAdvection::m_HalfMode
private

Definition at line 88 of file NavierStokesAdvection.h.

Referenced by v_InitObject().

bool Nektar::NavierStokesAdvection::m_homogen_dealiasing
private

Definition at line 86 of file NavierStokesAdvection.h.

Referenced by v_Advect(), and v_InitObject().

bool Nektar::NavierStokesAdvection::m_SingleMode
private

Definition at line 87 of file NavierStokesAdvection.h.

Referenced by v_InitObject().

bool Nektar::NavierStokesAdvection::m_specHP_dealiasing
private

Definition at line 85 of file NavierStokesAdvection.h.

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