Nektar++
Static Public Member Functions | Static Public Attributes | Protected Member Functions | Private Attributes | Friends | List of all members
Nektar::AlternateSkewAdvection Class Reference

#include <AlternateSkewAdvection.h>

Inheritance diagram for Nektar::AlternateSkewAdvection:
[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

 AlternateSkewAdvection ()
 
virtual ~AlternateSkewAdvection ()
 
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

int m_advectioncalls
 
bool m_SingleMode
 
bool m_HalfMode
 

Friends

class MemoryManager< AlternateSkewAdvection >
 

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 43 of file AlternateSkewAdvection.h.

Constructor & Destructor Documentation

◆ AlternateSkewAdvection()

Nektar::AlternateSkewAdvection::AlternateSkewAdvection ( )
protected

Constructor. Creates ...

Parameters

Definition at line 53 of file AlternateSkewAdvection.cpp.

53  : Advection()
54 {
55 }

◆ ~AlternateSkewAdvection()

Nektar::AlternateSkewAdvection::~AlternateSkewAdvection ( )
protectedvirtual

Definition at line 57 of file AlternateSkewAdvection.cpp.

58 {
59 }

Member Function Documentation

◆ create()

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

Creates an instance of this class.

Definition at line 50 of file AlternateSkewAdvection.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::AlternateSkewAdvection::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 69 of file AlternateSkewAdvection.cpp.

77 {
78  // use dimension of Velocity vector to dictate dimension of operation
79  int ndim = advVel.size();
80  int nPointsTot = fields[0]->GetNpoints();
81  Array<OneD, Array<OneD, NekDouble>> velocity(ndim);
82  for (int i = 0; i < ndim; ++i)
83  {
84  if (fields[i]->GetWaveSpace() && !m_SingleMode && !m_HalfMode)
85  {
86  velocity[i] = Array<OneD, NekDouble>(nPointsTot, 0.0);
87  fields[i]->HomogeneousBwdTrans(advVel[i], velocity[i]);
88  }
89  else
90  {
91  velocity[i] = advVel[i];
92  }
93  }
94  for (int n = 0; n < nConvectiveFields; ++n)
95  {
96  // ToDo: here we should add a check that V has right dimension
97  Array<OneD, NekDouble> gradV0, gradV1, gradV2, tmp, Up;
98 
99  gradV0 = Array<OneD, NekDouble>(nPointsTot);
100  tmp = Array<OneD, NekDouble>(nPointsTot);
101 
102  // Evaluate V\cdot Grad(u)
103  switch (ndim)
104  {
105  case 1:
106  if (m_advectioncalls % 2 == 0)
107  {
108  fields[0]->PhysDeriv(inarray[n], gradV0);
109  Vmath::Vmul(nPointsTot, gradV0, 1, velocity[0], 1,
110  outarray[n], 1);
111  }
112  else
113  {
114  Vmath::Vmul(nPointsTot, inarray[n], 1, velocity[0], 1,
115  gradV0, 1);
116  fields[0]->PhysDeriv(gradV0, outarray[n]);
117  }
118  Vmath::Smul(nPointsTot, 0.5, outarray[n], 1, outarray[n],
119  1); // must be mult by 0.5????
120  break;
121  case 2:
122  gradV1 = Array<OneD, NekDouble>(nPointsTot);
123  if (m_advectioncalls % 2 == 0)
124  {
125  fields[0]->PhysDeriv(inarray[n], gradV0, gradV1);
126  Vmath::Vmul(nPointsTot, gradV0, 1, velocity[0], 1,
127  outarray[n], 1);
128  Vmath::Vvtvp(nPointsTot, gradV1, 1, velocity[1], 1,
129  outarray[n], 1, outarray[n], 1);
130  }
131  else
132  {
133  Vmath::Vmul(nPointsTot, inarray[n], 1, velocity[0], 1,
134  gradV0, 1);
135  Vmath::Vmul(nPointsTot, inarray[n], 1, velocity[1], 1,
136  gradV1, 1);
137  fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[0],
138  gradV0, outarray[n]);
139  fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[1],
140  gradV1, tmp);
141  Vmath::Vadd(nPointsTot, tmp, 1, outarray[n], 1, outarray[n],
142  1);
143  }
144  Vmath::Smul(nPointsTot, 1.0, outarray[n], 1, outarray[n],
145  1); // must be mult by 0.5????
146  break;
147  case 3:
148  gradV1 = Array<OneD, NekDouble>(nPointsTot);
149  gradV2 = Array<OneD, NekDouble>(nPointsTot);
150 
151  // outarray[n] = 1/2(u*du/dx + v*du/dy + w*du/dz + duu/dx +
152  // duv/dy + duw/dz)
153 
154  if (fields[0]->GetWaveSpace() == true)
155  {
156  if (m_advectioncalls % 2 == 0)
157  {
158  // vector reused to avoid even more memory requirements
159  // names may be misleading
160  fields[0]->PhysDeriv(inarray[n], gradV0, gradV1,
161  gradV2);
162  fields[0]->HomogeneousBwdTrans(gradV0, tmp);
163  Vmath::Vmul(nPointsTot, tmp, 1, velocity[0], 1,
164  outarray[n], 1); // + u*du/dx
165  fields[0]->HomogeneousBwdTrans(gradV1, tmp);
166  Vmath::Vvtvp(nPointsTot, tmp, 1, velocity[1], 1,
167  outarray[n], 1, outarray[n],
168  1); // + v*du/dy
169  fields[0]->HomogeneousBwdTrans(gradV2, tmp);
170  Vmath::Vvtvp(nPointsTot, tmp, 1, velocity[2], 1,
171  outarray[n], 1, outarray[n],
172  1); // + w*du/dz
173  }
174  else
175  {
176  Up = Array<OneD, NekDouble>(nPointsTot);
177  fields[0]->HomogeneousBwdTrans(inarray[n], Up);
178  Vmath::Vmul(nPointsTot, Up, 1, velocity[0], 1, gradV0,
179  1);
180  Vmath::Vmul(nPointsTot, Up, 1, velocity[1], 1, gradV1,
181  1);
182  Vmath::Vmul(nPointsTot, Up, 1, velocity[2], 1, gradV2,
183  1);
184 
185  fields[0]->SetWaveSpace(false);
186  fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[0],
187  gradV0, outarray[n]); // duu/dx
188  fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[1],
189  gradV1, tmp); // duv/dy
190  Vmath::Vadd(nPointsTot, tmp, 1, outarray[n], 1,
191  outarray[n], 1);
192  fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[2],
193  gradV2, tmp); // duw/dz
194  Vmath::Vadd(nPointsTot, tmp, 1, outarray[n], 1,
195  outarray[n], 1);
196  fields[0]->SetWaveSpace(true);
197  }
198 
199  Vmath::Smul(nPointsTot, 1.0, outarray[n], 1, tmp,
200  1); // must be mult by 0.5????
201  fields[0]->HomogeneousFwdTrans(tmp, outarray[n]);
202  }
203  else
204  {
205  if (m_advectioncalls % 2 == 0)
206  {
207  fields[0]->PhysDeriv(inarray[n], gradV0, gradV1,
208  gradV2);
209  Vmath::Vmul(nPointsTot, gradV0, 1, velocity[0], 1,
210  outarray[n], 1);
211  Vmath::Vvtvp(nPointsTot, gradV1, 1, velocity[1], 1,
212  outarray[n], 1, outarray[n], 1);
213  Vmath::Vvtvp(nPointsTot, gradV2, 1, velocity[2], 1,
214  outarray[n], 1, outarray[n], 1);
215  }
216  else
217  {
218  Vmath::Vmul(nPointsTot, inarray[n], 1, velocity[0], 1,
219  gradV0, 1);
220  Vmath::Vmul(nPointsTot, inarray[n], 1, velocity[1], 1,
221  gradV1, 1);
222  Vmath::Vmul(nPointsTot, inarray[n], 1, velocity[2], 1,
223  gradV2, 1);
224  fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[0],
225  gradV0, outarray[n]);
226  fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[1],
227  gradV1, tmp);
228  Vmath::Vadd(nPointsTot, tmp, 1, outarray[n], 1,
229  outarray[n], 1);
230  fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[2],
231  gradV2, tmp);
232  Vmath::Vadd(nPointsTot, tmp, 1, outarray[n], 1,
233  outarray[n], 1);
234  }
235  Vmath::Smul(nPointsTot, 1.0, outarray[n], 1, outarray[n],
236  1); // must be mult by 0.5????
237  }
238  break;
239  default:
240  ASSERTL0(false, "dimension unknown");
241  }
242  }
243 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
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 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, Nektar::MultiRegions::DirCartesianMap, m_advectioncalls, m_HalfMode, m_SingleMode, Vmath::Smul(), Vmath::Vadd(), Vmath::Vmul(), and Vmath::Vvtvp().

◆ v_InitObject()

void Nektar::AlternateSkewAdvection::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 61 of file AlternateSkewAdvection.cpp.

64 {
65  pSession->MatchSolverInfo("ModeType", "SingleMode", m_SingleMode, false);
66  pSession->MatchSolverInfo("ModeType", "HalfMode", m_HalfMode, false);
67 }

References m_HalfMode, and m_SingleMode.

Friends And Related Function Documentation

◆ MemoryManager< AlternateSkewAdvection >

friend class MemoryManager< AlternateSkewAdvection >
friend

Definition at line 1 of file AlternateSkewAdvection.h.

Member Data Documentation

◆ className

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

Name of class.

Definition at line 55 of file AlternateSkewAdvection.h.

◆ className2

std::string Nektar::AlternateSkewAdvection::className2
static

Definition at line 56 of file AlternateSkewAdvection.h.

◆ m_advectioncalls

int Nektar::AlternateSkewAdvection::m_advectioncalls
private

Definition at line 79 of file AlternateSkewAdvection.h.

Referenced by v_Advect().

◆ m_HalfMode

bool Nektar::AlternateSkewAdvection::m_HalfMode
private

Definition at line 81 of file AlternateSkewAdvection.h.

Referenced by v_Advect(), and v_InitObject().

◆ m_SingleMode

bool Nektar::AlternateSkewAdvection::m_SingleMode
private

Definition at line 80 of file AlternateSkewAdvection.h.

Referenced by v_Advect(), and v_InitObject().