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:
[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, 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...
 

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

Constructor & Destructor Documentation

◆ SkewSymmetricAdvection()

Nektar::SkewSymmetricAdvection::SkewSymmetricAdvection ( )
protected

Definition at line 49 of file SkewSymmetricAdvection.cpp.

49  :
50  Advection()
51 
52 {
53 }

◆ ~SkewSymmetricAdvection()

Nektar::SkewSymmetricAdvection::~SkewSymmetricAdvection ( )
protectedvirtual

Definition at line 59 of file SkewSymmetricAdvection.cpp.

60 {
61 }

Member Function Documentation

◆ create()

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

Creates an instance of this class.

Definition at line 51 of file SkewSymmetricAdvection.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.

◆ v_Advect()

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,
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 83 of file SkewSymmetricAdvection.cpp.

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

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

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

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

friend class MemoryManager< SkewSymmetricAdvection >
friend

Definition at line 48 of file SkewSymmetricAdvection.h.

Member Data Documentation

◆ className

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

Name of class.

Definition at line 56 of file SkewSymmetricAdvection.h.

◆ className2

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

Definition at line 57 of file SkewSymmetricAdvection.h.

◆ m_CoeffState

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

Definition at line 80 of file SkewSymmetricAdvection.h.

Referenced by v_Advect(), and v_InitObject().

◆ m_HalfMode

bool Nektar::SkewSymmetricAdvection::m_HalfMode
private

Definition at line 83 of file SkewSymmetricAdvection.h.

Referenced by v_Advect(), and v_InitObject().

◆ m_homogen_dealiasing

bool Nektar::SkewSymmetricAdvection::m_homogen_dealiasing
private

Definition at line 81 of file SkewSymmetricAdvection.h.

Referenced by v_Advect(), and v_InitObject().

◆ m_SingleMode

bool Nektar::SkewSymmetricAdvection::m_SingleMode
private

Definition at line 82 of file SkewSymmetricAdvection.h.

Referenced by v_Advect(), and v_InitObject().