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

48  : Advection()
49 
50 {
51 }

◆ ~SkewSymmetricAdvection()

Nektar::SkewSymmetricAdvection::~SkewSymmetricAdvection ( )
protectedvirtual

Definition at line 56 of file SkewSymmetricAdvection.cpp.

57 {
58 }

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.

52  {
54  }
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 77 of file SkewSymmetricAdvection.cpp.

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

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

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

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", SkewSymmetricAdvection::create, "Skew Symmetric")
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
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 57 of file SkewSymmetricAdvection.h.

◆ className2

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

Definition at line 58 of file SkewSymmetricAdvection.h.

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