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_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

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

Constructor & Destructor Documentation

◆ AlternateSkewAdvection()

Nektar::AlternateSkewAdvection::AlternateSkewAdvection ( )
protected

Constructor. Creates ...

Parameters

Definition at line 52 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.

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

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

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

References m_HalfMode, and m_SingleMode.

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

Friends And Related Function Documentation

◆ MemoryManager< AlternateSkewAdvection >

friend class MemoryManager< AlternateSkewAdvection >
friend

Definition at line 47 of file AlternateSkewAdvection.h.

Member Data Documentation

◆ className

string Nektar::AlternateSkewAdvection::className
static
Initial value:

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