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

#include <UnsteadyViscousBurgers.h>

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

Public Member Functions

virtual ~UnsteadyViscousBurgers ()
 Destructor. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::AdvectionSystem
SOLVER_UTILS_EXPORT AdvectionSystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 
virtual SOLVER_UTILS_EXPORT ~AdvectionSystem ()
 
AdvectionSharedPtr GetAdvObject ()
 Returns the advection object held by this instance. More...
 
- 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 >
boost::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 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 void EvaluateFunction (Array< OneD, Array< OneD, NekDouble > > &pArray, std::string pFunctionName, const NekDouble pTime=0.0, const int domain=0)
 Evaluates a function as specified in the session file. More...
 
SOLVER_UTILS_EXPORT void EvaluateFunction (std::vector< std::string > pFieldNames, Array< OneD, Array< OneD, NekDouble > > &pFields, const std::string &pName, const NekDouble &pTime=0.0, const int domain=0)
 Populate given fields with the function from session. More...
 
SOLVER_UTILS_EXPORT void EvaluateFunction (std::vector< std::string > pFieldNames, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const std::string &pName, const NekDouble &pTime=0.0, const int domain=0)
 Populate given fields with the function from session. More...
 
SOLVER_UTILS_EXPORT void EvaluateFunction (std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, const NekDouble &pTime=0.0, const int domain=0)
 
SOLVER_UTILS_EXPORT std::string DescribeFunction (std::string pFieldName, const std::string &pFunctionName, const int domain)
 Provide a description of a function for a given field name. More...
 
SOLVER_UTILS_EXPORT void InitialiseBaseFlow (Array< OneD, Array< OneD, NekDouble > > &base)
 Perform initialisation of the base flow. 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, NekDouble
ErrorExtraPoints (unsigned int field)
 Compute error (L2 and L_inf) over an larger set of quadrature points return [L2 Linf]. More...
 
SOLVER_UTILS_EXPORT void WeakAdvectionGreensDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\nabla \phi \cdot F) $. More...
 
SOLVER_UTILS_EXPORT void WeakAdvectionDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\phi, \nabla \cdot F) $. More...
 
SOLVER_UTILS_EXPORT void WeakAdvectionNonConservativeForm (const Array< OneD, Array< OneD, NekDouble > > &V, const Array< OneD, const NekDouble > &u, Array< OneD, NekDouble > &outarray, bool UseContCoeffs=false)
 Compute the inner product $ (\phi, V\cdot \nabla u) $. More...
 
f SOLVER_UTILS_EXPORT void AdvectionNonConservativeForm (const Array< OneD, Array< OneD, NekDouble > > &V, const Array< OneD, const NekDouble > &u, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wk=NullNekDouble1DArray)
 Compute the non-conservative advection. More...
 
SOLVER_UTILS_EXPORT void WeakDGAdvection (const Array< OneD, Array< OneD, NekDouble > > &InField, Array< OneD, Array< OneD, NekDouble > > &OutField, bool NumericalFluxIncludesNormal=true, bool InFieldIsInPhysSpace=false, int nvariables=0)
 Calculate the weak discontinuous Galerkin advection. More...
 
SOLVER_UTILS_EXPORT void WeakDGDiffusion (const Array< OneD, Array< OneD, NekDouble > > &InField, Array< OneD, Array< OneD, NekDouble > > &OutField, bool NumericalFluxIncludesNormal=true, bool InFieldIsInPhysSpace=false)
 Calculate weak DG Diffusion in the LDG form. 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 ScanForHistoryPoints ()
 Builds map of which element holds each history point. More...
 
SOLVER_UTILS_EXPORT void WriteHistoryData (std::ostream &out)
 Probe each history point and write to 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::FieldMetaDataMap
UpdateFieldMetaDataMap ()
 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 GetNumElmVelocity ()
 
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 GetFluxVector (const int i, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
 
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &fluxX, Array< OneD, Array< OneD, NekDouble > > &fluxY)
 
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, const int j, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
 
SOLVER_UTILS_EXPORT void NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numflux)
 
SOLVER_UTILS_EXPORT void NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxX, Array< OneD, Array< OneD, NekDouble > > &numfluxY)
 
SOLVER_UTILS_EXPORT void NumFluxforScalar (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &uflux)
 
SOLVER_UTILS_EXPORT void NumFluxforVector (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
 
SOLVER_UTILS_EXPORT void SetModifiedBasis (const bool modbasis)
 
SOLVER_UTILS_EXPORT int NoCaseStringCompare (const string &s1, const string &s2)
 Perform a case-insensitive string comparison. More...
 
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)
 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)
 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 int i, const int j, const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &derivatives, Array< OneD, Array< OneD, NekDouble > > &flux)
 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::UnsteadySystem
SOLVER_UTILS_EXPORT UnsteadySystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 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 void v_NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numflux)
 
virtual SOLVER_UTILS_EXPORT void v_NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxX, Array< OneD, Array< OneD, NekDouble > > &numfluxY)
 
virtual SOLVER_UTILS_EXPORT void v_NumFluxforScalar (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &uflux)
 
virtual SOLVER_UTILS_EXPORT void v_NumFluxforVector (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
 
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_SteadyStateCheck (int step)
 
SOLVER_UTILS_EXPORT void CheckForRestartTime (NekDouble &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 time velocity. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 Initialises EquationSystem class members. More...
 
int nocase_cmp (const string &s1, const string &s2)
 
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)
 
SOLVER_UTILS_EXPORT void SetUpBaseFields (SpatialDomains::MeshGraphSharedPtr &mesh)
 
SOLVER_UTILS_EXPORT void ImportFldBase (std::string pInfile, SpatialDomains::MeshGraphSharedPtr pGraph)
 
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::ForcingSharedPtr
m_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...
 
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...
 
std::vector< int > m_intVariables
 
std::vector< FilterSharedPtrm_filters
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output. More...
 
map< std::string, Array< OneD,
Array< OneD, float > > > 
m_interpWeights
 Map of the interpolation weights for a specific filename. More...
 
map< std::string, Array< OneD,
Array< OneD, unsigned int > > > 
m_interpInds
 Map of the interpolation indices for a specific filename. More...
 
Array< OneD,
MultiRegions::ExpListSharedPtr
m_fields
 Array holding all dependent variables. More...
 
Array< OneD,
MultiRegions::ExpListSharedPtr
m_base
 Base fields. More...
 
Array< OneD,
MultiRegions::ExpListSharedPtr
m_derivedfields
 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...
 
std::set< std::string > m_loadedFields
 
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, Array< OneD,
Array< OneD, NekDouble > > > 
m_gradtan
 1 x nvariable x nq More...
 
Array< OneD, Array< OneD,
Array< OneD, NekDouble > > > 
m_tanbasis
 2 x m_spacedim x nq 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...
 
int m_NumMode
 Mode to use in case of single mode analysis. 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...
 

Detailed Description

Definition at line 49 of file UnsteadyViscousBurgers.h.

Constructor & Destructor Documentation

Nektar::UnsteadyViscousBurgers::~UnsteadyViscousBurgers ( )
virtual

Destructor.

Unsteady linear advection diffusion equation destructor.

Definition at line 172 of file UnsteadyViscousBurgers.cpp.

173  {
174  }
Nektar::UnsteadyViscousBurgers::UnsteadyViscousBurgers ( const LibUtilities::SessionReaderSharedPtr pSession)
protected

Session reader.

Definition at line 48 of file UnsteadyViscousBurgers.cpp.

References m_planeNumber.

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

Member Function Documentation

static SolverUtils::EquationSystemSharedPtr Nektar::UnsteadyViscousBurgers::create ( const LibUtilities::SessionReaderSharedPtr pSession)
inlinestatic

Creates an instance of this class.

Definition at line 55 of file UnsteadyViscousBurgers.h.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

56  {
59  AllocateSharedPtr(pSession);
60  p->InitObject();
61  return p;
62  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
boost::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.
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 315 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().

320  {
321  int nvariables = inarray.num_elements();
322  int nq = m_fields[0]->GetNpoints();
323 
325  factors[StdRegions::eFactorLambda] = 1.0/lambda/m_epsilon;
326 
327  if(m_useSpecVanVisc)
328  {
331  }
332 
333  Array<OneD, Array< OneD, NekDouble> > F(nvariables);
334  F[0] = Array<OneD, NekDouble> (nq*nvariables);
335 
336  for (int n = 1; n < nvariables; ++n)
337  {
338  F[n] = F[n-1] + nq;
339  }
340 
341  // We solve ( \nabla^2 - HHlambda ) Y[i] = rhs [i]
342  // inarray = input: \hat{rhs} -> output: \hat{Y}
343  // outarray = output: nabla^2 \hat{Y}
344  // where \hat = modal coeffs
345  for (int i = 0; i < nvariables; ++i)
346  {
347  // Multiply 1.0/timestep/lambda
348  Vmath::Smul(nq, -factors[StdRegions::eFactorLambda],
349  inarray[i], 1, F[i], 1);
350  }
351 
352  //Setting boundary conditions
353  SetBoundaryConditions(time);
354 
356  {
357  static int cnt = 0;
358 
359  if(cnt %10 == 0)
360  {
361  Array<OneD, Array<OneD, NekDouble> > vel(m_fields.num_elements());
362  for(int i = 0; i < m_fields.num_elements(); ++i)
363  {
364  m_fields[i]->ClearGlobalLinSysManager();
365  vel[i] = m_fields[i]->UpdatePhys();
366  }
368  }
369  ++cnt;
370  }
371  for (int i = 0; i < nvariables; ++i)
372  {
373  // Solve a system of equations with Helmholtz solver
374  m_fields[i]->HelmSolve(F[i], m_fields[i]->UpdateCoeffs(),
375  NullFlagList, factors, m_varCoeffLap);
376 
377  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(), outarray[i]);
378  }
379  }
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:251
StdRegions::VarCoeffMap m_varCoeffLap
Variable Coefficient map for the Laplacian which can be activated as part of SVV or otherwise...
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...
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:199
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.
static FlagList NullFlagList
An empty flag list.
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 273 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().

277  {
278  int i;
279  int nvariables = inarray.num_elements();
280  SetBoundaryConditions(time);
281  switch(m_projectionType)
282  {
284  {
285  // Just copy over array
286  int npoints = GetNpoints();
287 
288  for(i = 0; i < nvariables; ++i)
289  {
290  Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
291  }
292  break;
293  }
296  {
297  // Do nothing for the moment.
298  }
299  default:
300  {
301  ASSERTL0(false, "Unknown projection scheme");
302  break;
303  }
304  }
305  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
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:1047
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 216 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().

220  {
221  // Number of fields (variables of the problem)
222  int nVariables = inarray.num_elements();
223 
224  // Number of solution points
225  int nSolutionPts = GetNpoints();
226 
227  Array<OneD, Array<OneD, NekDouble> > outarrayDiff(nVariables);
228 
229  for (int i = 0; i < nVariables; ++i)
230  {
231  outarrayDiff[i] = Array<OneD, NekDouble>(nSolutionPts, 0.0);
232  }
233 
234  // RHS computation using the new advection base class
235  m_advObject->Advect(nVariables, m_fields, inarray,
236  inarray, outarray, time);
237 
238  // Negate the RHS
239  for (int i = 0; i < nVariables; ++i)
240  {
241  Vmath::Neg(nSolutionPts, outarray[i], 1);
242  }
243 
244  // No explicit diffusion for CG
246  {
247  m_diffusion->Diffuse(nVariables, m_fields, inarray, outarrayDiff);
248 
249  for (int i = 0; i < nVariables; ++i)
250  {
251  Vmath::Vadd(nSolutionPts, &outarray[i][0], 1,
252  &outarrayDiff[i][0], 1, &outarray[i][0], 1);
253  }
254  }
255 
256  // Add forcing terms
257  std::vector<SolverUtils::ForcingSharedPtr>::const_iterator x;
258  for (x = m_forcing.begin(); x != m_forcing.end(); ++x)
259  {
260  // set up non-linear terms
261  (*x)->Apply(m_fields, inarray, outarray, time);
262  }
263  }
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:382
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:285
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Forcing terms.
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 387 of file UnsteadyViscousBurgers.cpp.

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

Referenced by v_InitObject().

390  {
391 
392  const int nq = m_fields[0]->GetNpoints();
393 
394  for (int i = 0; i < flux.num_elements(); ++i)
395  {
396  for (int j = 0; j < flux[0].num_elements(); ++j)
397  {
398  Vmath::Vmul(nq, physfield[i], 1, physfield[j], 1,
399  flux[i][j], 1);
400  }
401  }
402  }
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:169
void Nektar::UnsteadyViscousBurgers::GetFluxVectorDiff ( const int  i,
const int  j,
const Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, NekDouble > > &  derivatives,
Array< OneD, Array< OneD, NekDouble > > &  flux 
)
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 413 of file UnsteadyViscousBurgers.cpp.

References Nektar::SolverUtils::EquationSystem::GetNpoints(), Vmath::Vcopy(), and Vmath::Zero().

Referenced by v_InitObject().

419  {
420  for (int k = 0; k < flux.num_elements(); ++k)
421  {
422  Vmath::Zero(GetNpoints(), flux[k], 1);
423  }
424  Vmath::Vcopy(GetNpoints(), physfield[i], 1, flux[j], 1);
425  }
SOLVER_UTILS_EXPORT int GetNpoints()
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:359
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
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 180 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().

182  {
183  // Number of trace (interface) points
184  int i;
185  int nTracePts = GetTraceNpoints();
186 
187  // Auxiliary variable to compute the normal velocity
188  Array<OneD, NekDouble> tmp(nTracePts);
189  m_traceVn = Array<OneD, NekDouble>(nTracePts, 0.0);
190 
191  // Reset the normal velocity
192  Vmath::Zero(nTracePts, m_traceVn, 1);
193 
194  for (i = 0; i < inarray.num_elements(); ++i)
195  {
196  m_fields[0]->ExtractTracePhys(inarray[i], tmp);
197 
198  Vmath::Vvtvp(nTracePts,
199  m_traceNormals[i], 1,
200  tmp, 1,
201  m_traceVn, 1,
202  m_traceVn, 1);
203  }
204 
205  return m_traceVn;
206  }
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:428
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:359
void Nektar::UnsteadyViscousBurgers::v_GenerateSummary ( SolverUtils::SummaryList s)
protectedvirtual

Print Summary.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 427 of file UnsteadyViscousBurgers.cpp.

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

429  {
430  AdvectionSystem::v_GenerateSummary(s);
431  if(m_useSpecVanVisc)
432  {
433  stringstream ss;
434  ss << "SVV (cut off = " << m_sVVCutoffRatio
435  << ", coeff = " << m_sVVDiffCoeff << ")";
436  AddSummaryItem(s, "Smoothing", ss.str());
437  }
438  }
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:50
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 61 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, and Nektar::SolverUtils::UnsteadySystem::SVVVarDiffCoeff().

62  {
63  AdvectionSystem::v_InitObject();
64 
65  m_session->LoadParameter("wavefreq", m_waveFreq, 0.0);
66  m_session->LoadParameter("epsilon", m_epsilon, 0.0);
67 
68  m_session->MatchSolverInfo(
69  "SpectralVanishingViscosity", "True", m_useSpecVanVisc, false);
70  m_session->MatchSolverInfo(
71  "SpectralVanishingViscosity", "VarDiff", m_useSpecVanViscVarDiff, false);
73  {
74  m_useSpecVanVisc = true;
75  }
76 
78  {
79  m_session->LoadParameter("SVVCutoffRatio",m_sVVCutoffRatio,0.75);
80  m_session->LoadParameter("SVVDiffCoeff",m_sVVDiffCoeff,0.1);
81  }
82 
83  // Type of advection and diffusion classes to be used
84  switch(m_projectionType)
85  {
86  // Discontinuous field
88  {
89  ASSERTL0(false,"Need to implement for DG");
90  // Do not forwards transform initial condition
91  m_homoInitialFwd = false;
92 
93  // Advection term
94  string advName;
95  string riemName;
96  m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
98  CreateInstance(advName, advName);
99  m_advObject->SetFluxVector(&UnsteadyViscousBurgers::
100  GetFluxVectorAdv, this);
101  m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
103  CreateInstance(riemName);
104  m_advObject->SetRiemannSolver(m_riemannSolver);
105  m_advObject->InitObject (m_session, m_fields);
106 
107  // Diffusion term
108  std::string diffName;
109  m_session->LoadSolverInfo("DiffusionType", diffName, "LDG");
111  CreateInstance(diffName, diffName);
112  m_diffusion->SetFluxVector(&UnsteadyViscousBurgers::
113  GetFluxVectorDiff, this);
114  m_diffusion->InitObject(m_session, m_fields);
115  break;
116  }
117  // Continuous field
120  {
121  // Advection term
122  std::string advName;
123  m_session->LoadSolverInfo("AdvectionType", advName,
124  "NonConservative");
126  CreateInstance(advName, advName);
127  m_advObject->SetFluxVector(&UnsteadyViscousBurgers::
128  GetFluxVectorAdv, this);
129 
131  {
132  Array<OneD, Array<OneD, NekDouble> > vel(m_fields.num_elements());
133  for(int i = 0; i < m_fields.num_elements(); ++i)
134  {
135  vel[i] = m_fields[i]->UpdatePhys();
136  }
138  }
139 
140  // In case of Galerkin explicit diffusion gives an error
142  {
143  ASSERTL0(false, "Explicit Galerkin diffusion not set up.");
144  }
145  // In case of Galerkin implicit diffusion: do nothing
146  break;
147  }
148  default:
149  {
150  ASSERTL0(false, "Unsupported projection type.");
151  break;
152  }
153  }
154 
155  // Forcing terms
157  m_fields.num_elements());
158 
161 
163  m_explicitDiffusion == 1)
164  {
166  }
167  }
void GetFluxVectorDiff(const int i, const int j, const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &derivatives, Array< OneD, Array< OneD, NekDouble > > &flux)
Evaluate the flux at each solution point for the diffusion part.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
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:42
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...
static SOLVER_UTILS_EXPORT std::vector< ForcingSharedPtr > Load(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields=0)
Definition: Forcing.cpp:84
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...
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
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:46
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.
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Forcing terms.

Friends And Related Function Documentation

friend class MemoryManager< UnsteadyViscousBurgers >
friend

Definition at line 52 of file UnsteadyViscousBurgers.h.

Member Data Documentation

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

Name of class.

Definition at line 64 of file UnsteadyViscousBurgers.h.

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

Definition at line 75 of file UnsteadyViscousBurgers.h.

Referenced by DoOdeRhs(), and v_InitObject().

NekDouble Nektar::UnsteadyViscousBurgers::m_epsilon
private

Definition at line 135 of file UnsteadyViscousBurgers.h.

Referenced by DoImplicitSolve(), and v_InitObject().

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

Forcing terms.

Definition at line 85 of file UnsteadyViscousBurgers.h.

Referenced by DoOdeRhs(), and v_InitObject().

int Nektar::UnsteadyViscousBurgers::m_planeNumber
protected

Definition at line 82 of file UnsteadyViscousBurgers.h.

Referenced by UnsteadyViscousBurgers().

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

Definition at line 74 of file UnsteadyViscousBurgers.h.

Referenced by v_InitObject().

NekDouble Nektar::UnsteadyViscousBurgers::m_sVVCutoffRatio
protected

Definition at line 72 of file UnsteadyViscousBurgers.h.

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

NekDouble Nektar::UnsteadyViscousBurgers::m_sVVDiffCoeff
protected

Definition at line 73 of file UnsteadyViscousBurgers.h.

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

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

Definition at line 76 of file UnsteadyViscousBurgers.h.

Referenced by GetNormalVelocity().

bool Nektar::UnsteadyViscousBurgers::m_useSpecVanVisc
protected

Definition at line 70 of file UnsteadyViscousBurgers.h.

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

bool Nektar::UnsteadyViscousBurgers::m_useSpecVanViscVarDiff
protected

Definition at line 71 of file UnsteadyViscousBurgers.h.

Referenced by DoImplicitSolve(), and v_InitObject().

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 79 of file UnsteadyViscousBurgers.h.

Referenced by DoImplicitSolve(), and v_InitObject().

NekDouble Nektar::UnsteadyViscousBurgers::m_waveFreq
private

Definition at line 134 of file UnsteadyViscousBurgers.h.

Referenced by v_InitObject().