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

#include <SkewSymmetricAdvection.h>

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

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

 SkewSymmetricAdvection ()
 
virtual ~SkewSymmetricAdvection ()
 
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_homogen_dealiasing
 
bool m_SingleMode
 
bool m_HalfMode
 

Friends

class MemoryManager< SkewSymmetricAdvection >
 

Additional Inherited Members

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

Constructor & Destructor Documentation

Nektar::SkewSymmetricAdvection::SkewSymmetricAdvection ( )
protected

Definition at line 48 of file SkewSymmetricAdvection.cpp.

48  :
49  Advection()
50 
51 {
52 }
Nektar::SkewSymmetricAdvection::~SkewSymmetricAdvection ( )
protectedvirtual

Definition at line 58 of file SkewSymmetricAdvection.cpp.

59 {
60 }

Member Function Documentation

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

Creates an instance of this class.

Definition at line 52 of file SkewSymmetricAdvection.h.

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

52  {
54  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
void Nektar::SkewSymmetricAdvection::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 82 of file SkewSymmetricAdvection.cpp.

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

89 {
90  int nqtot = fields[0]->GetTotPoints();
91  ASSERTL1(nConvectiveFields == inarray.num_elements(),"Number of convective fields and Inarray are not compatible");
92 
93  Array<OneD, NekDouble > Deriv = Array<OneD, NekDouble> (nqtot*nConvectiveFields);
94 
95  for(int n = 0; n < nConvectiveFields; ++n)
96  {
97  // use dimension of Velocity vector to dictate dimension of operation
98  int ndim = advVel.num_elements();
99 
100  // ToDo: here we should add a check that V has right dimension
101 
102  int nPointsTot = fields[0]->GetNpoints();
103  Array<OneD, NekDouble> gradV0,gradV1,gradV2, tmp, Up;
104 
105  gradV0 = Array<OneD, NekDouble> (nPointsTot);
106  tmp = Array<OneD, NekDouble> (nPointsTot);
107 
108  // Evaluate V\cdot Grad(u)
109  switch(ndim)
110  {
111  case 1:
112  fields[0]->PhysDeriv(inarray[n],gradV0);
113  Vmath::Vmul(nPointsTot,gradV0,1,advVel[0],1,outarray[n],1);
114  Vmath::Vmul(nPointsTot,inarray[n],1,advVel[0],1,gradV0,1);
115  fields[0]->PhysDeriv(gradV0,tmp);
116  Vmath::Vadd(nPointsTot,tmp,1,outarray[n],1,outarray[n],1);
117  Vmath::Smul(nPointsTot,0.5,outarray[n],1,outarray[n],1);
118  break;
119  case 2:
120  gradV1 = Array<OneD, NekDouble> (nPointsTot);
121  fields[0]->PhysDeriv(inarray[n],gradV0,gradV1);
122  Vmath::Vmul (nPointsTot,gradV0,1,advVel[0],1,outarray[n],1);
123  Vmath::Vvtvp(nPointsTot,gradV1,1,advVel[1],1,outarray[n],1,outarray[n],1);
124  Vmath::Vmul(nPointsTot,inarray[n],1,advVel[0],1,gradV0,1);
125  Vmath::Vmul(nPointsTot,inarray[n],1,advVel[1],1,gradV1,1);
126  fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[0],gradV0,tmp);
127  Vmath::Vadd(nPointsTot,tmp,1,outarray[n],1,outarray[n],1);
128  fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[1],gradV1,tmp);
129  Vmath::Vadd(nPointsTot,tmp,1,outarray[n],1,outarray[n],1);
130  Vmath::Smul(nPointsTot,0.5,outarray[n],1,outarray[n],1);
131  break;
132  case 3:
133  gradV1 = Array<OneD, NekDouble> (nPointsTot);
134  gradV2 = Array<OneD, NekDouble> (nPointsTot);
135 
136  fields[0]->PhysDeriv(inarray[n],gradV0,gradV1,gradV2);
137 
138  //outarray[n] = 1/2(u*du/dx + v*du/dy + w*du/dz + duu/dx + duv/dy + duw/dz)
139 
140  if(m_homogen_dealiasing == true && fields[0]->GetWaveSpace() == false)
141  {
142  fields[0]->DealiasedProd(advVel[0],gradV0,gradV0,m_CoeffState);
143  fields[0]->DealiasedProd(advVel[1],gradV1,gradV1,m_CoeffState);
144  fields[0]->DealiasedProd(advVel[2],gradV2,gradV2,m_CoeffState);
145  Vmath::Vadd(nPointsTot,gradV0,1,gradV1,1,outarray[n],1);
146  Vmath::Vadd(nPointsTot,gradV2,1,outarray[n],1,outarray[n],1);
147  fields[0]->DealiasedProd(inarray[n],advVel[0],gradV0,m_CoeffState);
148  fields[0]->DealiasedProd(inarray[n],advVel[1],gradV1,m_CoeffState);
149  fields[0]->DealiasedProd(inarray[n],advVel[2],gradV2,m_CoeffState);
150  fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[0],gradV0,tmp);
151  Vmath::Vadd(nPointsTot,tmp,1,outarray[n],1,outarray[n],1);
152  fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[1],gradV1,tmp);
153  Vmath::Vadd(nPointsTot,tmp,1,outarray[n],1,outarray[n],1);
154  fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[2],gradV2,tmp);
155  Vmath::Vadd(nPointsTot,tmp,1,outarray[n],1,outarray[n],1);
156  Vmath::Smul(nPointsTot,0.5,outarray[n],1,outarray[n],1);
157  }
158  else if(fields[0]->GetWaveSpace() == true && m_homogen_dealiasing == false)
159  {
160  Up = Array<OneD, NekDouble> (nPointsTot);
161  //vector reused to avoid even more memory requirements
162  //names may be misleading
163  fields[0]->HomogeneousBwdTrans(gradV0,tmp);
164  Vmath::Vmul(nPointsTot,tmp,1,advVel[0],1,outarray[n],1); // + u*du/dx
165  fields[0]->HomogeneousBwdTrans(gradV1,tmp);
166  Vmath::Vvtvp(nPointsTot,tmp,1,advVel[1],1,outarray[n],1,outarray[n],1);// + v*du/dy
167  fields[0]->HomogeneousBwdTrans(gradV2,tmp);
168  Vmath::Vvtvp(nPointsTot,tmp,1,advVel[2],1,outarray[n],1,outarray[n],1);// + w*du/dz
169 
170  fields[0]->HomogeneousBwdTrans(inarray[n],Up);
171  Vmath::Vmul(nPointsTot,Up,1,advVel[0],1,gradV0,1);
172  Vmath::Vmul(nPointsTot,Up,1,advVel[1],1,gradV1,1);
173  Vmath::Vmul(nPointsTot,Up,1,advVel[2],1,gradV2,1);
174 
175  fields[0]->SetWaveSpace(false);
176  fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[0],gradV0,tmp);//duu/dx
177  Vmath::Vadd(nPointsTot,tmp,1,outarray[n],1,outarray[n],1);
178  fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[1],gradV1,tmp);//duv/dy
179  Vmath::Vadd(nPointsTot,tmp,1,outarray[n],1,outarray[n],1);
180  fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[2],gradV2,tmp);//duw/dz
181  Vmath::Vadd(nPointsTot,tmp,1,outarray[n],1,outarray[n],1);
182  fields[0]->SetWaveSpace(true);
183 
184  Vmath::Smul(nPointsTot,0.5,outarray[n],1,tmp,1);
185  fields[0]->HomogeneousFwdTrans(tmp,outarray[n]);
186  }
187  else if(fields[0]->GetWaveSpace() == false && m_homogen_dealiasing == false)
188  {
189  Vmath::Vmul(nPointsTot,gradV0,1,advVel[0],1,outarray[n],1);
190  Vmath::Vvtvp(nPointsTot,gradV1,1,advVel[1],1,outarray[n],1,outarray[n],1);
191  Vmath::Vvtvp(nPointsTot,gradV2,1,advVel[2],1,outarray[n],1,outarray[n],1);
192  Vmath::Vmul(nPointsTot,inarray[n],1,advVel[0],1,gradV0,1);
193  Vmath::Vmul(nPointsTot,inarray[n],1,advVel[1],1,gradV1,1);
194  Vmath::Vmul(nPointsTot,inarray[n],1,advVel[2],1,gradV2,1);
195  fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[0],gradV0,tmp);
196  Vmath::Vadd(nPointsTot,tmp,1,outarray[n],1,outarray[n],1);
197  fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[1],gradV1,tmp);
198  Vmath::Vadd(nPointsTot,tmp,1,outarray[n],1,outarray[n],1);
199  fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[2],gradV2,tmp);
200  Vmath::Vadd(nPointsTot,tmp,1,outarray[n],1,outarray[n],1);
201  Vmath::Smul(nPointsTot,0.5,outarray[n],1,outarray[n],1);
202  }
203  else
204  {
205  ASSERTL0(false, "Dealiasing is not allowed in combination "
206  "with the Skew-Symmetric advection form for "
207  "efficiency reasons.");
208  }
209  break;
210  default:
211  ASSERTL0(false,"dimension unknown");
212  }
213 
214  Vmath::Neg(nqtot,outarray[n],1);
215  }
216 
217 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
MultiRegions::CoeffState m_CoeffState
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 Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:199
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:382
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 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::SkewSymmetricAdvection::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 66 of file SkewSymmetricAdvection.cpp.

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

69 {
70  Advection::v_InitObject(pSession, pFields);
71 
73  m_homogen_dealiasing = pSession->DefinesSolverInfo("dealiasing");
74  pSession->MatchSolverInfo("ModeType","SingleMode",m_SingleMode,false);
75  pSession->MatchSolverInfo("ModeType","HalfMode",m_HalfMode,false);
76 }
Local coefficients.
MultiRegions::CoeffState m_CoeffState

Friends And Related Function Documentation

friend class MemoryManager< SkewSymmetricAdvection >
friend

Definition at line 49 of file SkewSymmetricAdvection.h.

Member Data Documentation

string Nektar::SkewSymmetricAdvection::className
static
Initial value:

Name of class.

Definition at line 57 of file SkewSymmetricAdvection.h.

std::string Nektar::SkewSymmetricAdvection::className2
static

Definition at line 58 of file SkewSymmetricAdvection.h.

MultiRegions::CoeffState Nektar::SkewSymmetricAdvection::m_CoeffState
private

Definition at line 79 of file SkewSymmetricAdvection.h.

Referenced by v_Advect(), and v_InitObject().

bool Nektar::SkewSymmetricAdvection::m_HalfMode
private

Definition at line 82 of file SkewSymmetricAdvection.h.

Referenced by v_InitObject().

bool Nektar::SkewSymmetricAdvection::m_homogen_dealiasing
private

Definition at line 80 of file SkewSymmetricAdvection.h.

Referenced by v_Advect(), and v_InitObject().

bool Nektar::SkewSymmetricAdvection::m_SingleMode
private

Definition at line 81 of file SkewSymmetricAdvection.h.

Referenced by v_InitObject().