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

#include <UnsteadyViscousBurgers.h>

Inheritance diagram for Nektar::UnsteadyViscousBurgers:
[legend]

Public Member Functions

virtual ~UnsteadyViscousBurgers ()
 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 ()
 
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)
 
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT void SetUpTraceNormals (void)
 
SOLVER_UTILS_EXPORT void InitObject (bool DeclareField=true)
 Initialises the members of this object. More...
 
SOLVER_UTILS_EXPORT void DoInitialise ()
 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 NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation. 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 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 int GetPararealIterationNumber ()
 
SOLVER_UTILS_EXPORT void SetPararealIterationNumber (int num)
 
SOLVER_UTILS_EXPORT bool GetUseInitialCondition ()
 
SOLVER_UTILS_EXPORT void SetUseInitialCondition (bool 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...
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp ()
 Virtual function to 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

 UnsteadyViscousBurgers (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Session reader. More...
 
void GetFluxVectorAdv (const Array< OneD, Array< OneD, NekDouble >> &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble >>> &flux)
 Evaluate the flux at each solution point for the advection part. More...
 
void GetFluxVectorDiff (const Array< OneD, Array< OneD, NekDouble >> &inarray, const Array< OneD, Array< OneD, Array< OneD, NekDouble >>> &qfield, Array< OneD, Array< OneD, Array< OneD, NekDouble >>> &viscousTensor)
 Evaluate the flux at each solution point for the diffusion part. 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)
 Perform the projection. More...
 
void DoImplicitSolve (const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, NekDouble time, NekDouble lambda)
 Solve implicitly the diffusion term. More...
 
Array< OneD, NekDouble > & GetNormalVelocity (Array< OneD, Array< OneD, NekDouble >> &inarray)
 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...
 
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve () override
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise () override
 Sets up initial conditions. 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_RequireFwdTrans ()
 
virtual SOLVER_UTILS_EXPORT void v_SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
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...
 
virtual SOLVER_UTILS_EXPORT bool v_UpdateTimeStepCheck ()
 
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 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_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 void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble >> &fieldcoeffs, std::vector< std::string > &variables)
 

Protected Attributes

bool m_useSpecVanVisc
 
bool m_useSpecVanViscVarDiff
 
NekDouble m_sVVCutoffRatio
 
NekDouble m_sVVDiffCoeff
 
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
 
SolverUtils::DiffusionSharedPtr m_diffusion
 
Array< OneD, NekDoublem_traceVn
 
StdRegions::VarCoeffMap m_varCoeffLap
 Variable Coefficient map for the Laplacian which can be activated as part of SVV or otherwise. More...
 
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
int m_abortSteps
 Number of steps between checks for abort conditions. More...
 
int m_filtersInfosteps
 Number of time steps between outputting filters information. More...
 
int m_nanSteps
 
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
NekDouble m_epsilon
 
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...
 
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state. More...
 
NekDouble m_steadyStateTol
 Tolerance to which steady state should be evaluated at. More...
 
int m_steadyStateSteps
 Check for steady state at step interval. More...
 
NekDouble m_steadyStateRes = 1.0
 
NekDouble m_steadyStateRes0 = 1.0
 
Array< OneD, Array< OneD, NekDouble > > m_previousSolution
 Storage for previous solution for steady-state check. More...
 
std::ofstream m_errFile
 
std::vector< int > m_intVariables
 
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
 
NekDouble m_filterTimeWarning
 Number of time steps between outputting status information. More...
 
NekDouble m_TimeIntegLambda = 0.0
 coefff of spacial derivatives(rhs or m_F in GLM) in calculating the residual of the whole equation(used in unsteady time integrations) More...
 
bool m_flagImplicitItsStatistics
 
bool m_flagImplicitSolver = false
 
Array< OneD, NekDoublem_magnitdEstimat
 estimate the magnitude of each conserved varibles More...
 
Array< OneD, NekDoublem_locTimeStep
 local time step(notice only for jfnk other see m_cflSafetyFactor) More...
 
NekDouble m_inArrayNorm = -1.0
 
int m_TotLinItePerStep = 0
 
int m_StagesPerStep = 1
 
bool m_flagUpdatePreconMat
 
int m_maxLinItePerNewton
 
int m_TotNewtonIts = 0
 
int m_TotLinIts = 0
 
int m_TotImpStages = 0
 
bool m_CalcPhysicalAV = true
 flag to update artificial viscosity 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_timestepMax = -1.0
 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_pararealIter
 Number of parareal 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_useInitialCondition
 Flag to determine if IC are used. 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
 
NekDouble m_epsilon
 

Friends

class MemoryManager< UnsteadyViscousBurgers >
 

Additional Inherited Members

- Public Attributes inherited from Nektar::SolverUtils::UnsteadySystem
NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1). More...
 
NekDouble m_cflNonAcoustic
 
NekDouble m_CFLGrowth
 CFL growth rate. More...
 
NekDouble m_CFLEnd
 maximun cfl in cfl growth More...
 
- 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 []
 

Detailed Description

Definition at line 48 of file UnsteadyViscousBurgers.h.

Constructor & Destructor Documentation

◆ ~UnsteadyViscousBurgers()

Nektar::UnsteadyViscousBurgers::~UnsteadyViscousBurgers ( )
virtual

Destructor.

Unsteady linear advection diffusion equation destructor.

Definition at line 175 of file UnsteadyViscousBurgers.cpp.

176 {
177 }

◆ UnsteadyViscousBurgers()

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

Session reader.

Definition at line 50 of file UnsteadyViscousBurgers.cpp.

53  : UnsteadySystem(pSession, pGraph), AdvectionSystem(pSession, pGraph),
55 {
56  m_planeNumber = 0;
57 }
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.
StdRegions::VarCoeffMap m_varCoeffLap
Variable Coefficient map for the Laplacian which can be activated as part of SVV or otherwise.
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:344

References m_planeNumber.

Member Function Documentation

◆ create()

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

Creates an instance of this class.

Definition at line 54 of file UnsteadyViscousBurgers.h.

57  {
60  pGraph);
61  p->InitObject();
62  return p;
63  }
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.

◆ DoImplicitSolve()

void Nektar::UnsteadyViscousBurgers::DoImplicitSolve ( const Array< OneD, const Array< OneD, NekDouble >> &  inarray,
Array< OneD, Array< OneD, NekDouble >> &  outarray,
NekDouble  time,
NekDouble  lambda 
)
protected

Solve implicitly the diffusion term.

Definition at line 310 of file UnsteadyViscousBurgers.cpp.

314 {
315  int nvariables = inarray.size();
316  int nq = m_fields[0]->GetNpoints();
317 
319  factors[StdRegions::eFactorLambda] = 1.0 / lambda / m_epsilon;
320 
321  if (m_useSpecVanVisc)
322  {
325  }
326 
327  Array<OneD, Array<OneD, NekDouble>> F(nvariables);
328  F[0] = Array<OneD, NekDouble>(nq * nvariables);
329 
330  for (int n = 1; n < nvariables; ++n)
331  {
332  F[n] = F[n - 1] + nq;
333  }
334 
335  // We solve ( \nabla^2 - HHlambda ) Y[i] = rhs [i]
336  // inarray = input: \hat{rhs} -> output: \hat{Y}
337  // outarray = output: nabla^2 \hat{Y}
338  // where \hat = modal coeffs
339  for (int i = 0; i < nvariables; ++i)
340  {
341  // Multiply 1.0/timestep/lambda
342  Vmath::Smul(nq, -factors[StdRegions::eFactorLambda], inarray[i], 1,
343  F[i], 1);
344  }
345 
346  // Setting boundary conditions
347  SetBoundaryConditions(time);
348 
350  {
351  static int cnt = 0;
352 
353  if (cnt % 10 == 0)
354  {
355  Array<OneD, Array<OneD, NekDouble>> vel(m_fields.size());
356  for (int i = 0; i < m_fields.size(); ++i)
357  {
358  m_fields[i]->ClearGlobalLinSysManager();
359  vel[i] = m_fields[i]->UpdatePhys();
360  }
362  }
363  ++cnt;
364  }
365  for (int i = 0; i < nvariables; ++i)
366  {
367  // Solve a system of equations with Helmholtz solver
368  m_fields[i]->HelmSolve(F[i], m_fields[i]->UpdateCoeffs(), factors,
369  m_varCoeffLap);
370 
371  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(), outarray[i]);
372  }
373 }
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
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 t...
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:399
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:248

References Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::eFactorSVVCutoffRatio, Nektar::StdRegions::eFactorSVVDiffCoeff, m_epsilon, Nektar::SolverUtils::EquationSystem::m_fields, m_sVVCutoffRatio, m_sVVDiffCoeff, m_useSpecVanVisc, m_useSpecVanViscVarDiff, m_varCoeffLap, Nektar::SolverUtils::EquationSystem::SetBoundaryConditions(), Vmath::Smul(), and Nektar::SolverUtils::UnsteadySystem::SVVVarDiffCoeff().

Referenced by v_InitObject().

◆ DoOdeProjection()

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

Perform the projection.

Compute the projection for the unsteady advection diffusion problem.

Parameters
inarrayGiven fields.
outarrayCalculated solution.
timeTime.

Definition at line 270 of file UnsteadyViscousBurgers.cpp.

273 {
274  int i;
275  int nvariables = inarray.size();
276  SetBoundaryConditions(time);
277  switch (m_projectionType)
278  {
280  {
281  // Just copy over array
282  int npoints = GetNpoints();
283 
284  for (i = 0; i < nvariables; ++i)
285  {
286  Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
287  }
288  break;
289  }
292  {
293  // Do nothing for the moment.
294  }
295  default:
296  {
297  ASSERTL0(false, "Unknown projection scheme");
298  break;
299  }
300  }
301 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
SOLVER_UTILS_EXPORT int GetNpoints()
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255

References ASSERTL0, Nektar::MultiRegions::eDiscontinuous, Nektar::MultiRegions::eGalerkin, Nektar::MultiRegions::eMixed_CG_Discontinuous, Nektar::SolverUtils::EquationSystem::GetNpoints(), Nektar::SolverUtils::EquationSystem::m_projectionType, Nektar::SolverUtils::EquationSystem::SetBoundaryConditions(), and Vmath::Vcopy().

Referenced by v_InitObject().

◆ DoOdeRhs()

void Nektar::UnsteadyViscousBurgers::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 unsteady linear advection diffusion problem.

Parameters
inarrayGiven fields.
outarrayCalculated solution.
timeTime.

Definition at line 216 of file UnsteadyViscousBurgers.cpp.

219 {
220  // Number of fields (variables of the problem)
221  int nVariables = inarray.size();
222 
223  // Number of solution points
224  int nSolutionPts = GetNpoints();
225 
226  Array<OneD, Array<OneD, NekDouble>> outarrayDiff(nVariables);
227 
228  for (int i = 0; i < nVariables; ++i)
229  {
230  outarrayDiff[i] = Array<OneD, NekDouble>(nSolutionPts, 0.0);
231  }
232 
233  // RHS computation using the new advection base class
234  m_advObject->Advect(nVariables, m_fields, inarray, inarray, outarray, time);
235 
236  // Negate the RHS
237  for (int i = 0; i < nVariables; ++i)
238  {
239  Vmath::Neg(nSolutionPts, outarray[i], 1);
240  }
241 
242  // No explicit diffusion for CG
244  {
245  m_diffusion->Diffuse(nVariables, m_fields, inarray, outarrayDiff);
246 
247  for (int i = 0; i < nVariables; ++i)
248  {
249  Vmath::Vadd(nSolutionPts, &outarray[i][0], 1, &outarrayDiff[i][0],
250  1, &outarray[i][0], 1);
251  }
252  }
253 
254  // Add forcing terms
255  for (auto &x : m_forcing)
256  {
257  // set up non-linear terms
258  x->Apply(m_fields, inarray, outarray, time);
259  }
260 }
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
SolverUtils::DiffusionSharedPtr m_diffusion
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Forcing terms.
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:518
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:359

References Nektar::MultiRegions::eDiscontinuous, Nektar::SolverUtils::EquationSystem::GetNpoints(), Nektar::SolverUtils::AdvectionSystem::m_advObject, m_diffusion, Nektar::SolverUtils::EquationSystem::m_fields, m_forcing, Nektar::SolverUtils::EquationSystem::m_projectionType, Vmath::Neg(), and Vmath::Vadd().

Referenced by v_InitObject().

◆ GetFluxVectorAdv()

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

Evaluate the flux at each solution point for the advection part.

Return the flux vector for the advection part.

Parameters
physfieldFields.
fluxResulting flux.

Definition at line 381 of file UnsteadyViscousBurgers.cpp.

384 {
385 
386  const int nq = m_fields[0]->GetNpoints();
387 
388  for (int i = 0; i < flux.size(); ++i)
389  {
390  for (int j = 0; j < flux[0].size(); ++j)
391  {
392  Vmath::Vmul(nq, physfield[i], 1, physfield[j], 1, flux[i][j], 1);
393  }
394  }
395 }
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:209

References Nektar::SolverUtils::EquationSystem::m_fields, and Vmath::Vmul().

Referenced by v_InitObject().

◆ GetFluxVectorDiff()

void Nektar::UnsteadyViscousBurgers::GetFluxVectorDiff ( const Array< OneD, Array< OneD, NekDouble >> &  inarray,
const Array< OneD, Array< OneD, Array< OneD, NekDouble >>> &  qfield,
Array< OneD, Array< OneD, Array< OneD, NekDouble >>> &  viscousTensor 
)
protected

Evaluate the flux at each solution point for the diffusion part.

Return the flux vector for the diffusion part.

Parameters
iEquation number.
jSpatial direction.
physfieldFields.
derivativesFirst order derivatives.
fluxResulting flux.

Definition at line 406 of file UnsteadyViscousBurgers.cpp.

410 {
411  boost::ignore_unused(inarray);
412 
413  unsigned int nDim = qfield.size();
414  unsigned int nConvectiveFields = qfield[0].size();
415  unsigned int nPts = qfield[0][0].size();
416 
417  for (unsigned int j = 0; j < nDim; ++j)
418  {
419  for (unsigned int i = 0; i < nConvectiveFields; ++i)
420  {
421  Vmath::Smul(nPts, m_epsilon, qfield[j][i], 1, viscousTensor[j][i],
422  1);
423  }
424  }
425 }

References m_epsilon, and Vmath::Smul().

Referenced by v_InitObject().

◆ GetNormalVelocity()

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

Get the normal velocity.

Get the normal velocity for the unsteady linear advection diffusion equation.

Definition at line 183 of file UnsteadyViscousBurgers.cpp.

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

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

◆ v_GenerateSummary()

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

Print Summary.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 427 of file UnsteadyViscousBurgers.cpp.

428 {
430  if (m_useSpecVanVisc)
431  {
432  stringstream ss;
433  ss << "SVV (cut off = " << m_sVVCutoffRatio
434  << ", coeff = " << m_sVVDiffCoeff << ")";
435  AddSummaryItem(s, "Smoothing", ss.str());
436  }
437 }
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_sVVCutoffRatio, m_sVVDiffCoeff, m_useSpecVanVisc, and Nektar::SolverUtils::UnsteadySystem::v_GenerateSummary().

◆ v_InitObject()

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

Initialise the object.

Initialisation object for the unsteady linear advection diffusion equation.

Reimplemented from Nektar::SolverUtils::AdvectionSystem.

Definition at line 63 of file UnsteadyViscousBurgers.cpp.

64 {
65  AdvectionSystem::v_InitObject(DeclareFields);
66 
67  m_session->LoadParameter("wavefreq", m_waveFreq, 0.0);
68  m_session->LoadParameter("epsilon", m_epsilon, 0.0);
69 
70  m_session->MatchSolverInfo("SpectralVanishingViscosity", "True",
71  m_useSpecVanVisc, false);
72  m_session->MatchSolverInfo("SpectralVanishingViscosity", "VarDiff",
75  {
76  m_useSpecVanVisc = true;
77  }
78 
79  if (m_useSpecVanVisc)
80  {
81  m_session->LoadParameter("SVVCutoffRatio", m_sVVCutoffRatio, 0.75);
82  m_session->LoadParameter("SVVDiffCoeff", m_sVVDiffCoeff, 0.1);
83  }
84 
85  // Type of advection and diffusion classes to be used
86  switch (m_projectionType)
87  {
88  // Discontinuous field
90  {
91  ASSERTL0(false, "Need to implement for DG");
92  // Do not forwards transform initial condition
93  m_homoInitialFwd = false;
94 
95  // Advection term
96  string advName;
97  string riemName;
98  m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
100  advName, advName);
101  m_advObject->SetFluxVector(
103  m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
106  riemName, m_session);
107  m_advObject->SetRiemannSolver(m_riemannSolver);
108  m_advObject->InitObject(m_session, m_fields);
109 
110  // Diffusion term
111  std::string diffName;
112  m_session->LoadSolverInfo("DiffusionType", diffName, "LDG");
114  diffName, diffName);
115  m_diffusion->SetFluxVector(
117  m_diffusion->InitObject(m_session, m_fields);
118  break;
119  }
120  // Continuous field
123  {
124  // Advection term
125  std::string advName;
126  m_session->LoadSolverInfo("AdvectionType", advName,
127  "NonConservative");
129  advName, advName);
130  m_advObject->SetFluxVector(
132 
134  {
135  Array<OneD, Array<OneD, NekDouble>> vel(m_fields.size());
136  for (int i = 0; i < m_fields.size(); ++i)
137  {
138  vel[i] = m_fields[i]->UpdatePhys();
139  }
141  }
142 
143  // In case of Galerkin explicit diffusion gives an error
145  {
146  ASSERTL0(false, "Explicit Galerkin diffusion not set up.");
147  }
148  // In case of Galerkin implicit diffusion: do nothing
149  break;
150  }
151  default:
152  {
153  ASSERTL0(false, "Unsupported projection type.");
154  break;
155  }
156  }
157 
158  // Forcing terms
159  m_forcing = SolverUtils::Forcing::Load(m_session, shared_from_this(),
160  m_fields, m_fields.size());
161 
164 
166  m_explicitDiffusion == 1)
167  {
169  }
170 }
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)
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
virtual SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareField=true) override
Init object for UnsteadySystem class.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
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_explicitDiffusion
Indicates if explicit or implicit treatment of diffusion is used.
bool m_homoInitialFwd
Flag to determine if simulation should start in homogeneous forward transformed state.
void DoImplicitSolve(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, NekDouble time, NekDouble lambda)
Solve implicitly the diffusion term.
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)
Perform the projection.
void GetFluxVectorAdv(const Array< OneD, Array< OneD, NekDouble >> &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble >>> &flux)
Evaluate the flux at each solution point for the advection part.
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
void GetFluxVectorDiff(const Array< OneD, Array< OneD, NekDouble >> &inarray, const Array< OneD, Array< OneD, Array< OneD, NekDouble >>> &qfield, Array< OneD, Array< OneD, Array< OneD, NekDouble >>> &viscousTensor)
Evaluate the flux at each solution point for the diffusion part.
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:41
RiemannSolverFactory & GetRiemannSolverFactory()

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineImplicitSolve(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineOdeRhs(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineProjection(), DoImplicitSolve(), DoOdeProjection(), DoOdeRhs(), Nektar::MultiRegions::eDiscontinuous, Nektar::MultiRegions::eGalerkin, Nektar::MultiRegions::eMixed_CG_Discontinuous, Nektar::SolverUtils::GetAdvectionFactory(), Nektar::SolverUtils::GetDiffusionFactory(), GetFluxVectorAdv(), GetFluxVectorDiff(), Nektar::SolverUtils::GetRiemannSolverFactory(), Nektar::SolverUtils::Forcing::Load(), Nektar::SolverUtils::AdvectionSystem::m_advObject, m_diffusion, m_epsilon, Nektar::SolverUtils::UnsteadySystem::m_explicitDiffusion, Nektar::SolverUtils::EquationSystem::m_fields, m_forcing, Nektar::SolverUtils::UnsteadySystem::m_homoInitialFwd, Nektar::SolverUtils::UnsteadySystem::m_ode, Nektar::SolverUtils::EquationSystem::m_projectionType, m_riemannSolver, Nektar::SolverUtils::EquationSystem::m_session, m_sVVCutoffRatio, m_sVVDiffCoeff, m_useSpecVanVisc, m_useSpecVanViscVarDiff, m_varCoeffLap, m_waveFreq, Nektar::SolverUtils::UnsteadySystem::SVVVarDiffCoeff(), and Nektar::SolverUtils::AdvectionSystem::v_InitObject().

Friends And Related Function Documentation

◆ MemoryManager< UnsteadyViscousBurgers >

friend class MemoryManager< UnsteadyViscousBurgers >
friend

Definition at line 1 of file UnsteadyViscousBurgers.h.

Member Data Documentation

◆ className

string Nektar::UnsteadyViscousBurgers::className
static
Initial value:
=
"UnsteadyViscousBurgers", UnsteadyViscousBurgers::create)
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 65 of file UnsteadyViscousBurgers.h.

◆ m_diffusion

SolverUtils::DiffusionSharedPtr Nektar::UnsteadyViscousBurgers::m_diffusion
protected

Definition at line 78 of file UnsteadyViscousBurgers.h.

Referenced by DoOdeRhs(), and v_InitObject().

◆ m_epsilon

NekDouble Nektar::UnsteadyViscousBurgers::m_epsilon
private

Definition at line 135 of file UnsteadyViscousBurgers.h.

Referenced by DoImplicitSolve(), GetFluxVectorDiff(), and v_InitObject().

◆ m_forcing

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

Forcing terms.

Definition at line 90 of file UnsteadyViscousBurgers.h.

Referenced by DoOdeRhs(), and v_InitObject().

◆ m_planeNumber

int Nektar::UnsteadyViscousBurgers::m_planeNumber
protected

Definition at line 87 of file UnsteadyViscousBurgers.h.

Referenced by UnsteadyViscousBurgers().

◆ m_riemannSolver

SolverUtils::RiemannSolverSharedPtr Nektar::UnsteadyViscousBurgers::m_riemannSolver
protected

Definition at line 77 of file UnsteadyViscousBurgers.h.

Referenced by v_InitObject().

◆ m_sVVCutoffRatio

NekDouble Nektar::UnsteadyViscousBurgers::m_sVVCutoffRatio
protected

Definition at line 75 of file UnsteadyViscousBurgers.h.

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

◆ m_sVVDiffCoeff

NekDouble Nektar::UnsteadyViscousBurgers::m_sVVDiffCoeff
protected

Definition at line 76 of file UnsteadyViscousBurgers.h.

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

◆ m_traceVn

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

Definition at line 79 of file UnsteadyViscousBurgers.h.

Referenced by GetNormalVelocity().

◆ m_useSpecVanVisc

bool Nektar::UnsteadyViscousBurgers::m_useSpecVanVisc
protected

Definition at line 71 of file UnsteadyViscousBurgers.h.

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

◆ m_useSpecVanViscVarDiff

bool Nektar::UnsteadyViscousBurgers::m_useSpecVanViscVarDiff
protected

Definition at line 72 of file UnsteadyViscousBurgers.h.

Referenced by DoImplicitSolve(), and v_InitObject().

◆ m_varCoeffLap

StdRegions::VarCoeffMap Nektar::UnsteadyViscousBurgers::m_varCoeffLap
protected

Variable Coefficient map for the Laplacian which can be activated as part of SVV or otherwise.

Definition at line 83 of file UnsteadyViscousBurgers.h.

Referenced by DoImplicitSolve(), and v_InitObject().

◆ m_waveFreq

NekDouble Nektar::UnsteadyViscousBurgers::m_waveFreq
private

Definition at line 134 of file UnsteadyViscousBurgers.h.

Referenced by v_InitObject().