Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator 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.
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)
 Advects a vector field.
- 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.

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.
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.
template<typename FuncPointerT , typename ObjectPointerT >
void SetFluxVector (FuncPointerT func, ObjectPointerT obj)
 Set the flux vector callback function.
void SetRiemannSolver (RiemannSolverSharedPtr riemann)
 Set a Riemann solver object for this advection object.
void SetFluxVector (AdvectionFluxVecCB fluxVector)
 Set the flux vector callback function.
void SetBaseFlow (const Array< OneD, Array< OneD, NekDouble > > &inarray)
 Set the base flow used for linearised advection objects.
- Protected Attributes inherited from Nektar::SolverUtils::Advection
AdvectionFluxVecCB m_fluxVector
 Callback function to the flux vector (set when advection is in conservative form).
RiemannSolverSharedPtr m_riemann
 Riemann solver for DG-type schemes.
int m_spaceDim
 Storage for space dimension. Used for homogeneous extension.

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

{
}

Member Function Documentation

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

Definition at line 48 of file AdvectionNonConservative.h.

References AdvectionNonConservative().

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

Advects a vector field.

Implements Nektar::SolverUtils::Advection.

Definition at line 65 of file AdvectionNonConservative.cpp.

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

{
int nDim = advVel.num_elements();
int nPointsTot = fields[0]->GetNpoints();
Array<OneD, NekDouble> grad0,grad1,grad2;
grad0 = Array<OneD, NekDouble> (nPointsTot);
if (nDim > 1)
{
grad1 = Array<OneD,NekDouble>(nPointsTot);
}
if (nDim > 2)
{
grad2 = Array<OneD,NekDouble>(nPointsTot);
}
for (int i = 0; i < nConvectiveFields; ++i)
{
// Evaluate V \cdot Grad(u)
switch(nDim)
{
case 1:
fields[0]->PhysDeriv(inarray[i], grad0);
Vmath::Vmul(nPointsTot,
grad0, 1,
advVel[0], 1,
outarray[i], 1);
break;
case 2:
fields[0]->PhysDeriv(inarray[i], grad0, grad1);
// Calculate advection terms
Vmath::Vmul (nPointsTot,
grad0, 1,
advVel[0], 1,
outarray[i], 1);
Vmath::Vvtvp(nPointsTot,
grad1, 1,
advVel[1], 1,
outarray[i], 1,
outarray[i], 1);
break;
case 3:
fields[0]->PhysDeriv(inarray[i], grad0, grad1, grad2);
// Calculate advection terms
Vmath::Vmul (nPointsTot,
grad0, 1,
advVel[0], 1,
outarray[i], 1);
Vmath::Vvtvp(nPointsTot,
grad1, 1,
advVel[1], 1,
outarray[i], 1,
outarray[i], 1);
Vmath::Vvtvp(nPointsTot,
grad2, 1,
advVel[2], 1,
outarray[i], 1,
outarray[i], 1);
break;
default:
ASSERTL0(false,"dimension unknown");
}
}
}
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.

{
Advection::v_InitObject(pSession, pFields);
}

Member Data Documentation

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

Definition at line 53 of file AdvectionNonConservative.h.