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 (void)
 
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...
 
- 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 ()
 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, Array< OneD, NekDouble > &output)
 
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 void SetTime (const NekDouble time)
 
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...
 

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

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...
 
virtual 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...
 
virtual 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 ()
 Initialise the object. More...
 
virtual void v_GenerateSummary (SolverUtils::SummaryList &s)
 Print Summary. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::AdvectionSystem
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT Array< OneD, NekDoublev_GetMaxStdVelocity ()
 
- 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 ()
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise ()
 Sets up initial conditions. More...
 
virtual SOLVER_UTILS_EXPORT void v_AppendOutput1D (Array< OneD, Array< OneD, NekDouble >> &solution1D)
 Print the solution at each solution point in a txt file. 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 ()
 
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...
 
- 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_infosteps
 Number of time steps between outputting status information. More...
 
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::TimeIntegrationWrapperSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
LibUtilities::TimeIntegrationSolutionSharedPtr m_intSoln
 
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...
 
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...
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
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...
 
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_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...
 
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...
 
- 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 176 of file UnsteadyViscousBurgers.cpp.

177  {
178  }

◆ UnsteadyViscousBurgers()

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

Session reader.

Definition at line 51 of file UnsteadyViscousBurgers.cpp.

References m_planeNumber.

54  : UnsteadySystem(pSession, pGraph),
55  AdvectionSystem(pSession, pGraph),
57  {
58  m_planeNumber = 0;
59  }
StdRegions::VarCoeffMap m_varCoeffLap
Variable Coefficient map for the Laplacian which can be activated as part of SVV or otherwise...
SOLVER_UTILS_EXPORT UnsteadySystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises UnsteadySystem class members.
SOLVER_UTILS_EXPORT AdvectionSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:265

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.

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

57  {
60  AllocateSharedPtr(pSession, 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.

◆ DoImplicitSolve()

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

Solve implicitly the diffusion term.

Definition at line 318 of file UnsteadyViscousBurgers.cpp.

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::NullFlagList, Nektar::SolverUtils::EquationSystem::SetBoundaryConditions(), Vmath::Smul(), and Nektar::SolverUtils::UnsteadySystem::SVVVarDiffCoeff().

Referenced by v_InitObject().

323  {
324  int nvariables = inarray.num_elements();
325  int nq = m_fields[0]->GetNpoints();
326 
328  factors[StdRegions::eFactorLambda] = 1.0/lambda/m_epsilon;
329 
330  if(m_useSpecVanVisc)
331  {
334  }
335 
336  Array<OneD, Array< OneD, NekDouble> > F(nvariables);
337  F[0] = Array<OneD, NekDouble> (nq*nvariables);
338 
339  for (int n = 1; n < nvariables; ++n)
340  {
341  F[n] = F[n-1] + nq;
342  }
343 
344  // We solve ( \nabla^2 - HHlambda ) Y[i] = rhs [i]
345  // inarray = input: \hat{rhs} -> output: \hat{Y}
346  // outarray = output: nabla^2 \hat{Y}
347  // where \hat = modal coeffs
348  for (int i = 0; i < nvariables; ++i)
349  {
350  // Multiply 1.0/timestep/lambda
352  inarray[i], 1, F[i], 1);
353  }
354 
355  //Setting boundary conditions
356  SetBoundaryConditions(time);
357 
359  {
360  static int cnt = 0;
361 
362  if(cnt %10 == 0)
363  {
364  Array<OneD, Array<OneD, NekDouble> > vel(m_fields.num_elements());
365  for(int i = 0; i < m_fields.num_elements(); ++i)
366  {
367  m_fields[i]->ClearGlobalLinSysManager();
368  vel[i] = m_fields[i]->UpdatePhys();
369  }
371  }
372  ++cnt;
373  }
374  for (int i = 0; i < nvariables; ++i)
375  {
376  // Solve a system of equations with Helmholtz solver
377  m_fields[i]->HelmSolve(F[i], m_fields[i]->UpdateCoeffs(),
378  NullFlagList, factors, m_varCoeffLap);
379 
380  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(), outarray[i]);
381  }
382  }
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:294
StdRegions::VarCoeffMap m_varCoeffLap
Variable Coefficient map for the Laplacian which can be activated as part of SVV or otherwise...
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:216
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT void SVVVarDiffCoeff(const Array< OneD, Array< OneD, NekDouble >> vel, StdRegions::VarCoeffMap &varCoeffMap)
Evaluate the SVV diffusion coefficient according to Moura&#39;s paper where it should proportional to h t...
static FlagList NullFlagList
An empty flag list.

◆ 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 276 of file UnsteadyViscousBurgers.cpp.

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

280  {
281  int i;
282  int nvariables = inarray.num_elements();
283  SetBoundaryConditions(time);
284  switch(m_projectionType)
285  {
287  {
288  // Just copy over array
289  int npoints = GetNpoints();
290 
291  for(i = 0; i < nvariables; ++i)
292  {
293  Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
294  }
295  break;
296  }
299  {
300  // Do nothing for the moment.
301  }
302  default:
303  {
304  ASSERTL0(false, "Unknown projection scheme");
305  break;
306  }
307  }
308  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
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.
SOLVER_UTILS_EXPORT int GetNpoints()
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064

◆ DoOdeRhs()

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

Compute the RHS.

Compute the right-hand side for the unsteady linear advection diffusion problem.

Parameters
inarrayGiven fields.
outarrayCalculated solution.
timeTime.

Definition at line 220 of file UnsteadyViscousBurgers.cpp.

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

224  {
225  // Number of fields (variables of the problem)
226  int nVariables = inarray.num_elements();
227 
228  // Number of solution points
229  int nSolutionPts = GetNpoints();
230 
231  Array<OneD, Array<OneD, NekDouble> > outarrayDiff(nVariables);
232 
233  for (int i = 0; i < nVariables; ++i)
234  {
235  outarrayDiff[i] = Array<OneD, NekDouble>(nSolutionPts, 0.0);
236  }
237 
238  // RHS computation using the new advection base class
239  m_advObject->Advect(nVariables, m_fields, inarray,
240  inarray, outarray, time);
241 
242  // Negate the RHS
243  for (int i = 0; i < nVariables; ++i)
244  {
245  Vmath::Neg(nSolutionPts, outarray[i], 1);
246  }
247 
248  // No explicit diffusion for CG
250  {
251  m_diffusion->Diffuse(nVariables, m_fields, inarray, outarrayDiff);
252 
253  for (int i = 0; i < nVariables; ++i)
254  {
255  Vmath::Vadd(nSolutionPts, &outarray[i][0], 1,
256  &outarrayDiff[i][0], 1, &outarray[i][0], 1);
257  }
258  }
259 
260  // Add forcing terms
261  for (auto &x : m_forcing)
262  {
263  // set up non-linear terms
264  x->Apply(m_fields, inarray, outarray, time);
265  }
266  }
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:399
SolverUtils::DiffusionSharedPtr m_diffusion
SOLVER_UTILS_EXPORT int GetNpoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
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:302
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Forcing terms.

◆ 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 390 of file UnsteadyViscousBurgers.cpp.

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

Referenced by v_InitObject().

393  {
394 
395  const int nq = m_fields[0]->GetNpoints();
396 
397  for (int i = 0; i < flux.num_elements(); ++i)
398  {
399  for (int j = 0; j < flux[0].num_elements(); ++j)
400  {
401  Vmath::Vmul(nq, physfield[i], 1, physfield[j], 1,
402  flux[i][j], 1);
403  }
404  }
405  }
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
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:186

◆ 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 416 of file UnsteadyViscousBurgers.cpp.

References m_epsilon, and Vmath::Smul().

Referenced by v_InitObject().

420  {
421  boost::ignore_unused(inarray);
422 
423  unsigned int nDim = qfield.num_elements();
424  unsigned int nConvectiveFields = qfield[0].num_elements();
425  unsigned int nPts = qfield[0][0].num_elements();
426 
427  for (unsigned int j = 0; j < nDim; ++j)
428  {
429  for (unsigned int i = 0; i < nConvectiveFields; ++i)
430  {
431  Vmath::Smul(nPts, m_epsilon, qfield[j][i], 1,
432  viscousTensor[j][i], 1 );
433  }
434  }
435  }
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:216

◆ 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 184 of file UnsteadyViscousBurgers.cpp.

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

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

◆ v_GenerateSummary()

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

Print Summary.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 437 of file UnsteadyViscousBurgers.cpp.

References Nektar::SolverUtils::AddSummaryItem(), m_sVVCutoffRatio, m_sVVDiffCoeff, m_useSpecVanVisc, and Nektar::SolverUtils::UnsteadySystem::v_GenerateSummary().

439  {
441  if(m_useSpecVanVisc)
442  {
443  stringstream ss;
444  ss << "SVV (cut off = " << m_sVVCutoffRatio
445  << ", coeff = " << m_sVVDiffCoeff << ")";
446  AddSummaryItem(s, "Smoothing", ss.str());
447  }
448  }
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s)
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

◆ v_InitObject()

void Nektar::UnsteadyViscousBurgers::v_InitObject ( )
protectedvirtual

Initialise the object.

Initialisation object for the unsteady linear advection diffusion equation.

Reimplemented from Nektar::SolverUtils::AdvectionSystem.

Definition at line 65 of file UnsteadyViscousBurgers.cpp.

References ASSERTL0, 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().

66  {
68 
69  m_session->LoadParameter("wavefreq", m_waveFreq, 0.0);
70  m_session->LoadParameter("epsilon", m_epsilon, 0.0);
71 
72  m_session->MatchSolverInfo(
73  "SpectralVanishingViscosity", "True", m_useSpecVanVisc, false);
74  m_session->MatchSolverInfo(
75  "SpectralVanishingViscosity", "VarDiff", m_useSpecVanViscVarDiff, false);
77  {
78  m_useSpecVanVisc = true;
79  }
80 
82  {
83  m_session->LoadParameter("SVVCutoffRatio",m_sVVCutoffRatio,0.75);
84  m_session->LoadParameter("SVVDiffCoeff",m_sVVDiffCoeff,0.1);
85  }
86 
87  // Type of advection and diffusion classes to be used
88  switch(m_projectionType)
89  {
90  // Discontinuous field
92  {
93  ASSERTL0(false,"Need to implement for DG");
94  // Do not forwards transform initial condition
95  m_homoInitialFwd = false;
96 
97  // Advection term
98  string advName;
99  string riemName;
100  m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
102  CreateInstance(advName, advName);
103  m_advObject->SetFluxVector(&UnsteadyViscousBurgers::
104  GetFluxVectorAdv, this);
105  m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
107  CreateInstance(riemName, m_session);
108  m_advObject->SetRiemannSolver(m_riemannSolver);
109  m_advObject->InitObject (m_session, m_fields);
110 
111  // Diffusion term
112  std::string diffName;
113  m_session->LoadSolverInfo("DiffusionType", diffName, "LDG");
115  CreateInstance(diffName, diffName);
116  m_diffusion->SetFluxVector(&UnsteadyViscousBurgers::
117  GetFluxVectorDiff, this);
118  m_diffusion->InitObject(m_session, m_fields);
119  break;
120  }
121  // Continuous field
124  {
125  // Advection term
126  std::string advName;
127  m_session->LoadSolverInfo("AdvectionType", advName,
128  "NonConservative");
130  CreateInstance(advName, advName);
131  m_advObject->SetFluxVector(&UnsteadyViscousBurgers::
132  GetFluxVectorAdv, this);
133 
135  {
136  Array<OneD, Array<OneD, NekDouble> > vel(m_fields.num_elements());
137  for(int i = 0; i < m_fields.num_elements(); ++i)
138  {
139  vel[i] = m_fields[i]->UpdatePhys();
140  }
142  }
143 
144  // In case of Galerkin explicit diffusion gives an error
146  {
147  ASSERTL0(false, "Explicit Galerkin diffusion not set up.");
148  }
149  // In case of Galerkin implicit diffusion: do nothing
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.num_elements());
162 
165 
167  m_explicitDiffusion == 1)
168  {
170  }
171  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
bool m_homoInitialFwd
Flag to determine if simulation should start in homogeneous forward transformed state.
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
bool m_explicitDiffusion
Indicates if explicit or implicit treatment of diffusion is used.
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:41
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
StdRegions::VarCoeffMap m_varCoeffLap
Variable Coefficient map for the Laplacian which can be activated as part of SVV or otherwise...
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
UnsteadyViscousBurgers(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
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:85
RiemannSolverFactory & GetRiemannSolverFactory()
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.
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Perform the projection.
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47
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.
virtual void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the RHS.
SolverUtils::DiffusionSharedPtr m_diffusion
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
virtual 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.
SOLVER_UTILS_EXPORT void SVVVarDiffCoeff(const Array< OneD, Array< OneD, NekDouble >> vel, StdRegions::VarCoeffMap &varCoeffMap)
Evaluate the SVV diffusion coefficient according to Moura&#39;s paper where it should proportional to h t...
virtual SOLVER_UTILS_EXPORT void v_InitObject()
Init object for UnsteadySystem class.
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Forcing terms.

Friends And Related Function Documentation

◆ MemoryManager< UnsteadyViscousBurgers >

friend class MemoryManager< UnsteadyViscousBurgers >
friend

Definition at line 51 of file UnsteadyViscousBurgers.h.

Member Data Documentation

◆ className

string Nektar::UnsteadyViscousBurgers::className
static
Initial value:

Name of class.

Definition at line 65 of file UnsteadyViscousBurgers.h.

◆ m_diffusion

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

Definition at line 76 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 86 of file UnsteadyViscousBurgers.h.

Referenced by DoOdeRhs(), and v_InitObject().

◆ m_planeNumber

int Nektar::UnsteadyViscousBurgers::m_planeNumber
protected

Definition at line 83 of file UnsteadyViscousBurgers.h.

Referenced by UnsteadyViscousBurgers().

◆ m_riemannSolver

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

Definition at line 75 of file UnsteadyViscousBurgers.h.

Referenced by v_InitObject().

◆ m_sVVCutoffRatio

NekDouble Nektar::UnsteadyViscousBurgers::m_sVVCutoffRatio
protected

Definition at line 73 of file UnsteadyViscousBurgers.h.

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

◆ m_sVVDiffCoeff

NekDouble Nektar::UnsteadyViscousBurgers::m_sVVDiffCoeff
protected

Definition at line 74 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 77 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 80 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().