Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | Private Attributes | Friends | List of all members
Nektar::UnsteadyAdvection Class Reference

#include <UnsteadyAdvection.h>

Inheritance diagram for Nektar::UnsteadyAdvection:
[legend]

Public Member Functions

virtual ~UnsteadyAdvection ()
 Destructor. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::AdvectionSystem
SOLVER_UTILS_EXPORT AdvectionSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
virtual SOLVER_UTILS_EXPORT ~AdvectionSystem ()
 
virtual SOLVER_UTILS_EXPORT void v_InitObject (bool DeclareField=true) override
 Init object for UnsteadySystem class. More...
 
SOLVER_UTILS_EXPORT AdvectionSharedPtr GetAdvObject ()
 Returns the advection object held by this instance. More...
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleGetElmtCFLVals (const bool FlagAcousticCFL=true)
 
SOLVER_UTILS_EXPORT NekDouble GetCFLEstimate (int &elmtid)
 
- Public Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
virtual SOLVER_UTILS_EXPORT ~UnsteadySystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT NekDouble GetTimeStep (const Array< OneD, const Array< OneD, NekDouble > > &inarray)
 Calculate the larger time-step mantaining the problem stable. More...
 
SOLVER_UTILS_EXPORT void SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
SOLVER_UTILS_EXPORT LibUtilities::TimeIntegrationSchemeSharedPtrGetTimeIntegrationScheme ()
 Returns the time integration scheme. More...
 
SOLVER_UTILS_EXPORT LibUtilities::TimeIntegrationSchemeOperatorsGetTimeIntegrationSchemeOperators ()
 Returns the time integration scheme operators. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT void InitObject (bool DeclareField=true)
 Initialises the members of this object. More...
 
SOLVER_UTILS_EXPORT void DoInitialise (bool dumpInitialConditions=true)
 Perform any initialisation necessary before solving the problem. More...
 
SOLVER_UTILS_EXPORT void DoSolve ()
 Solve the problem. More...
 
SOLVER_UTILS_EXPORT void TransCoeffToPhys ()
 Transform from coefficient to physical space. More...
 
SOLVER_UTILS_EXPORT void TransPhysToCoeff ()
 Transform from physical to coefficient space. More...
 
SOLVER_UTILS_EXPORT void Output ()
 Perform output operations after solve. More...
 
SOLVER_UTILS_EXPORT std::string GetSessionName ()
 Get Session name. More...
 
template<class T >
std::shared_ptr< T > as ()
 
SOLVER_UTILS_EXPORT void ResetSessionName (std::string newname)
 Reset Session name. More...
 
SOLVER_UTILS_EXPORT LibUtilities::SessionReaderSharedPtr GetSession ()
 Get Session name. More...
 
SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr GetPressure ()
 Get pressure field if available. More...
 
SOLVER_UTILS_EXPORT void ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 
SOLVER_UTILS_EXPORT void PrintSummary (std::ostream &out)
 Print a summary of parameters and solver characteristics. More...
 
SOLVER_UTILS_EXPORT void SetLambda (NekDouble lambda)
 Set parameter m_lambda. More...
 
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction (std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
 Get a SessionFunction by name. More...
 
SOLVER_UTILS_EXPORT void SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 Initialise the data in the dependent fields. More...
 
SOLVER_UTILS_EXPORT void EvaluateExactSolution (int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 Evaluates an exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln, bool Normalised=false)
 Compute the L2 error between fields and a given exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, bool Normalised=false)
 Compute the L2 error of the fields. More...
 
SOLVER_UTILS_EXPORT NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation. More...
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleErrorExtraPoints (unsigned int field)
 Compute error (L2 and L_inf) over an larger set of quadrature points return [L2 Linf]. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n)
 Write checkpoint file of m_fields. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 Write checkpoint file of custom data fields. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_BaseFlow (const int n)
 Write base flow file of m_fields. More...
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname)
 Write field data to the given filename. More...
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 Write input fields to the given filename. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 Input field data from the given file. More...
 
SOLVER_UTILS_EXPORT void ImportFldToMultiDomains (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int ndomains)
 Input field data from the given file to multiple domains. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, std::vector< std::string > &fieldStr, Array< OneD, Array< OneD, NekDouble > > &coeffs)
 Output a field. Input field data into array from the given file. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, MultiRegions::ExpListSharedPtr &pField, std::string &pFieldName)
 Output a field. Input field data into ExpList from the given file. More...
 
SOLVER_UTILS_EXPORT void SessionSummary (SummaryList &vSummary)
 Write out a session summary. More...
 
SOLVER_UTILS_EXPORT Array< OneD, MultiRegions::ExpListSharedPtr > & UpdateFields ()
 
SOLVER_UTILS_EXPORT LibUtilities::FieldMetaDataMapUpdateFieldMetaDataMap ()
 Get hold of FieldInfoMap so it can be updated. More...
 
SOLVER_UTILS_EXPORT NekDouble GetFinalTime ()
 Return final time. More...
 
SOLVER_UTILS_EXPORT int GetNcoeffs ()
 
SOLVER_UTILS_EXPORT int GetNcoeffs (const int eid)
 
SOLVER_UTILS_EXPORT int GetNumExpModes ()
 
SOLVER_UTILS_EXPORT const Array< OneD, int > GetNumExpModesPerExp ()
 
SOLVER_UTILS_EXPORT int GetNvariables ()
 
SOLVER_UTILS_EXPORT const std::string GetVariable (unsigned int i)
 
SOLVER_UTILS_EXPORT int GetTraceTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTraceNpoints ()
 
SOLVER_UTILS_EXPORT int GetExpSize ()
 
SOLVER_UTILS_EXPORT int GetPhys_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetCoeff_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTotPoints (int n)
 
SOLVER_UTILS_EXPORT int GetNpoints ()
 
SOLVER_UTILS_EXPORT int GetSteps ()
 
SOLVER_UTILS_EXPORT NekDouble GetTimeStep ()
 
SOLVER_UTILS_EXPORT void CopyFromPhysField (const int i, Array< OneD, NekDouble > &output)
 
SOLVER_UTILS_EXPORT void CopyToPhysField (const int i, const Array< OneD, const NekDouble > &input)
 
SOLVER_UTILS_EXPORT void SetSteps (const int steps)
 
SOLVER_UTILS_EXPORT void ZeroPhysFields ()
 
SOLVER_UTILS_EXPORT void FwdTransFields ()
 
SOLVER_UTILS_EXPORT void SetModifiedBasis (const bool modbasis)
 
SOLVER_UTILS_EXPORT int GetCheckpointNumber ()
 
SOLVER_UTILS_EXPORT void SetCheckpointNumber (int num)
 
SOLVER_UTILS_EXPORT int GetCheckpointSteps ()
 
SOLVER_UTILS_EXPORT void SetCheckpointSteps (int num)
 
SOLVER_UTILS_EXPORT int GetInfoSteps ()
 
SOLVER_UTILS_EXPORT void SetInfoSteps (int num)
 
SOLVER_UTILS_EXPORT void SetIterationNumberPIT (int num)
 
SOLVER_UTILS_EXPORT void SetWindowNumberPIT (int num)
 
SOLVER_UTILS_EXPORT Array< OneD, const Array< OneD, NekDouble > > GetTraceNormals ()
 
SOLVER_UTILS_EXPORT void SetTime (const NekDouble time)
 
SOLVER_UTILS_EXPORT void SetTimeStep (const NekDouble timestep)
 
SOLVER_UTILS_EXPORT void SetInitialStep (const int step)
 
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time. More...
 
SOLVER_UTILS_EXPORT bool NegatedOp ()
 Identify if operator is negated in DoSolve. More...
 
SOLVER_UTILS_EXPORT bool ParallelInTime ()
 Check if solver use Parallel-in-Time. More...
 

Static Public Member Functions

static SolverUtils::EquationSystemSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Creates an instance of this class. More...
 

Static Public Attributes

static std::string className
 Name of class. More...
 
- Static Public Attributes inherited from Nektar::SolverUtils::UnsteadySystem
static std::string cmdSetStartTime
 
static std::string cmdSetStartChkNum
 

Protected Member Functions

 UnsteadyAdvection (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Session reader. More...
 
void GetFluxVector (const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
 Evaluate the flux at each solution point. More...
 
void GetFluxVectorDeAlias (const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
 Evaluate the flux at each solution point using dealiasing. More...
 
void DoOdeRhs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 Compute the RHS. More...
 
void DoOdeProjection (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 Compute the projection. More...
 
Array< OneD, NekDouble > & GetNormalVelocity ()
 Get the normal velocity. More...
 
virtual void v_InitObject (bool DeclareFields=true) override
 Initialise the object. More...
 
virtual void v_GenerateSummary (SolverUtils::SummaryList &s) override
 Print Summary. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::AdvectionSystem
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step) override
 
virtual SOLVER_UTILS_EXPORT Array< OneD, NekDoublev_GetMaxStdVelocity (const NekDouble SpeedSoundFactor=1.0)
 
- Protected Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
SOLVER_UTILS_EXPORT UnsteadySystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises UnsteadySystem class members. More...
 
virtual SOLVER_UTILS_EXPORT void v_InitObject (bool DeclareField=true) override
 Init object for UnsteadySystem class. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve () override
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_PrintStatusInformation (const int step, const NekDouble cpuTime)
 Print Status Information. More...
 
virtual SOLVER_UTILS_EXPORT void v_PrintSummaryStatistics (const NekDouble intTime)
 Print Summary Statistics. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise (bool dumpInitialConditions=true) override
 Sets up initial conditions. More...
 
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary (SummaryList &s) override
 Print a summary of time stepping parameters. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_GetTimeStep (const Array< OneD, const Array< OneD, NekDouble > > &inarray)
 Return the timestep to be used for the next step in the time-marching loop. More...
 
virtual SOLVER_UTILS_EXPORT bool v_PreIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT bool v_RequireFwdTrans ()
 
virtual SOLVER_UTILS_EXPORT void v_SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
virtual SOLVER_UTILS_EXPORT bool v_UpdateTimeStepCheck ()
 
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control. More...
 
SOLVER_UTILS_EXPORT void CheckForRestartTime (NekDouble &time, int &nchk)
 
SOLVER_UTILS_EXPORT void SVVVarDiffCoeff (const Array< OneD, Array< OneD, NekDouble > > vel, StdRegions::VarCoeffMap &varCoeffMap)
 Evaluate the SVV diffusion coefficient according to Moura's paper where it should proportional to h time velocity. More...
 
SOLVER_UTILS_EXPORT void DoDummyProjection (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 Perform dummy projection. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises EquationSystem class members. More...
 
virtual SOLVER_UTILS_EXPORT void v_InitObject (bool DeclareFeld=true)
 Initialisation object for EquationSystem. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise (bool dumpInitialConditions=true)
 Virtual function for initialisation implementation. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve ()
 Virtual function for solve implementation. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Virtual function for the L_inf error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false)
 Virtual function for the L_2 error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransCoeffToPhys ()
 Virtual function for transformation to physical space. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransPhysToCoeff ()
 Virtual function for transformation to coefficient space. More...
 
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary (SummaryList &l)
 Virtual function for generating summary information. More...
 
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 
virtual SOLVER_UTILS_EXPORT void v_EvaluateExactSolution (unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 
virtual SOLVER_UTILS_EXPORT void v_Output (void)
 
virtual SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr v_GetPressure (void)
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp (void)
 Virtual function to identify if operator is negated in DoSolve. More...
 
virtual SOLVER_UTILS_EXPORT void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 

Protected Attributes

bool m_useGJPStabilisation
 
NekDouble m_GJPJumpScale
 
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
 
Array< OneD, Array< OneD, NekDouble > > m_velocity
 Advection velocity. More...
 
Array< OneD, NekDoublem_traceVn
 
int m_planeNumber
 
std::vector< SolverUtils::ForcingSharedPtrm_forcing
 Forcing terms. More...
 
- Protected Attributes inherited from Nektar::SolverUtils::AdvectionSystem
SolverUtils::AdvectionSharedPtr m_advObject
 Advection term. More...
 
- Protected Attributes inherited from Nektar::SolverUtils::UnsteadySystem
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
Array< OneD, Array< OneD, NekDouble > > m_previousSolution
 Storage for previous solution for steady-state check. More...
 
std::vector< int > m_intVariables
 
NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1). More...
 
NekDouble m_CFLGrowth
 CFL growth rate. More...
 
NekDouble m_CFLEnd
 Maximun cfl in cfl growth. More...
 
int m_abortSteps
 Number of steps between checks for abort conditions. More...
 
bool m_explicitDiffusion
 Indicates if explicit or implicit treatment of diffusion is used. More...
 
bool m_explicitAdvection
 Indicates if explicit or implicit treatment of advection is used. More...
 
bool m_explicitReaction
 Indicates if explicit or implicit treatment of reaction is used. More...
 
int m_steadyStateSteps
 Check for steady state at step interval. More...
 
NekDouble m_steadyStateTol
 Tolerance to which steady state should be evaluated at. More...
 
int m_filtersInfosteps
 Number of time steps between outputting filters information. More...
 
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
 
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state. More...
 
std::ofstream m_errFile
 
NekDouble m_epsilon
 Diffusion coefficient. More...
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
bool m_verbose
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader. More...
 
std::map< std::string, SolverUtils::SessionFunctionSharedPtrm_sessionFunctions
 Map of known SessionFunctions. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 Array holding all dependent variables. More...
 
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
 Pointer to boundary conditions object. More...
 
SpatialDomains::MeshGraphSharedPtr m_graph
 Pointer to graph defining mesh. More...
 
std::string m_sessionName
 Name of the session. More...
 
NekDouble m_time
 Current time of simulation. More...
 
int m_initialStep
 Number of the step where the simulation should begin. More...
 
NekDouble m_fintime
 Finish time of the simulation. More...
 
NekDouble m_timestep
 Time step size. More...
 
NekDouble m_lambda
 Lambda constant in real system if one required. More...
 
NekDouble m_checktime
 Time between checkpoints. More...
 
NekDouble m_lastCheckTime
 
NekDouble m_TimeIncrementFactor
 
int m_nchk
 Number of checkpoints written so far. More...
 
int m_steps
 Number of steps to take. More...
 
int m_checksteps
 Number of steps between checkpoints. More...
 
int m_infosteps
 Number of time steps between outputting status information. More...
 
int m_iterPIT = 0
 Number of parallel-in-time time iteration. More...
 
int m_windowPIT = 0
 Index of windows for parallel-in-time time iteration. More...
 
int m_spacedim
 Spatial dimension (>= expansion dim). More...
 
int m_expdim
 Expansion dimension. More...
 
bool m_singleMode
 Flag to determine if single homogeneous mode is used. More...
 
bool m_halfMode
 Flag to determine if half homogeneous mode is used. More...
 
bool m_multipleModes
 Flag to determine if use multiple homogenenous modes are used. More...
 
bool m_useFFT
 Flag to determine if FFT is used for homogeneous transform. More...
 
bool m_homogen_dealiasing
 Flag to determine if dealiasing is used for homogeneous simulations. More...
 
bool m_specHP_dealiasing
 Flag to determine if dealisising is usde for the Spectral/hp element discretisation. More...
 
enum MultiRegions::ProjectionType m_projectionType
 Type of projection; e.g continuous or discontinuous. More...
 
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
 Array holding trace normals for DG simulations in the forwards direction. More...
 
Array< OneD, bool > m_checkIfSystemSingular
 Flag to indicate if the fields should be checked for singularity. More...
 
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
 Map to identify relevant solver info to dump in output fields. More...
 
Array< OneD, NekDoublem_movingFrameVelsxyz
 Moving frame of reference velocities. More...
 
Array< OneD, NekDoublem_movingFrameTheta
 Moving frame of reference angles with respect to the. More...
 
boost::numeric::ublas::matrix< NekDoublem_movingFrameProjMat
 Projection matrix for transformation between inertial and moving. More...
 
int m_NumQuadPointsError
 Number of Quadrature points used to work out the error. More...
 
enum HomogeneousType m_HomogeneousType
 
NekDouble m_LhomX
 physical length in X direction (if homogeneous) More...
 
NekDouble m_LhomY
 physical length in Y direction (if homogeneous) More...
 
NekDouble m_LhomZ
 physical length in Z direction (if homogeneous) More...
 
int m_npointsX
 number of points in X direction (if homogeneous) More...
 
int m_npointsY
 number of points in Y direction (if homogeneous) More...
 
int m_npointsZ
 number of points in Z direction (if homogeneous) More...
 
int m_HomoDirec
 number of homogenous directions More...
 

Private Attributes

NekDouble m_waveFreq
 

Friends

class MemoryManager< UnsteadyAdvection >
 

Additional Inherited Members

- Protected Types inherited from Nektar::SolverUtils::EquationSystem
enum  HomogeneousType { eHomogeneous1D , eHomogeneous2D , eHomogeneous3D , eNotHomogeneous }
 Parameter for homogeneous expansions. More...
 
- Static Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
static std::string equationSystemTypeLookupIds []
 
static std::string projectionTypeLookupIds []
 

Detailed Description

Definition at line 44 of file UnsteadyAdvection.h.

Constructor & Destructor Documentation

◆ ~UnsteadyAdvection()

Nektar::UnsteadyAdvection::~UnsteadyAdvection ( )
virtual

Destructor.

Unsteady linear advection equation destructor.

Definition at line 179 of file UnsteadyAdvection.cpp.

180{
181}

◆ UnsteadyAdvection()

Nektar::UnsteadyAdvection::UnsteadyAdvection ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr pGraph 
)
protected

Session reader.

Definition at line 49 of file UnsteadyAdvection.cpp.

52 : UnsteadySystem(pSession, pGraph), AdvectionSystem(pSession, pGraph)
53{
54 m_planeNumber = 0;
55}
SOLVER_UTILS_EXPORT AdvectionSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
SOLVER_UTILS_EXPORT UnsteadySystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises UnsteadySystem class members.

References m_planeNumber.

Member Function Documentation

◆ create()

static SolverUtils::EquationSystemSharedPtr Nektar::UnsteadyAdvection::create ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr pGraph 
)
inlinestatic

Creates an instance of this class.

Definition at line 50 of file UnsteadyAdvection.h.

53 {
56 pGraph);
57 p->InitObject();
58 return p;
59 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and CellMLToNektar.cellml_metadata::p.

◆ DoOdeProjection()

void Nektar::UnsteadyAdvection::DoOdeProjection ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time 
)
protected

Compute the projection.

Compute the projection for the linear advection equation.

Parameters
inarrayGiven fields.
outarrayCalculated solution.
timeTime.

Definition at line 258 of file UnsteadyAdvection.cpp.

261{
262 // Counter variable
263 int i;
264
265 // Number of fields (variables of the problem)
266 int nVariables = inarray.size();
267
268 // Set the boundary conditions
270
271 // Switch on the projection type (Discontinuous or Continuous)
272 switch (m_projectionType)
273 {
274 // Discontinuous projection
276 {
277 // Number of quadrature points
278 int nQuadraturePts = GetNpoints();
279
280 // Just copy over array
281 if (inarray != outarray)
282 {
283 for (i = 0; i < nVariables; ++i)
284 {
285 Vmath::Vcopy(nQuadraturePts, inarray[i], 1, outarray[i], 1);
286 }
287 }
288 break;
289 }
290
291 // Continuous projection
294 {
295 int ncoeffs = m_fields[0]->GetNcoeffs();
296 Array<OneD, NekDouble> coeffs(ncoeffs, 0.0);
297
298#if 0
299 for(i = 0; i < nVariables; ++i)
300 {
301 m_fields[i]->FwdTrans(inarray[i], coeffs);
302 m_fields[i]->BwdTrans_IterPerExp(coeffs, outarray[i]);
303 }
304#else
307
308 Array<OneD, NekDouble> wsp(ncoeffs);
309
310 for (i = 0; i < nVariables; ++i)
311 {
313 std::dynamic_pointer_cast<MultiRegions::ContField>(
314 m_fields[i]);
315
316 // copy inarray
317 Array<OneD, NekDouble> in = inarray[i];
318
319 m_fields[i]->IProductWRTBase(in, wsp);
320
322 {
324 cfield->GetGJPForcing();
325
328
329 if (GJPData->IsSemiImplicit())
330 {
331 mtype = StdRegions::eMassGJP;
332 }
333
334 // to set up forcing need initial guess in
335 // physical space
337
338 GJPData->Apply(inarray[i], wsp, NullNekDouble1DArray,
339 scale);
340 }
341
342 // Solve the system
343 MultiRegions::GlobalLinSysKey key(
344 mtype, cfield->GetLocalToGlobalMap(), factors);
345
346 cfield->GlobalSolve(key, wsp, coeffs, NullNekDouble1DArray);
347
348 m_fields[i]->BwdTrans(coeffs, outarray[i]);
349 }
350#endif
351 break;
352 }
353 default:
354 ASSERTL0(false, "Unknown projection scheme");
355 break;
356 }
357}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
NekDouble m_timestep
Time step size.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetNpoints()
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
std::shared_ptr< GJPStabilisation > GJPStabilisationSharedPtr
std::shared_ptr< ContField > ContFieldSharedPtr
Definition: ContField.h:270
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:408
StdRegions::ConstFactorMap factors
static Array< OneD, NekDouble > NullNekDouble1DArray
double NekDouble
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191

References ASSERTL0, Nektar::MultiRegions::eDiscontinuous, Nektar::StdRegions::eFactorGJP, Nektar::MultiRegions::eGalerkin, Nektar::StdRegions::eMass, Nektar::StdRegions::eMassGJP, Nektar::MultiRegions::eMixed_CG_Discontinuous, Nektar::VarcoeffHashingTest::factors, Nektar::SolverUtils::EquationSystem::GetNpoints(), Nektar::SolverUtils::EquationSystem::m_fields, m_GJPJumpScale, Nektar::SolverUtils::EquationSystem::m_projectionType, Nektar::SolverUtils::EquationSystem::m_timestep, m_useGJPStabilisation, Nektar::NullNekDouble1DArray, Nektar::SolverUtils::EquationSystem::SetBoundaryConditions(), and Vmath::Vcopy().

Referenced by v_InitObject().

◆ DoOdeRhs()

void Nektar::UnsteadyAdvection::DoOdeRhs ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time 
)
protected

Compute the RHS.

Compute the right-hand side for the linear advection equation.

Parameters
inarrayGiven fields.
outarrayCalculated solution.
timeTime.

Definition at line 216 of file UnsteadyAdvection.cpp.

219{
220 // Counter variable
221 int i;
222
223 // Number of fields (variables of the problem)
224 int nVariables = inarray.size();
225
226 // Number of solution points
227 int nSolutionPts = GetNpoints();
228
229 LibUtilities::Timer timer;
230 timer.Start();
231 // RHS computation using the new advection base class
232 m_advObject->Advect(nVariables, m_fields, m_velocity, inarray, outarray,
233 time);
234 timer.Stop();
235 // Elapsed time
236 timer.AccumulateRegion("Advect");
237
238 // Negate the RHS
239 for (i = 0; i < nVariables; ++i)
240 {
241 Vmath::Neg(nSolutionPts, outarray[i], 1);
242 }
243
244 for (auto &x : m_forcing)
245 {
246 // set up non-linear terms
247 x->Apply(m_fields, inarray, outarray, time);
248 }
249}
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Forcing terms.
Array< OneD, Array< OneD, NekDouble > > m_velocity
Advection velocity.
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:513

References Nektar::LibUtilities::Timer::AccumulateRegion(), Nektar::SolverUtils::EquationSystem::GetNpoints(), Nektar::SolverUtils::AdvectionSystem::m_advObject, Nektar::SolverUtils::EquationSystem::m_fields, m_forcing, m_velocity, Vmath::Neg(), Nektar::LibUtilities::Timer::Start(), and Nektar::LibUtilities::Timer::Stop().

Referenced by v_InitObject().

◆ GetFluxVector()

void Nektar::UnsteadyAdvection::GetFluxVector ( const Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &  flux 
)
protected

Evaluate the flux at each solution point.

Return the flux vector for the linear advection equation.

Parameters
iComponent of the flux vector to calculate.
physfieldFields.
fluxResulting flux.

Definition at line 366 of file UnsteadyAdvection.cpp.

369{
370 ASSERTL1(flux[0].size() == m_velocity.size(),
371 "Dimension of flux array and velocity array do not match");
372
373 int i, j;
374 int nq = physfield[0].size();
375
376 for (i = 0; i < flux.size(); ++i)
377 {
378 for (j = 0; j < flux[0].size(); ++j)
379 {
380 Vmath::Vmul(nq, physfield[i], 1, m_velocity[j], 1, flux[i][j], 1);
381 }
382 }
383}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
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:207

References ASSERTL1, m_velocity, and Vmath::Vmul().

Referenced by v_InitObject().

◆ GetFluxVectorDeAlias()

void Nektar::UnsteadyAdvection::GetFluxVectorDeAlias ( const Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &  flux 
)
protected

Evaluate the flux at each solution point using dealiasing.

Return the flux vector for the linear advection equation using the dealiasing technique.

Parameters
iComponent of the flux vector to calculate.
physfieldFields.
fluxResulting flux.

Definition at line 393 of file UnsteadyAdvection.cpp.

396{
397 ASSERTL1(flux[0].size() == m_velocity.size(),
398 "Dimension of flux array and velocity array do not match");
399
400 int i, j;
401 int nq = physfield[0].size();
402 int nVariables = physfield.size();
403
404 // Factor to rescale 1d points in dealiasing
405 NekDouble OneDptscale = 2;
406
407 Array<OneD, Array<OneD, NekDouble>> advVel_plane(m_velocity.size());
408
409 // Get number of points to dealias a cubic non-linearity
410 nq = m_fields[0]->Get1DScaledTotPoints(OneDptscale);
411
412 // Initialisation of higher-space variables
413 Array<OneD, Array<OneD, NekDouble>> physfieldInterp(nVariables);
414 Array<OneD, Array<OneD, NekDouble>> velocityInterp(m_expdim);
415 Array<OneD, Array<OneD, Array<OneD, NekDouble>>> fluxInterp(nVariables);
416
417 // Interpolation to higher space of physfield
418 for (i = 0; i < nVariables; ++i)
419 {
420 physfieldInterp[i] = Array<OneD, NekDouble>(nq);
421 fluxInterp[i] = Array<OneD, Array<OneD, NekDouble>>(m_expdim);
422 for (j = 0; j < m_expdim; ++j)
423 {
424 fluxInterp[i][j] = Array<OneD, NekDouble>(nq);
425 }
426
427 m_fields[0]->PhysInterp1DScaled(OneDptscale, physfield[i],
428 physfieldInterp[i]);
429 }
430
431 // Interpolation to higher space of velocity
432 for (j = 0; j < m_expdim; ++j)
433 {
434 velocityInterp[j] = Array<OneD, NekDouble>(nq);
435
436 m_fields[0]->PhysInterp1DScaled(OneDptscale, m_velocity[j],
437 velocityInterp[j]);
438 }
439
440 // Evaluation of flux vector in the higher space
441 for (i = 0; i < flux.size(); ++i)
442 {
443 for (j = 0; j < flux[0].size(); ++j)
444 {
445 Vmath::Vmul(nq, physfieldInterp[i], 1, velocityInterp[j], 1,
446 fluxInterp[i][j], 1);
447 }
448 }
449
450 // Galerkin project solution back to original space
451 for (i = 0; i < nVariables; ++i)
452 {
453 for (j = 0; j < m_spacedim; ++j)
454 {
455 m_fields[0]->PhysGalerkinProjection1DScaled(
456 OneDptscale, fluxInterp[i][j], flux[i][j]);
457 }
458 }
459}
int m_spacedim
Spatial dimension (>= expansion dim).
int m_expdim
Expansion dimension.

References ASSERTL1, Nektar::SolverUtils::EquationSystem::m_expdim, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_spacedim, m_velocity, and Vmath::Vmul().

Referenced by v_InitObject().

◆ GetNormalVelocity()

Array< OneD, NekDouble > & Nektar::UnsteadyAdvection::GetNormalVelocity ( )
protected

Get the normal velocity.

Get the normal velocity for the linear advection equation.

Definition at line 186 of file UnsteadyAdvection.cpp.

187{
188 // Number of trace (interface) points
189 int i;
190 int nTracePts = GetTraceNpoints();
191
192 // Auxiliary variable to compute the normal velocity
193 Array<OneD, NekDouble> tmp(nTracePts);
194
195 // Reset the normal velocity
196 Vmath::Zero(nTracePts, m_traceVn, 1);
197
198 for (i = 0; i < m_velocity.size(); ++i)
199 {
200 m_fields[0]->ExtractTracePhys(m_velocity[i], tmp);
201
202 Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, tmp, 1, m_traceVn, 1,
203 m_traceVn, 1);
204 }
205
206 return m_traceVn;
207}
SOLVER_UTILS_EXPORT int GetTraceNpoints()
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
Array< OneD, NekDouble > m_traceVn
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:569
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:487

References Nektar::SolverUtils::EquationSystem::GetTraceNpoints(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_traceNormals, m_traceVn, m_velocity, Vmath::Vvtvp(), and Vmath::Zero().

Referenced by v_InitObject().

◆ v_GenerateSummary()

void Nektar::UnsteadyAdvection::v_GenerateSummary ( SolverUtils::SummaryList s)
overrideprotectedvirtual

Print Summary.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 461 of file UnsteadyAdvection.cpp.

462{
465 {
467 s, "GJP Stab. Impl. ",
468 m_session->GetSolverInfo("GJPStabilisation"));
469 SolverUtils::AddSummaryItem(s, "GJP Stab. JumpScale", m_GJPJumpScale);
470 }
471}
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s) override
Print a summary of time stepping parameters.
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:49

References Nektar::SolverUtils::AddSummaryItem(), m_GJPJumpScale, Nektar::SolverUtils::EquationSystem::m_session, m_useGJPStabilisation, and Nektar::SolverUtils::UnsteadySystem::v_GenerateSummary().

◆ v_InitObject()

void Nektar::UnsteadyAdvection::v_InitObject ( bool  DeclareFields = true)
overrideprotectedvirtual

Initialise the object.

Initialisation object for the unsteady linear advection equation.

Reimplemented from Nektar::SolverUtils::AdvectionSystem.

Definition at line 60 of file UnsteadyAdvection.cpp.

61{
62 // Call to the initialisation object of UnsteadySystem
63 AdvectionSystem::v_InitObject(DeclareFields);
64
65 // Read the advection velocities from session file
66 m_session->LoadParameter("wavefreq", m_waveFreq, 0.0);
67
68 // check to see if it is explicity turned off
69 m_session->MatchSolverInfo("GJPStabilisation", "False",
71
72 // if GJPStabilisation set to False bool will be true and
73 // if not false so negate/revese bool
75
76 m_session->LoadParameter("GJPJumpScale", m_GJPJumpScale, 1.0);
77
78 std::vector<std::string> vel;
79 vel.push_back("Vx");
80 vel.push_back("Vy");
81 vel.push_back("Vz");
82
83 // Resize the advection velocities vector to dimension of the problem
84 vel.resize(m_spacedim);
85
86 // Store in the global variable m_velocity the advection velocities
87 m_velocity = Array<OneD, Array<OneD, NekDouble>>(m_spacedim);
88 GetFunction("AdvectionVelocity")->Evaluate(vel, m_velocity);
89
90 // Type of advection class to be used
91 switch (m_projectionType)
92 {
93 // Continuous field
96 {
97 string advName;
98 m_session->LoadSolverInfo("AdvectionType", advName,
99 "NonConservative");
101 advName, advName);
103 {
104 m_advObject->SetFluxVector(
106 }
107 else
108 {
110 this);
111 }
112 break;
113 }
114 // Discontinuous field
116 {
117 // Do not forwards transform initial condition
118 m_homoInitialFwd = false;
119
120 // Define the normal velocity fields
121 if (m_fields[0]->GetTrace())
122 {
123 m_traceVn = Array<OneD, NekDouble>(GetTraceNpoints());
124 }
125
126 string advName;
127 string riemName;
128 m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
130 advName, advName);
132 {
133 m_advObject->SetFluxVector(
135 }
136 else
137 {
139 this);
140 }
141 m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
144 riemName, m_session);
145 m_riemannSolver->SetScalar(
147
148 m_advObject->SetRiemannSolver(m_riemannSolver);
149 m_advObject->InitObject(m_session, m_fields);
150 break;
151 }
152 default:
153 {
154 ASSERTL0(false, "Unsupported projection type.");
155 break;
156 }
157 }
158
159 // Forcing terms
160 m_forcing = SolverUtils::Forcing::Load(m_session, shared_from_this(),
161 m_fields, m_fields.size());
162
163 // If explicit it computes RHS and PROJECTION for the time integration
165 {
168 }
169 // Otherwise it gives an error (no implicit integration)
170 else
171 {
172 ASSERTL0(false, "Implicit unsteady Advection not set up.");
173 }
174}
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
virtual SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareField=true) override
Init object for UnsteadySystem class.
bool m_specHP_dealiasing
Flag to determine if dealisising is usde for the Spectral/hp element discretisation.
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
Get a SessionFunction by name.
static SOLVER_UTILS_EXPORT std::vector< ForcingSharedPtr > Load(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields=0)
Definition: Forcing.cpp:120
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
bool m_explicitAdvection
Indicates if explicit or implicit treatment of advection is used.
bool m_homoInitialFwd
Flag to determine if simulation should start in homogeneous forward transformed state.
Array< OneD, NekDouble > & GetNormalVelocity()
Get the normal velocity.
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the RHS.
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the projection.
void GetFluxVectorDeAlias(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
Evaluate the flux at each solution point using dealiasing.
void GetFluxVector(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
Evaluate the flux at each solution point.
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47
RiemannSolverFactory & GetRiemannSolverFactory()

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineOdeRhs(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineProjection(), DoOdeProjection(), DoOdeRhs(), Nektar::MultiRegions::eDiscontinuous, Nektar::MultiRegions::eGalerkin, Nektar::MultiRegions::eMixed_CG_Discontinuous, Nektar::SolverUtils::GetAdvectionFactory(), GetFluxVector(), GetFluxVectorDeAlias(), Nektar::SolverUtils::EquationSystem::GetFunction(), GetNormalVelocity(), Nektar::SolverUtils::GetRiemannSolverFactory(), Nektar::SolverUtils::EquationSystem::GetTraceNpoints(), Nektar::SolverUtils::Forcing::Load(), Nektar::SolverUtils::AdvectionSystem::m_advObject, Nektar::SolverUtils::UnsteadySystem::m_explicitAdvection, Nektar::SolverUtils::EquationSystem::m_fields, m_forcing, m_GJPJumpScale, Nektar::SolverUtils::UnsteadySystem::m_homoInitialFwd, Nektar::SolverUtils::UnsteadySystem::m_ode, Nektar::SolverUtils::EquationSystem::m_projectionType, m_riemannSolver, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_spacedim, Nektar::SolverUtils::EquationSystem::m_specHP_dealiasing, m_traceVn, m_useGJPStabilisation, m_velocity, m_waveFreq, and Nektar::SolverUtils::AdvectionSystem::v_InitObject().

Friends And Related Function Documentation

◆ MemoryManager< UnsteadyAdvection >

friend class MemoryManager< UnsteadyAdvection >
friend

Definition at line 1 of file UnsteadyAdvection.h.

Member Data Documentation

◆ className

string Nektar::UnsteadyAdvection::className
static
Initial value:
=
"UnsteadyAdvection", UnsteadyAdvection::create,
"Unsteady Advection equation.")
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
EquationSystemFactory & GetEquationSystemFactory()

Name of class.

Definition at line 61 of file UnsteadyAdvection.h.

◆ m_forcing

std::vector<SolverUtils::ForcingSharedPtr> Nektar::UnsteadyAdvection::m_forcing
protected

Forcing terms.

Definition at line 81 of file UnsteadyAdvection.h.

Referenced by DoOdeRhs(), and v_InitObject().

◆ m_GJPJumpScale

NekDouble Nektar::UnsteadyAdvection::m_GJPJumpScale
protected

Definition at line 69 of file UnsteadyAdvection.h.

Referenced by DoOdeProjection(), v_GenerateSummary(), and v_InitObject().

◆ m_planeNumber

int Nektar::UnsteadyAdvection::m_planeNumber
protected

Definition at line 78 of file UnsteadyAdvection.h.

Referenced by UnsteadyAdvection().

◆ m_riemannSolver

SolverUtils::RiemannSolverSharedPtr Nektar::UnsteadyAdvection::m_riemannSolver
protected

Definition at line 70 of file UnsteadyAdvection.h.

Referenced by v_InitObject().

◆ m_traceVn

Array<OneD, NekDouble> Nektar::UnsteadyAdvection::m_traceVn
protected

Definition at line 74 of file UnsteadyAdvection.h.

Referenced by GetNormalVelocity(), and v_InitObject().

◆ m_useGJPStabilisation

bool Nektar::UnsteadyAdvection::m_useGJPStabilisation
protected

Definition at line 67 of file UnsteadyAdvection.h.

Referenced by DoOdeProjection(), v_GenerateSummary(), and v_InitObject().

◆ m_velocity

Array<OneD, Array<OneD, NekDouble> > Nektar::UnsteadyAdvection::m_velocity
protected

Advection velocity.

Definition at line 73 of file UnsteadyAdvection.h.

Referenced by DoOdeRhs(), GetFluxVector(), GetFluxVectorDeAlias(), GetNormalVelocity(), and v_InitObject().

◆ m_waveFreq

NekDouble Nektar::UnsteadyAdvection::m_waveFreq
private

Definition at line 116 of file UnsteadyAdvection.h.

Referenced by v_InitObject().