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) override
 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) override
 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 
)
overrideprotectedvirtual

Advects a vector field.

Implements Nektar::SolverUtils::Advection.

Definition at line 77 of file SkewSymmetricAdvection.cpp.

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

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