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

50  :
51  Advection()
52 
53  {
54 
55  }
Nektar::NavierStokesAdvection::~NavierStokesAdvection ( )
protectedvirtual

Definition at line 57 of file NavierStokesAdvection.cpp.

58  {
59  }

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 81 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().

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

References Nektar::MultiRegions::eLocal, m_CoeffState, m_HalfMode, m_homogen_dealiasing, m_SingleMode, and m_specHP_dealiasing.

65  {
67  m_homogen_dealiasing = pSession->DefinesSolverInfo("dealiasing");
68 
69  pSession->MatchSolverInfo("SPECTRALHPDEALIASING","True",m_specHP_dealiasing,false);
70  if(m_specHP_dealiasing == false)
71  {
72  pSession->MatchSolverInfo("SPECTRALHPDEALIASING","On",m_specHP_dealiasing,false);
73  }
74  pSession->MatchSolverInfo("ModeType","SingleMode",m_SingleMode,false);
75  pSession->MatchSolverInfo("ModeType","HalfMode",m_HalfMode,false);
76 
77  Advection::v_InitObject(pSession, pFields);
78  }
Local coefficients.
MultiRegions::CoeffState m_CoeffState

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