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

#include <NavierStokesAdvection.h>

Inheritance diagram for Nektar::NavierStokesAdvection:
Inheritance graph
[legend]
Collaboration diagram for Nektar::NavierStokesAdvection:
Collaboration graph
[legend]

Public Member Functions

void SetSpecHPDealiasing (bool value)
- 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.

Static Public Member Functions

static
SolverUtils::AdvectionSharedPtr 
create (std::string)
 Creates an instance of this class.

Static Public Attributes

static std::string className = SolverUtils::GetAdvectionFactory().RegisterCreatorFunction("Convective", NavierStokesAdvection::create)
 Name of class.
static std::string className2 = SolverUtils::GetAdvectionFactory().RegisterCreatorFunction("NonConservative", NavierStokesAdvection::create)

Protected Member Functions

 NavierStokesAdvection ()
virtual ~NavierStokesAdvection ()
virtual void v_InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
 Initialises the advection object.
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.
- 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.

Private Attributes

MultiRegions::CoeffState m_CoeffState
bool m_specHP_dealiasing
bool m_homogen_dealiasing
bool m_SingleMode
bool m_HalfMode

Friends

class MemoryManager< NavierStokesAdvection >

Additional Inherited Members

- 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 NavierStokesAdvection.h.

Constructor & Destructor Documentation

Nektar::NavierStokesAdvection::NavierStokesAdvection ( )
protected

Constructor. Creates ...

Parameters
\param

Definition at line 50 of file NavierStokesAdvection.cpp.

:
Advection()
{
}
Nektar::NavierStokesAdvection::~NavierStokesAdvection ( )
protectedvirtual

Definition at line 57 of file NavierStokesAdvection.cpp.

{
}

Member Function Documentation

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

Creates an instance of this class.

Definition at line 52 of file NavierStokesAdvection.h.

void Nektar::NavierStokesAdvection::SetSpecHPDealiasing ( bool  value)
inline

Definition at line 60 of file NavierStokesAdvection.h.

References m_specHP_dealiasing.

{
}
void Nektar::NavierStokesAdvection::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 81 of file NavierStokesAdvection.cpp.

References ASSERTL0, ASSERTL1, Nektar::MultiRegions::DirCartesianMap, m_CoeffState, m_homogen_dealiasing, m_specHP_dealiasing, Vmath::Neg(), Vmath::Vadd(), Vmath::Vcopy(), Vmath::Vmul(), and Vmath::Vvtvp().

{
int nqtot = fields[0]->GetTotPoints();
ASSERTL1(nConvectiveFields == inarray.num_elements(),"Number of convective fields and Inarray are not compatible");
Array<OneD, NekDouble > Deriv(nqtot*nConvectiveFields);
for(int n = 0; n < nConvectiveFields; ++n)
{
// use dimension of Velocity vector to dictate dimension of operation
int ndim = advVel.num_elements();
Array<OneD, Array<OneD, NekDouble> > AdvVel (advVel.num_elements());
Array<OneD, NekDouble> Outarray;
int nPointsTot = fields[0]->GetNpoints();
Array<OneD, NekDouble> grad0,grad1,grad2,wkSp;
NekDouble OneDptscale = 1.5; // factor to rescale 1d points in dealiasing
{
// Get number of points to dealias a quadratic non-linearity
nPointsTot = fields[0]->Get1DScaledTotPoints(OneDptscale);
}
grad0 = Array<OneD, NekDouble> (nPointsTot);
// interpolate Advection velocity
int nadv = advVel.num_elements();
if(m_specHP_dealiasing) // interpolate advection field to higher space.
{
AdvVel[0] = Array<OneD, NekDouble> (nPointsTot*(nadv+1));
for(int i = 0; i < nadv; ++i)
{
if(i)
{
AdvVel[i] = AdvVel[i-1]+nPointsTot;
}
// interpolate infield to 3/2 dimension
fields[0]->PhysInterp1DScaled(OneDptscale,advVel[i],AdvVel[i]);
}
Outarray = AdvVel[nadv-1] + nPointsTot;
}
else
{
for(int i = 0; i < nadv; ++i)
{
AdvVel[i] = advVel[i];
}
Outarray = outarray[n];
}
wkSp = Array<OneD, NekDouble> (nPointsTot);
// Evaluate V\cdot Grad(u)
switch(ndim)
{
case 1:
fields[0]->PhysDeriv(inarray[n],grad0);
Vmath::Vmul(nPointsTot,grad0,1,advVel[0],1,outarray[n],1);
break;
case 2:
{
grad1 = Array<OneD, NekDouble> (nPointsTot);
fields[0]->PhysDeriv(inarray[n],grad0,grad1);
if(m_specHP_dealiasing) // interpolate gradient field
{
fields[0]->PhysInterp1DScaled(OneDptscale,grad0,wkSp);
Vmath::Vcopy(nPointsTot,wkSp,1,grad0,1);
fields[0]->PhysInterp1DScaled(OneDptscale,grad1,wkSp);
Vmath::Vcopy(nPointsTot,wkSp,1,grad1,1);
}
Vmath::Vmul (nPointsTot,grad0,1,AdvVel[0],1,Outarray,1);
Vmath::Vvtvp(nPointsTot,grad1,1,AdvVel[1],1,Outarray,1,Outarray,1);
if(m_specHP_dealiasing) // Galerkin project solution back to origianl space
{
fields[0]->PhysGalerkinProjection1DScaled(OneDptscale,Outarray,outarray[n]);
}
}
break;
case 3:
grad1 = Array<OneD, NekDouble> (fields[0]->GetNpoints());
grad2 = Array<OneD, NekDouble> (fields[0]->GetNpoints());
if(fields[0]->GetWaveSpace() == false && m_homogen_dealiasing == true )
{
ASSERTL0(m_specHP_dealiasing == false,"Spectral/hp element dealaising is not set up for this option");
fields[0]->PhysDeriv(inarray[n],grad0,grad1,grad2);
fields[0]->DealiasedProd(advVel[0],grad0,grad0,m_CoeffState);
fields[0]->DealiasedProd(advVel[1],grad1,grad1,m_CoeffState);
fields[0]->DealiasedProd(advVel[2],grad2,grad2,m_CoeffState);
Vmath::Vadd(nPointsTot,grad0,1,grad1,1,outarray[n],1);
Vmath::Vadd(nPointsTot,grad2,1,outarray[n],1,outarray[n],1);
}
else if(fields[0]->GetWaveSpace() == true && m_homogen_dealiasing == false)
{
// take d/dx, d/dy gradients in physical Fourier space
fields[0]->PhysDeriv(advVel[n],grad0,grad1);
// Take d/dz derivative using wave space field
fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[2],inarray[n],
outarray[n]);
fields[0]->HomogeneousBwdTrans(outarray[n],grad2);
if(m_specHP_dealiasing) //interpolate spectral/hp gradient field
{
fields[0]->PhysInterp1DScaled(OneDptscale,grad0,wkSp);
Vmath::Vmul(nPointsTot,wkSp,1,AdvVel[0],1,Outarray,1);
}
else
{
Vmath::Vmul(nPointsTot,grad0,1,AdvVel[0],1,Outarray,1);
}
if(m_specHP_dealiasing) //interpolate spectral/hp gradient field
{
fields[0]->PhysInterp1DScaled(OneDptscale,grad1,wkSp);
Vmath::Vvtvp(nPointsTot,wkSp,1,AdvVel[1],1,Outarray,1,
Outarray,1);
}
else
{
Vmath::Vvtvp(nPointsTot,grad1,1,AdvVel[1],1,Outarray,1,
Outarray,1);
}
if(m_specHP_dealiasing) //interpolate spectral/hp gradient field
{
fields[0]->PhysInterp1DScaled(OneDptscale,grad2,wkSp);
Vmath::Vvtvp(nPointsTot,wkSp,1,AdvVel[2],1,Outarray,1,Outarray,1);
fields[0]->PhysGalerkinProjection1DScaled(OneDptscale,Outarray,grad2);
fields[0]->HomogeneousFwdTrans(grad2,outarray[n]);
}
else
{
Vmath::Vvtvp(nPointsTot,grad2,1,AdvVel[2],1,Outarray,1,grad0,1);
fields[0]->HomogeneousFwdTrans(grad0,outarray[n]);
}
}
else if(fields[0]->GetWaveSpace() == false && m_homogen_dealiasing == false)
{
fields[0]->PhysDeriv(inarray[n],grad0,grad1,grad2);
if(m_specHP_dealiasing) //interpolate spectral/hp gradient field
{
fields[0]->PhysInterp1DScaled(OneDptscale,grad0,wkSp);
Vmath::Vmul(nPointsTot,wkSp,1,AdvVel[0],1,Outarray,1);
}
else
{
Vmath::Vmul(nPointsTot,grad0,1,AdvVel[0],1,Outarray,1);
}
if(m_specHP_dealiasing) //interpolate spectral/hp gradient field
{
fields[0]->PhysInterp1DScaled(OneDptscale,grad1,wkSp);
Vmath::Vvtvp(nPointsTot,wkSp,1,AdvVel[1],1,Outarray,1,
Outarray,1);
}
else
{
Vmath::Vvtvp(nPointsTot,grad1,1,AdvVel[1],1,Outarray,1,
Outarray,1);
}
if(m_specHP_dealiasing) //interpolate spectral/hp gradient field
{
fields[0]->PhysInterp1DScaled(OneDptscale,grad2,wkSp);
Vmath::Vvtvp(nPointsTot,wkSp,1,AdvVel[2],1,Outarray,1,Outarray,1);
fields[0]->PhysGalerkinProjection1DScaled(OneDptscale,Outarray,outarray[n]);
}
else
{
Vmath::Vvtvp(nPointsTot,grad2,1,AdvVel[2],1,Outarray,1,outarray[n],1);
}
}
else if(fields[0]->GetWaveSpace() == true && m_homogen_dealiasing == true)
{
ASSERTL0(m_specHP_dealiasing == false,"Spectral/hp element dealaising is not set up for this option");
fields[0]->PhysDeriv(inarray[n],grad0,grad1,grad2);
fields[0]->HomogeneousBwdTrans(grad0, outarray[n]);
fields[0]->DealiasedProd(advVel[0], outarray[n], grad0,
fields[0]->HomogeneousBwdTrans(grad1,outarray[n]);
fields[0]->DealiasedProd(advVel[1], outarray[n], grad1,
fields[0]->HomogeneousBwdTrans(grad2,outarray[n]);
fields[0]->DealiasedProd(advVel[2], outarray[n], grad2,
Vmath::Vadd(nPointsTot, grad0, 1, grad1, 1, grad0, 1);
Vmath::Vadd(nPointsTot, grad0, 1, grad2, 1, grad0, 1);
fields[0]->HomogeneousFwdTrans(grad0,outarray[n]);
}
else
{
ASSERTL0(false, "Advection term calculation not implented or "
"possible with the current problem set up");
}
break;
default:
ASSERTL0(false,"dimension unknown");
}
Vmath::Neg(nqtot,outarray[n],1);
}
}
void Nektar::NavierStokesAdvection::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 NavierStokesAdvection.cpp.

References Nektar::MultiRegions::eLocal, m_CoeffState, m_HalfMode, m_homogen_dealiasing, m_SingleMode, and m_specHP_dealiasing.

{
m_homogen_dealiasing = pSession->DefinesSolverInfo("dealiasing");
pSession->MatchSolverInfo("SPECTRALHPDEALIASING","True",m_specHP_dealiasing,false);
if(m_specHP_dealiasing == false)
{
pSession->MatchSolverInfo("SPECTRALHPDEALIASING","On",m_specHP_dealiasing,false);
}
pSession->MatchSolverInfo("ModeType","SingleMode",m_SingleMode,false);
pSession->MatchSolverInfo("ModeType","HalfMode",m_HalfMode,false);
Advection::v_InitObject(pSession, pFields);
}

Friends And Related Function Documentation

friend class MemoryManager< NavierStokesAdvection >
friend

Definition at line 49 of file NavierStokesAdvection.h.

Member Data Documentation

string Nektar::NavierStokesAdvection::className = SolverUtils::GetAdvectionFactory().RegisterCreatorFunction("Convective", NavierStokesAdvection::create)
static

Name of class.

Definition at line 57 of file NavierStokesAdvection.h.

string Nektar::NavierStokesAdvection::className2 = SolverUtils::GetAdvectionFactory().RegisterCreatorFunction("NonConservative", NavierStokesAdvection::create)
static

Definition at line 58 of file NavierStokesAdvection.h.

MultiRegions::CoeffState Nektar::NavierStokesAdvection::m_CoeffState
private

Definition at line 84 of file NavierStokesAdvection.h.

Referenced by v_Advect(), and v_InitObject().

bool Nektar::NavierStokesAdvection::m_HalfMode
private

Definition at line 88 of file NavierStokesAdvection.h.

Referenced by v_InitObject().

bool Nektar::NavierStokesAdvection::m_homogen_dealiasing
private

Definition at line 86 of file NavierStokesAdvection.h.

Referenced by v_Advect(), and v_InitObject().

bool Nektar::NavierStokesAdvection::m_SingleMode
private

Definition at line 87 of file NavierStokesAdvection.h.

Referenced by v_InitObject().

bool Nektar::NavierStokesAdvection::m_specHP_dealiasing
private

Definition at line 85 of file NavierStokesAdvection.h.

Referenced by SetSpecHPDealiasing(), v_Advect(), and v_InitObject().