Nektar++
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:
[legend]

Static Public Member Functions

static AdvectionSharedPtr create (std::string advType)
 

Static Public Attributes

static std::string type
 

Protected Member Functions

 AdvectionNonConservative ()
 
void v_InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields) override
 Initialise AdvectionNonConservative objects and store them before starting the time-stepping. More...
 
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) override
 Advects a vector field. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::Advection
virtual SOLVER_UTILS_EXPORT void v_InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
 Initialises the advection object. More...
 
virtual SOLVER_UTILS_EXPORT 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)=0
 Advects a vector field. More...
 
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
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...
 
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)
 
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 42 of file AdvectionNonConservative.h.

Constructor & Destructor Documentation

◆ AdvectionNonConservative()

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

Definition at line 44 of file AdvectionNonConservative.cpp.

45{
46}

Referenced by create().

Member Function Documentation

◆ create()

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

Definition at line 45 of file AdvectionNonConservative.h.

46 {
48 }
std::shared_ptr< Advection > AdvectionSharedPtr
A shared pointer to an Advection object.
Definition: Advection.h:54

References AdvectionNonConservative().

◆ v_Advect()

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

Advects a vector field.

Implements Nektar::SolverUtils::Advection.

Definition at line 65 of file AdvectionNonConservative.cpp.

74{
75 int nDim = advVel.size();
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 for (int i = 0; i < nConvectiveFields; ++i)
92 {
93 // Evaluate V \cdot Grad(u)
94 switch (nDim)
95 {
96 case 1:
97 fields[0]->PhysDeriv(inarray[i], grad0);
98
99 Vmath::Vmul(nPointsTot, grad0, 1, advVel[0], 1, outarray[i], 1);
100 break;
101 case 2:
102 fields[0]->PhysDeriv(inarray[i], grad0, grad1);
103
104 // Calculate advection terms
105 Vmath::Vmul(nPointsTot, grad0, 1, advVel[0], 1, outarray[i], 1);
106
107 Vmath::Vvtvp(nPointsTot, grad1, 1, advVel[1], 1, outarray[i], 1,
108 outarray[i], 1);
109
110 break;
111 case 3:
112 fields[0]->PhysDeriv(inarray[i], grad0, grad1, grad2);
113
114 // Calculate advection terms
115 Vmath::Vmul(nPointsTot, grad0, 1, advVel[0], 1, outarray[i], 1);
116
117 Vmath::Vvtvp(nPointsTot, grad1, 1, advVel[1], 1, outarray[i], 1,
118 outarray[i], 1);
119
120 Vmath::Vvtvp(nPointsTot, grad2, 1, advVel[2], 1, outarray[i], 1,
121 outarray[i], 1);
122 break;
123 default:
124 ASSERTL0(false, "dimension unknown");
125 }
126 }
127}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
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.hpp:72
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.hpp:366

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

◆ v_InitObject()

void Nektar::SolverUtils::AdvectionNonConservative::v_InitObject ( LibUtilities::SessionReaderSharedPtr  pSession,
Array< OneD, MultiRegions::ExpListSharedPtr pFields 
)
overrideprotectedvirtual

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 55 of file AdvectionNonConservative.cpp.

58{
59 Advection::v_InitObject(pSession, pFields);
60}
virtual SOLVER_UTILS_EXPORT void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
Initialises the advection object.
Definition: Advection.cpp:295

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

Member Data Documentation

◆ type

std::string Nektar::SolverUtils::AdvectionNonConservative::type
static
Initial value:
=
"Non Conservative")
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
static AdvectionSharedPtr create(std::string advType)
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:43

Definition at line 50 of file AdvectionNonConservative.h.