Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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:
Inheritance graph
[legend]
Collaboration diagram for Nektar::AlternateSkewAdvection:
Collaboration graph
[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)
 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)
 Overrides the base flow used during linearised advection. More...
 

Private Attributes

int m_advectioncalls
 

Friends

class MemoryManager< AlternateSkewAdvection >
 

Additional Inherited Members

- Public Member Functions inherited from Nektar::SolverUtils::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)
 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)
 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 44 of file AlternateSkewAdvection.h.

Constructor & Destructor Documentation

Nektar::AlternateSkewAdvection::AlternateSkewAdvection ( )
protected

Constructor. Creates ...

Parameters
param

Definition at line 53 of file AlternateSkewAdvection.cpp.

54  : Advection()
55 {
56 }
Nektar::AlternateSkewAdvection::~AlternateSkewAdvection ( )
protectedvirtual

Definition at line 58 of file AlternateSkewAdvection.cpp.

59 {
60 }

Member Function Documentation

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

Creates an instance of this class.

Definition at line 51 of file AlternateSkewAdvection.h.

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

52  {
54  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
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 
)
protectedvirtual

Advects a vector field.

Implements Nektar::SolverUtils::Advection.

Definition at line 68 of file AlternateSkewAdvection.cpp.

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

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

65 {
66 }

Friends And Related Function Documentation

friend class MemoryManager< AlternateSkewAdvection >
friend

Definition at line 48 of file AlternateSkewAdvection.h.

Member Data Documentation

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

Name of class.

Definition at line 56 of file AlternateSkewAdvection.h.

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

Definition at line 57 of file AlternateSkewAdvection.h.

int Nektar::AlternateSkewAdvection::m_advectioncalls
private

Definition at line 78 of file AlternateSkewAdvection.h.

Referenced by v_Advect().