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::AdvectionWeakDG Class Reference

#include <AdvectionWeakDG.h>

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

Static Public Member Functions

static AdvectionSharedPtr create (std::string advType)
 

Static Public Attributes

static std::string type
 

Protected Member Functions

 AdvectionWeakDG ()
 
virtual void v_InitObject (LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields)
 Initialise AdvectionWeakDG 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)
 Compute the advection term at each time-step using the Discontinuous Glaerkin approach (DG). 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...
 

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

Constructor & Destructor Documentation

Nektar::SolverUtils::AdvectionWeakDG::AdvectionWeakDG ( )
protected

Definition at line 47 of file AdvectionWeakDG.cpp.

Referenced by create().

48  {
49  }

Member Function Documentation

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

Definition at line 48 of file AdvectionWeakDG.h.

References AdvectionWeakDG().

49  {
50  return AdvectionSharedPtr(new AdvectionWeakDG());
51  }
boost::shared_ptr< Advection > AdvectionSharedPtr
A shared pointer to an Advection object.
Definition: Advection.h:158
void Nektar::SolverUtils::AdvectionWeakDG::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

Compute the advection term at each time-step using the Discontinuous Glaerkin approach (DG).

Parameters
nConvectiveFieldsNumber of fields.
fieldsPointer to fields.
advVelAdvection velocities.
inarraySolution at the previous time-step.
outarrayAdvection term to be passed at the time integration class.

Implements Nektar::SolverUtils::Advection.

Definition at line 76 of file AdvectionWeakDG.cpp.

References ASSERTL1, Nektar::SolverUtils::Advection::m_fluxVector, Nektar::SolverUtils::Advection::m_riemann, Nektar::SolverUtils::Advection::m_spaceDim, and Vmath::Neg().

83  {
84  int nPointsTot = fields[0]->GetTotPoints();
85  int nCoeffs = fields[0]->GetNcoeffs();
86  int nTracePointsTot = fields[0]->GetTrace()->GetTotPoints();
87  int i, j;
88 
89  Array<OneD, Array<OneD, NekDouble> > tmp(nConvectiveFields);
91  nConvectiveFields);
92 
93  // Allocate storage for flux vector F(u).
94  for (i = 0; i < nConvectiveFields; ++i)
95  {
96  fluxvector[i] =
98  for (j = 0; j < m_spaceDim; ++j)
99  {
100  fluxvector[i][j] = Array<OneD, NekDouble>(nPointsTot);
101  }
102  }
103 
105  "Riemann solver must be provided for AdvectionWeakDG.");
106 
107  m_fluxVector(inarray, fluxvector);
108 
109  // Get the advection part (without numerical flux)
110  for(i = 0; i < nConvectiveFields; ++i)
111  {
112  tmp[i] = Array<OneD, NekDouble>(nCoeffs, 0.0);
113 
114  fields[i]->IProductWRTDerivBase(fluxvector[i],tmp[i]);
115  }
116 
117  // Store forwards/backwards space along trace space
118  Array<OneD, Array<OneD, NekDouble> > Fwd (nConvectiveFields);
119  Array<OneD, Array<OneD, NekDouble> > Bwd (nConvectiveFields);
120  Array<OneD, Array<OneD, NekDouble> > numflux(nConvectiveFields);
121 
122  for(i = 0; i < nConvectiveFields; ++i)
123  {
124  Fwd[i] = Array<OneD, NekDouble>(nTracePointsTot, 0.0);
125  Bwd[i] = Array<OneD, NekDouble>(nTracePointsTot, 0.0);
126  numflux[i] = Array<OneD, NekDouble>(nTracePointsTot, 0.0);
127  fields[i]->GetFwdBwdTracePhys(inarray[i], Fwd[i], Bwd[i]);
128  }
129 
130  m_riemann->Solve(m_spaceDim, Fwd, Bwd, numflux);
131 
132  // Evaulate <\phi, \hat{F}\cdot n> - OutField[i]
133  for(i = 0; i < nConvectiveFields; ++i)
134  {
135  Vmath::Neg (nCoeffs, tmp[i], 1);
136  fields[i]->AddTraceIntegral (numflux[i], tmp[i]);
137  fields[i]->MultiplyByElmtInvMass(tmp[i], tmp[i]);
138  fields[i]->BwdTrans (tmp[i], outarray[i]);
139  }
140  }
RiemannSolverSharedPtr m_riemann
Riemann solver for DG-type schemes.
Definition: Advection.h:134
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:382
AdvectionFluxVecCB m_fluxVector
Callback function to the flux vector (set when advection is in conservative form).
Definition: Advection.h:132
int m_spaceDim
Storage for space dimension. Used for homogeneous extension.
Definition: Advection.h:136
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
void Nektar::SolverUtils::AdvectionWeakDG::v_InitObject ( LibUtilities::SessionReaderSharedPtr  pSession,
Array< OneD, MultiRegions::ExpListSharedPtr pFields 
)
protectedvirtual

Initialise AdvectionWeakDG 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 AdvectionWeakDG.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:97

Member Data Documentation

std::string Nektar::SolverUtils::AdvectionWeakDG::type
static
Initial value:
RegisterCreatorFunction("WeakDG", AdvectionWeakDG::create)

Definition at line 53 of file AdvectionWeakDG.h.