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_AdvectVolumeFlux (const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &advVel, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &pVolumeFlux, const NekDouble &time)
 Advects Volume Flux. More...
 
virtual SOLVER_UTILS_EXPORT void v_AdvectTraceFlux (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 >> &pTraceFlux, const NekDouble &time, const Array< OneD, Array< OneD, NekDouble > > &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayOfArray)
 Advects Trace Flux. More...
 
virtual SOLVER_UTILS_EXPORT void v_AdvectCoeffs (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)
 
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...
 
virtual SOLVER_UTILS_EXPORT void v_AddVolumJacToMat (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int &nConvectiveFields, const TensorOfArray5D< NekDouble > &ElmtJacArray, Array< OneD, Array< OneD, SNekBlkMatSharedPtr >> &gmtxarray)
 

Private Attributes

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...
 
SOLVER_UTILS_EXPORT void AdvectVolumeFlux (const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble >> &pAdvVel, const Array< OneD, Array< OneD, NekDouble >> &pInarray, TensorOfArray3D< NekDouble > &pVolumeFlux, const NekDouble &pTime)
 Interface function to advect the Volume field. More...
 
SOLVER_UTILS_EXPORT void AdvectTraceFlux (const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, Array< OneD, NekDouble >> &pAdvVel, const Array< OneD, Array< OneD, NekDouble >> &pInarray, Array< OneD, Array< OneD, NekDouble >> &pTraceFlux, const NekDouble &pTime, const Array< OneD, Array< OneD, NekDouble >> &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble >> &pBwd=NullNekDoubleArrayOfArray)
 Interface function to advect the Trace field. More...
 
SOLVER_UTILS_EXPORT void AdvectCoeffs (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)
 Similar with Advection::Advect(): calculate the advection flux The difference is in the outarray: it is the coefficients of basis for AdvectCoeffs() it is the physics on quadrature points for Advect() 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...
 
template<typename DataType , typename TypeNekBlkMatSharedPtr >
SOLVER_UTILS_EXPORT void AddTraceJacToMat (const int nConvectiveFields, const int nSpaceDim, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, TypeNekBlkMatSharedPtr > &TracePntJacCons, Array< OneD, Array< OneD, TypeNekBlkMatSharedPtr > > &gmtxarray, const Array< OneD, TypeNekBlkMatSharedPtr > &TracePntJacGrad, const Array< OneD, Array< OneD, DataType > > &TracePntJacGradSign)
 
template<typename DataType , typename TypeNekBlkMatSharedPtr >
void CalcJacobTraceInteg (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int m, const int n, const Array< OneD, const TypeNekBlkMatSharedPtr > &PntJac, const Array< OneD, const Array< OneD, DataType > > &PntJacSign, Array< OneD, DNekMatSharedPtr > &TraceJacFwd, Array< OneD, DNekMatSharedPtr > &TraceJacBwd)
 
SOLVER_UTILS_EXPORT void AddVolumJacToMat (const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int &nConvectiveFields, const TensorOfArray5D< NekDouble > &ElmtJacArray, Array< OneD, Array< OneD, SNekBlkMatSharedPtr >> &gmtxarray)
 
template<typename DataType , typename TypeNekBlkMatSharedPtr >
void AddTraceJacToMat (const int nConvectiveFields, const int nSpaceDim, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const Array< OneD, TypeNekBlkMatSharedPtr > &TracePntJacCons, Array< OneD, Array< OneD, TypeNekBlkMatSharedPtr > > &gmtxarray, const Array< OneD, TypeNekBlkMatSharedPtr > &TracePntJacGrad, const Array< OneD, Array< OneD, DataType > > &TracePntJacGradSign)
 
- 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.

51  {
53  }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

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

◆ 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 82 of file SkewSymmetricAdvection.cpp.

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

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

◆ 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.

70 {
71  Advection::v_InitObject(pSession, pFields);
72 
73  m_homogen_dealiasing = pSession->DefinesSolverInfo("dealiasing");
74  pSession->MatchSolverInfo("ModeType","SingleMode",m_SingleMode,false);
75  pSession->MatchSolverInfo("ModeType","HalfMode",m_HalfMode,false);
76 }
virtual SOLVER_UTILS_EXPORT void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Initialises the advection object.
Definition: Advection.cpp:366

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

Friends And Related Function Documentation

◆ MemoryManager< SkewSymmetricAdvection >

friend class MemoryManager< SkewSymmetricAdvection >
friend

Definition at line 1 of file SkewSymmetricAdvection.h.

Member Data Documentation

◆ className

string Nektar::SkewSymmetricAdvection::className
static
Initial value:
"SkewSymmetric",
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:200
static SolverUtils::AdvectionSharedPtr create(std::string)
Creates an instance of this class.
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47

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_HalfMode

bool Nektar::SkewSymmetricAdvection::m_HalfMode
private

Definition at line 82 of file SkewSymmetricAdvection.h.

Referenced by v_Advect(), and v_InitObject().

◆ m_homogen_dealiasing

bool Nektar::SkewSymmetricAdvection::m_homogen_dealiasing
private

Definition at line 80 of file SkewSymmetricAdvection.h.

Referenced by v_Advect(), and v_InitObject().

◆ m_SingleMode

bool Nektar::SkewSymmetricAdvection::m_SingleMode
private

Definition at line 81 of file SkewSymmetricAdvection.h.

Referenced by v_Advect(), and v_InitObject().