Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Static Public Member Functions | Static Public Attributes | Protected Member Functions | List of all members
Nektar::SolverUtils::AdvectionNonConservative Class Reference

#include <AdvectionNonConservative.h>

Inheritance diagram for Nektar::SolverUtils::AdvectionNonConservative:
Inheritance graph
[legend]
Collaboration diagram for Nektar::SolverUtils::AdvectionNonConservative:
Collaboration graph
[legend]

Static Public Member Functions

static AdvectionSharedPtr create (std::string advType)
 

Static Public Attributes

static std::string type
 

Protected Member Functions

 AdvectionNonConservative ()
 
virtual void v_InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
 Initialise AdvectionNonConservative objects and store them before starting the time-stepping. More...
 
virtual void v_Advect (const int nConvective, 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...
 

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, 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 45 of file AdvectionNonConservative.h.

Constructor & Destructor Documentation

Nektar::SolverUtils::AdvectionNonConservative::AdvectionNonConservative ( )
protected

Definition at line 46 of file AdvectionNonConservative.cpp.

Referenced by create().

47  {
48 
49  }

Member Function Documentation

static AdvectionSharedPtr Nektar::SolverUtils::AdvectionNonConservative::create ( std::string  advType)
inlinestatic

Definition at line 48 of file AdvectionNonConservative.h.

References AdvectionNonConservative().

49  {
51  }
boost::shared_ptr< Advection > AdvectionSharedPtr
A shared pointer to an Advection object.
Definition: Advection.h:165
void Nektar::SolverUtils::AdvectionNonConservative::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 65 of file AdvectionNonConservative.cpp.

References ASSERTL0, Vmath::Vmul(), and Vmath::Vvtvp().

74  {
75  int nDim = advVel.num_elements();
76  int nPointsTot = fields[0]->GetNpoints();
77  Array<OneD, NekDouble> grad0,grad1,grad2;
78 
79  grad0 = Array<OneD, NekDouble> (nPointsTot);
80 
81  if (nDim > 1)
82  {
83  grad1 = Array<OneD,NekDouble>(nPointsTot);
84  }
85 
86  if (nDim > 2)
87  {
88  grad2 = Array<OneD,NekDouble>(nPointsTot);
89  }
90 
91 
92  for (int i = 0; i < nConvectiveFields; ++i)
93  {
94  // Evaluate V \cdot Grad(u)
95  switch(nDim)
96  {
97  case 1:
98  fields[0]->PhysDeriv(inarray[i], grad0);
99 
100  Vmath::Vmul(nPointsTot,
101  grad0, 1,
102  advVel[0], 1,
103  outarray[i], 1);
104  break;
105  case 2:
106  fields[0]->PhysDeriv(inarray[i], grad0, grad1);
107 
108 
109  // Calculate advection terms
110  Vmath::Vmul (nPointsTot,
111  grad0, 1,
112  advVel[0], 1,
113  outarray[i], 1);
114 
115  Vmath::Vvtvp(nPointsTot,
116  grad1, 1,
117  advVel[1], 1,
118  outarray[i], 1,
119  outarray[i], 1);
120 
121  break;
122  case 3:
123  fields[0]->PhysDeriv(inarray[i], grad0, grad1, grad2);
124 
125  // Calculate advection terms
126  Vmath::Vmul (nPointsTot,
127  grad0, 1,
128  advVel[0], 1,
129  outarray[i], 1);
130 
131  Vmath::Vvtvp(nPointsTot,
132  grad1, 1,
133  advVel[1], 1,
134  outarray[i], 1,
135  outarray[i], 1);
136 
137  Vmath::Vvtvp(nPointsTot,
138  grad2, 1,
139  advVel[2], 1,
140  outarray[i], 1,
141  outarray[i], 1);
142  break;
143  default:
144  ASSERTL0(false,"dimension unknown");
145  }
146  }
147  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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:442
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:183
void Nektar::SolverUtils::AdvectionNonConservative::v_InitObject ( LibUtilities::SessionReaderSharedPtr  pSession,
Array< OneD, MultiRegions::ExpListSharedPtr pFields 
)
protectedvirtual

Initialise AdvectionNonConservative objects and store them before starting the time-stepping.

Parameters
pSessionPointer to session reader.
pFieldsPointer to fields.

Reimplemented from Nektar::SolverUtils::Advection.

Definition at line 58 of file AdvectionNonConservative.cpp.

References Nektar::SolverUtils::Advection::v_InitObject().

61  {
62  Advection::v_InitObject(pSession, pFields);
63  }
virtual SOLVER_UTILS_EXPORT void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Initialises the advection object.
Definition: Advection.cpp:100

Member Data Documentation

std::string Nektar::SolverUtils::AdvectionNonConservative::type
static
Initial value:
RegisterCreatorFunction("NonConservative",

Definition at line 53 of file AdvectionNonConservative.h.