Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 50 of file SkewSymmetricAdvection.cpp.

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

Definition at line 60 of file SkewSymmetricAdvection.cpp.

61 {
62 }

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

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

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

71 {
72  Advection::v_InitObject(pSession, pFields);
73 
75  m_homogen_dealiasing = pSession->DefinesSolverInfo("dealiasing");
76  pSession->MatchSolverInfo("ModeType","SingleMode",m_SingleMode,false);
77  pSession->MatchSolverInfo("ModeType","HalfMode",m_HalfMode,false);
78 }
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< 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().