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

#include <UnsteadyAdvectionDiffusion.h>

Inheritance diagram for Nektar::UnsteadyAdvectionDiffusion:
[legend]

Public Member Functions

 ~UnsteadyAdvectionDiffusion () override=default
 Destructor. More...
 
void v_ALEInitObject (int spaceDim, Array< OneD, MultiRegions::ExpListSharedPtr > &fields) override
 
- Public Member Functions inherited from Nektar::UnsteadyAdvection
 ~UnsteadyAdvection () override=default
 Destructor. More...
 
void v_ALEInitObject (int spaceDim, Array< OneD, MultiRegions::ExpListSharedPtr > &fields) override
 
- Public Member Functions inherited from Nektar::SolverUtils::AdvectionSystem
SOLVER_UTILS_EXPORT AdvectionSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
SOLVER_UTILS_EXPORT ~AdvectionSystem () override
 
SOLVER_UTILS_EXPORT void v_InitObject (bool DeclareField=true) override
 Initialisation object for EquationSystem. More...
 
SOLVER_UTILS_EXPORT AdvectionSharedPtr GetAdvObject ()
 Returns the advection object held by this instance. More...
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleGetElmtCFLVals (const bool FlagAcousticCFL=true)
 
SOLVER_UTILS_EXPORT NekDouble GetCFLEstimate (int &elmtid)
 
- Public Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
SOLVER_UTILS_EXPORT ~UnsteadySystem () override
 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 NekDouble GetTimeStep ()
 
SOLVER_UTILS_EXPORT void SetTimeStep (const NekDouble timestep)
 
SOLVER_UTILS_EXPORT void SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
SOLVER_UTILS_EXPORT LibUtilities::TimeIntegrationSchemeSharedPtrGetTimeIntegrationScheme ()
 Returns the time integration scheme. More...
 
SOLVER_UTILS_EXPORT LibUtilities::TimeIntegrationSchemeOperatorsGetTimeIntegrationSchemeOperators ()
 Returns the time integration scheme operators. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT void InitObject (bool DeclareField=true)
 Initialises the members of this object. More...
 
SOLVER_UTILS_EXPORT void DoInitialise (bool dumpInitialConditions=true)
 Perform any initialisation necessary before solving the problem. More...
 
SOLVER_UTILS_EXPORT void DoSolve ()
 Solve the problem. More...
 
SOLVER_UTILS_EXPORT void TransCoeffToPhys ()
 Transform from coefficient to physical space. More...
 
SOLVER_UTILS_EXPORT void TransPhysToCoeff ()
 Transform from physical to coefficient space. More...
 
SOLVER_UTILS_EXPORT void Output ()
 Perform output operations after solve. More...
 
SOLVER_UTILS_EXPORT std::string GetSessionName ()
 Get Session name. More...
 
template<class T >
std::shared_ptr< T > as ()
 
SOLVER_UTILS_EXPORT void ResetSessionName (std::string newname)
 Reset Session name. More...
 
SOLVER_UTILS_EXPORT LibUtilities::SessionReaderSharedPtr GetSession ()
 Get Session name. More...
 
SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr GetPressure ()
 Get pressure field if available. More...
 
SOLVER_UTILS_EXPORT void ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 
SOLVER_UTILS_EXPORT void PrintSummary (std::ostream &out)
 Print a summary of parameters and solver characteristics. More...
 
SOLVER_UTILS_EXPORT void SetLambda (NekDouble lambda)
 Set parameter m_lambda. More...
 
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction (std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
 Get a SessionFunction by name. More...
 
SOLVER_UTILS_EXPORT void SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 Initialise the data in the dependent fields. More...
 
SOLVER_UTILS_EXPORT void EvaluateExactSolution (int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 Evaluates an exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln, bool Normalised=false)
 Compute the L2 error between fields and a given exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, bool Normalised=false)
 Compute the L2 error of the fields. More...
 
SOLVER_UTILS_EXPORT NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation. More...
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleErrorExtraPoints (unsigned int field)
 Compute error (L2 and L_inf) over an larger set of quadrature points return [L2 Linf]. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n)
 Write checkpoint file of m_fields. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 Write checkpoint file of custom data fields. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_BaseFlow (const int n)
 Write base flow file of m_fields. More...
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname)
 Write field data to the given filename. More...
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 Write input fields to the given filename. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 Input field data from the given file. More...
 
SOLVER_UTILS_EXPORT void ImportFldToMultiDomains (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int ndomains)
 Input field data from the given file to multiple domains. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, std::vector< std::string > &fieldStr, Array< OneD, Array< OneD, NekDouble > > &coeffs)
 Output a field. Input field data into array from the given file. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, MultiRegions::ExpListSharedPtr &pField, std::string &pFieldName)
 Output a field. Input field data into ExpList from the given file. More...
 
SOLVER_UTILS_EXPORT void SessionSummary (SummaryList &vSummary)
 Write out a session summary. More...
 
SOLVER_UTILS_EXPORT Array< OneD, MultiRegions::ExpListSharedPtr > & UpdateFields ()
 
SOLVER_UTILS_EXPORT LibUtilities::FieldMetaDataMapUpdateFieldMetaDataMap ()
 Get hold of FieldInfoMap so it can be updated. More...
 
SOLVER_UTILS_EXPORT NekDouble GetTime ()
 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 Array< OneD, NekDouble > & UpdatePhysField (const int i)
 
SOLVER_UTILS_EXPORT void SetSteps (const int steps)
 
SOLVER_UTILS_EXPORT void ZeroPhysFields ()
 
SOLVER_UTILS_EXPORT void FwdTransFields ()
 
SOLVER_UTILS_EXPORT void SetModifiedBasis (const bool modbasis)
 
SOLVER_UTILS_EXPORT int GetCheckpointNumber ()
 
SOLVER_UTILS_EXPORT void SetCheckpointNumber (int num)
 
SOLVER_UTILS_EXPORT int GetCheckpointSteps ()
 
SOLVER_UTILS_EXPORT void SetCheckpointSteps (int num)
 
SOLVER_UTILS_EXPORT int GetInfoSteps ()
 
SOLVER_UTILS_EXPORT void SetInfoSteps (int num)
 
SOLVER_UTILS_EXPORT void SetIterationNumberPIT (int num)
 
SOLVER_UTILS_EXPORT void SetWindowNumberPIT (int num)
 
SOLVER_UTILS_EXPORT Array< OneD, const Array< OneD, NekDouble > > GetTraceNormals ()
 
SOLVER_UTILS_EXPORT void SetTime (const NekDouble time)
 
SOLVER_UTILS_EXPORT void SetTimeStep (const NekDouble timestep)
 
SOLVER_UTILS_EXPORT void SetInitialStep (const int step)
 
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time. More...
 
SOLVER_UTILS_EXPORT bool NegatedOp ()
 Identify if operator is negated in DoSolve. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::ALEHelper
virtual ~ALEHelper ()=default
 
virtual SOLVER_UTILS_EXPORT void v_ALEInitObject (int spaceDim, Array< OneD, MultiRegions::ExpListSharedPtr > &fields)
 
SOLVER_UTILS_EXPORT void InitObject (int spaceDim, Array< OneD, MultiRegions::ExpListSharedPtr > &fields)
 
virtual SOLVER_UTILS_EXPORT void v_UpdateGridVelocity (const NekDouble &time)
 
virtual SOLVER_UTILS_EXPORT void v_ALEPreMultiplyMass (Array< OneD, Array< OneD, NekDouble > > &fields)
 
SOLVER_UTILS_EXPORT void ALEDoElmtInvMass (Array< OneD, Array< OneD, NekDouble > > &traceNormals, Array< OneD, Array< OneD, NekDouble > > &fields, NekDouble time)
 Update m_fields with u^n by multiplying by inverse mass matrix. That's then used in e.g. checkpoint output and L^2 error calculation. More...
 
SOLVER_UTILS_EXPORT void ALEDoElmtInvMassBwdTrans (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
SOLVER_UTILS_EXPORT void MoveMesh (const NekDouble &time, Array< OneD, Array< OneD, NekDouble > > &traceNormals)
 
const Array< OneD, const Array< OneD, NekDouble > > & GetGridVelocity ()
 
SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & GetGridVelocityTrace ()
 
SOLVER_UTILS_EXPORT void ExtraFldOutputGridVelocity (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 

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 Member Functions inherited from Nektar::UnsteadyAdvection
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::UnsteadyAdvection
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

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

Protected Attributes

bool m_subSteppingScheme
 
bool m_useSpecVanVisc
 
NekDouble m_sVVCutoffRatio
 
NekDouble m_sVVDiffCoeff
 
SolverUtils::DiffusionSharedPtr m_diffusion
 
LibUtilities::TimeIntegrationSchemeSharedPtr m_subStepIntegrationScheme
 
LibUtilities::TimeIntegrationSchemeOperators m_subStepIntegrationOps
 
int m_intSteps
 
int m_minsubsteps
 
- Protected Attributes inherited from Nektar::UnsteadyAdvection
bool m_useGJPStabilisation
 
NekDouble m_GJPJumpScale
 
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
 
Array< OneD, Array< OneD, NekDouble > > m_velocity
 Advection velocity. More...
 
Array< OneD, NekDoublem_traceVn
 
std::vector< SolverUtils::ForcingSharedPtrm_forcing
 Forcing terms. More...
 
- Protected Attributes inherited from Nektar::SolverUtils::AdvectionSystem
SolverUtils::AdvectionSharedPtr m_advObject
 Advection term. More...
 
- Protected Attributes inherited from Nektar::SolverUtils::UnsteadySystem
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
Array< OneD, Array< OneD, NekDouble > > m_previousSolution
 Storage for previous solution for steady-state check. More...
 
std::vector< int > m_intVariables
 
NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1). More...
 
NekDouble m_CFLGrowth
 CFL growth rate. More...
 
NekDouble m_CFLEnd
 Maximun cfl in cfl growth. More...
 
int m_abortSteps
 Number of steps between checks for abort conditions. More...
 
bool m_explicitDiffusion
 Indicates if explicit or implicit treatment of diffusion is used. More...
 
bool m_explicitAdvection
 Indicates if explicit or implicit treatment of advection is used. More...
 
bool m_explicitReaction
 Indicates if explicit or implicit treatment of reaction is used. More...
 
int m_steadyStateSteps
 Check for steady state at step interval. More...
 
NekDouble m_steadyStateTol
 Tolerance to which steady state should be evaluated at. More...
 
int m_filtersInfosteps
 Number of time steps between outputting filters information. More...
 
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
 
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state. More...
 
std::ofstream m_errFile
 
NekDouble m_epsilon
 Diffusion coefficient. More...
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
bool m_verbose
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader. More...
 
std::map< std::string, SolverUtils::SessionFunctionSharedPtrm_sessionFunctions
 Map of known SessionFunctions. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 Array holding all dependent variables. More...
 
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
 Pointer to boundary conditions object. More...
 
SpatialDomains::MeshGraphSharedPtr m_graph
 Pointer to graph defining mesh. More...
 
std::string m_sessionName
 Name of the session. More...
 
NekDouble m_time
 Current time of simulation. More...
 
int m_initialStep
 Number of the step where the simulation should begin. More...
 
NekDouble m_fintime
 Finish time of the simulation. More...
 
NekDouble m_timestep
 Time step size. More...
 
NekDouble m_lambda
 Lambda constant in real system if one required. More...
 
NekDouble m_checktime
 Time between checkpoints. More...
 
NekDouble m_lastCheckTime
 
NekDouble m_TimeIncrementFactor
 
int m_nchk
 Number of checkpoints written so far. More...
 
int m_steps
 Number of steps to take. More...
 
int m_checksteps
 Number of steps between checkpoints. More...
 
int m_infosteps
 Number of time steps between outputting status information. More...
 
int m_iterPIT = 0
 Number of parallel-in-time time iteration. More...
 
int m_windowPIT = 0
 Index of windows for parallel-in-time time iteration. More...
 
int m_spacedim
 Spatial dimension (>= expansion dim). More...
 
int m_expdim
 Expansion dimension. More...
 
bool m_singleMode
 Flag to determine if single homogeneous mode is used. More...
 
bool m_halfMode
 Flag to determine if half homogeneous mode is used. More...
 
bool m_multipleModes
 Flag to determine if use multiple homogenenous modes are used. More...
 
bool m_useFFT
 Flag to determine if FFT is used for homogeneous transform. More...
 
bool m_homogen_dealiasing
 Flag to determine if dealiasing is used for homogeneous simulations. More...
 
bool m_specHP_dealiasing
 Flag to determine if dealisising is usde for the Spectral/hp element discretisation. More...
 
enum MultiRegions::ProjectionType m_projectionType
 Type of projection; e.g continuous or discontinuous. More...
 
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
 Array holding trace normals for DG simulations in the forwards direction. More...
 
Array< OneD, bool > m_checkIfSystemSingular
 Flag to indicate if the fields should be checked for singularity. More...
 
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
 Map to identify relevant solver info to dump in output fields. More...
 
Array< OneD, NekDoublem_movingFrameData
 Moving reference frame status in the inertial frame X, Y, Z, Theta_x, Theta_y, Theta_z, U, V, W, Omega_x, Omega_y, Omega_z, A_x, A_y, A_z, DOmega_x, DOmega_y, DOmega_z, pivot_x, pivot_y, pivot_z. More...
 
std::vector< std::string > m_strFrameData
 variable name in m_movingFrameData 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...
 
- Protected Attributes inherited from Nektar::SolverUtils::ALEHelper
Array< OneD, MultiRegions::ExpListSharedPtrm_fieldsALE
 
Array< OneD, Array< OneD, NekDouble > > m_gridVelocity
 
Array< OneD, Array< OneD, NekDouble > > m_gridVelocityTrace
 
std::vector< ALEBaseShPtrm_ALEs
 
bool m_ALESolver = false
 
bool m_ImplicitALESolver = false
 
NekDouble m_prevStageTime = 0.0
 
int m_spaceDim
 

Private Attributes

NekDouble m_epsilon
 

Friends

class MemoryManager< UnsteadyAdvectionDiffusion >
 

Additional Inherited Members

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

Detailed Description

Definition at line 43 of file UnsteadyAdvectionDiffusion.h.

Constructor & Destructor Documentation

◆ ~UnsteadyAdvectionDiffusion()

Nektar::UnsteadyAdvectionDiffusion::~UnsteadyAdvectionDiffusion ( )
overridedefault

Destructor.

◆ UnsteadyAdvectionDiffusion()

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

Session reader.

Definition at line 48 of file UnsteadyAdvectionDiffusion.cpp.

51 : UnsteadySystem(pSession, pGraph), UnsteadyAdvection(pSession, pGraph)
52{
53}
SOLVER_UTILS_EXPORT UnsteadySystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises UnsteadySystem class members.
UnsteadyAdvection(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Session reader.

Member Function Documentation

◆ AddAdvectionPenaltyFlux()

void Nektar::UnsteadyAdvectionDiffusion::AddAdvectionPenaltyFlux ( const Array< OneD, const Array< OneD, NekDouble > > &  velfield,
const Array< OneD, const Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, NekDouble > > &  Outarray 
)
protected

Number of trace points

Forward state array

Backward state array

upwind numerical flux state array

Normal velocity array

Extract forwards/backwards trace spaces Note: Needs to have correct i value to get boundary conditions

Upwind between elements

Construct difference between numflux and Fwd,Bwd

Calculate the numerical fluxes multipling Fwd, Bwd and numflux by the normal advection velocity

Definition at line 514 of file UnsteadyAdvectionDiffusion.cpp.

518{
519 ASSERTL1(physfield.size() == Outarray.size(),
520 "Physfield and outarray are of different dimensions");
521
522 int i;
523
524 /// Number of trace points
525 int nTracePts = m_fields[0]->GetTrace()->GetNpoints();
526
527 /// Forward state array
528 Array<OneD, NekDouble> Fwd(3 * nTracePts);
529
530 /// Backward state array
531 Array<OneD, NekDouble> Bwd = Fwd + nTracePts;
532
533 /// upwind numerical flux state array
534 Array<OneD, NekDouble> numflux = Bwd + nTracePts;
535
536 /// Normal velocity array
537 Array<OneD, NekDouble> Vn = GetNormalVel(velfield);
538
539 for (i = 0; i < physfield.size(); ++i)
540 {
541 /// Extract forwards/backwards trace spaces
542 /// Note: Needs to have correct i value to get boundary conditions
543 m_fields[i]->GetFwdBwdTracePhys(physfield[i], Fwd, Bwd);
544
545 /// Upwind between elements
546 m_fields[0]->GetTrace()->Upwind(Vn, Fwd, Bwd, numflux);
547
548 /// Construct difference between numflux and Fwd,Bwd
549 Vmath::Vsub(nTracePts, numflux, 1, Fwd, 1, Fwd, 1);
550 Vmath::Vsub(nTracePts, numflux, 1, Bwd, 1, Bwd, 1);
551
552 /// Calculate the numerical fluxes multipling Fwd, Bwd and
553 /// numflux by the normal advection velocity
554 Vmath::Vmul(nTracePts, Fwd, 1, Vn, 1, Fwd, 1);
555 Vmath::Vmul(nTracePts, Bwd, 1, Vn, 1, Bwd, 1);
556
557 m_fields[0]->AddFwdBwdTraceIntegral(Fwd, Bwd, Outarray[i]);
558 }
559}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
Array< OneD, NekDouble > & GetNormalVel(const Array< OneD, const Array< OneD, NekDouble > > &velfield)
Get the normal velocity based on input velfield.
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.hpp:72
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.hpp:220

References ASSERTL1, Nektar::UnsteadyAdvection::GetNormalVel(), Nektar::SolverUtils::EquationSystem::m_fields, Vmath::Vmul(), and Vmath::Vsub().

Referenced by SubStepAdvection().

◆ create()

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

Creates an instance of this class.

Definition at line 49 of file UnsteadyAdvectionDiffusion.h.

52 {
55 pSession, pGraph);
56 p->InitObject();
57 return p;
58 }
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::UnsteadyAdvectionDiffusion::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 209 of file UnsteadyAdvectionDiffusion.cpp.

213{
214 int nvariables = inarray.size();
215 int nq = m_fields[0]->GetNpoints();
216
219
220 // Add factor is GJP turned on
222 {
224 }
225
227 {
230 }
231
233 {
235 }
236
237 Array<OneD, Array<OneD, NekDouble>> F(nvariables);
238 for (int n = 0; n < nvariables; ++n)
239 {
240 F[n] = Array<OneD, NekDouble>(nq);
241 }
242
243 // We solve ( \nabla^2 - HHlambda ) Y[i] = rhs [i]
244 // inarray = input: \hat{rhs} -> output: \hat{Y}
245 // outarray = output: nabla^2 \hat{Y}
246 // where \hat = modal coeffs
247 for (int i = 0; i < nvariables; ++i)
248 {
249 // Multiply 1.0/timestep/lambda
250 Vmath::Smul(nq, -factors[StdRegions::eFactorLambda], inarray[i], 1,
251 F[i], 1);
252 }
253
254 // Setting boundary conditions
256
257 for (int i = 0; i < nvariables; ++i)
258 {
259 // Solve a system of equations with Helmholtz solver
260 m_fields[i]->HelmSolve(F[i], m_fields[i]->UpdateCoeffs(), factors);
261
262 m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(), outarray[i]);
263 }
264}
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:430
StdRegions::ConstFactorMap factors
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.hpp:100

References Nektar::MultiRegions::eDiscontinuous, Nektar::StdRegions::eFactorGJP, Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::eFactorSVVCutoffRatio, Nektar::StdRegions::eFactorSVVDiffCoeff, Nektar::StdRegions::eFactorTau, Nektar::VarcoeffHashingTest::factors, m_epsilon, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::UnsteadyAdvection::m_GJPJumpScale, Nektar::SolverUtils::EquationSystem::m_projectionType, m_sVVCutoffRatio, m_sVVDiffCoeff, Nektar::UnsteadyAdvection::m_useGJPStabilisation, m_useSpecVanVisc, Nektar::SolverUtils::EquationSystem::SetBoundaryConditions(), and Vmath::Smul().

Referenced by v_InitObject().

◆ DoOdeRhs()

void Nektar::UnsteadyAdvectionDiffusion::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 164 of file UnsteadyAdvectionDiffusion.cpp.

167{
168 // Number of fields (variables of the problem)
169 int nVariables = inarray.size();
170
171 UnsteadyAdvection::DoOdeRhs(inarray, outarray, time);
172
174 {
175 Array<OneD, Array<OneD, NekDouble>> outarrayDiff(nVariables);
176 for (int i = 0; i < nVariables; ++i)
177 {
178 outarrayDiff[i] = Array<OneD, NekDouble>(outarray[i].size(), 0.0);
179 }
180
181 if (m_ALESolver)
182 {
183 Array<OneD, Array<OneD, NekDouble>> tmpIn(nVariables);
184 // If ALE we must take Mu coefficient space to u physical space
186 m_diffusion->DiffuseCoeffs(nVariables, m_fields, tmpIn,
187 outarrayDiff); // Using tmpIn from above
188 }
189 else
190 {
191 m_diffusion->Diffuse(nVariables, m_fields, inarray, outarrayDiff);
192 }
193
194 for (int i = 0; i < nVariables; ++i)
195 {
196 Vmath::Vadd(outarray[i].size(), &outarrayDiff[i][0], 1,
197 &outarray[i][0], 1, &outarray[i][0], 1);
198 }
199 }
200}
SOLVER_UTILS_EXPORT void ALEDoElmtInvMassBwdTrans(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: ALEHelper.cpp:149
bool m_explicitDiffusion
Indicates if explicit or implicit treatment of diffusion is used.
SolverUtils::DiffusionSharedPtr m_diffusion
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the RHS.
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.hpp:180

References Nektar::SolverUtils::ALEHelper::ALEDoElmtInvMassBwdTrans(), Nektar::UnsteadyAdvection::DoOdeRhs(), Nektar::SolverUtils::ALEHelper::m_ALESolver, m_diffusion, Nektar::SolverUtils::UnsteadySystem::m_explicitDiffusion, Nektar::SolverUtils::EquationSystem::m_fields, and Vmath::Vadd().

Referenced by v_InitObject().

◆ GetFluxVectorDiff()

void Nektar::UnsteadyAdvectionDiffusion::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.

Definition at line 270 of file UnsteadyAdvectionDiffusion.cpp.

274{
275 unsigned int nDim = qfield.size();
276 unsigned int nConvectiveFields = qfield[0].size();
277 unsigned int nPts = qfield[0][0].size();
278 for (unsigned int j = 0; j < nDim; ++j)
279 {
280 for (unsigned int i = 0; i < nConvectiveFields; ++i)
281 {
282 Vmath::Smul(nPts, m_epsilon, qfield[j][i], 1, viscousTensor[j][i],
283 1);
284 }
285 }
286}

References m_epsilon, and Vmath::Smul().

Referenced by v_InitObject().

◆ GetMaxStdVelocity()

Array< OneD, NekDouble > Nektar::UnsteadyAdvectionDiffusion::GetMaxStdVelocity ( const Array< OneD, Array< OneD, NekDouble > >  inarray)
protected

Definition at line 561 of file UnsteadyAdvectionDiffusion.cpp.

563{
564
565 int n_points_0 = m_fields[0]->GetExp(0)->GetTotPoints();
566 int n_element = m_fields[0]->GetExpSize();
567 int nvel = inarray.size();
568 int cnt;
569
570 ASSERTL0(nvel >= 2, "Method not implemented for 1D");
571
572 NekDouble pntVelocity;
573
574 // Getting the standard velocity vector on the 2D normal space
575 Array<OneD, Array<OneD, NekDouble>> stdVelocity(nvel);
576 Array<OneD, NekDouble> maxV(n_element, 0.0);
578
579 for (int i = 0; i < nvel; ++i)
580 {
581 stdVelocity[i] = Array<OneD, NekDouble>(n_points_0);
582 }
583
584 if (nvel == 2)
585 {
586 cnt = 0.0;
587 for (int el = 0; el < n_element; ++el)
588 {
589 int n_points = m_fields[0]->GetExp(el)->GetTotPoints();
590 ptsKeys = m_fields[0]->GetExp(el)->GetPointsKeys();
591
592 // reset local space if necessary
593 if (n_points != n_points_0)
594 {
595 for (int j = 0; j < nvel; ++j)
596 {
597 stdVelocity[j] = Array<OneD, NekDouble>(n_points);
598 }
599 n_points_0 = n_points;
600 }
601
602 Array<TwoD, const NekDouble> gmat = m_fields[0]
603 ->GetExp(el)
604 ->GetGeom()
605 ->GetMetricInfo()
606 ->GetDerivFactors(ptsKeys);
607
608 if (m_fields[0]
609 ->GetExp(el)
610 ->GetGeom()
611 ->GetMetricInfo()
612 ->GetGtype() == SpatialDomains::eDeformed)
613 {
614 for (int i = 0; i < n_points; i++)
615 {
616 stdVelocity[0][i] = gmat[0][i] * inarray[0][i + cnt] +
617 gmat[2][i] * inarray[1][i + cnt];
618
619 stdVelocity[1][i] = gmat[1][i] * inarray[0][i + cnt] +
620 gmat[3][i] * inarray[1][i + cnt];
621 }
622 }
623 else
624 {
625 for (int i = 0; i < n_points; i++)
626 {
627 stdVelocity[0][i] = gmat[0][0] * inarray[0][i + cnt] +
628 gmat[2][0] * inarray[1][i + cnt];
629
630 stdVelocity[1][i] = gmat[1][0] * inarray[0][i + cnt] +
631 gmat[3][0] * inarray[1][i + cnt];
632 }
633 }
634
635 cnt += n_points;
636
637 for (int i = 0; i < n_points; i++)
638 {
639 pntVelocity = stdVelocity[0][i] * stdVelocity[0][i] +
640 stdVelocity[1][i] * stdVelocity[1][i];
641
642 if (pntVelocity > maxV[el])
643 {
644 maxV[el] = pntVelocity;
645 }
646 }
647 maxV[el] = sqrt(maxV[el]);
648 }
649 }
650 else
651 {
652 cnt = 0;
653 for (int el = 0; el < n_element; ++el)
654 {
655
656 int n_points = m_fields[0]->GetExp(el)->GetTotPoints();
657 ptsKeys = m_fields[0]->GetExp(el)->GetPointsKeys();
658
659 // reset local space if necessary
660 if (n_points != n_points_0)
661 {
662 for (int j = 0; j < nvel; ++j)
663 {
664 stdVelocity[j] = Array<OneD, NekDouble>(n_points);
665 }
666 n_points_0 = n_points;
667 }
668
669 Array<TwoD, const NekDouble> gmat = m_fields[0]
670 ->GetExp(el)
671 ->GetGeom()
672 ->GetMetricInfo()
673 ->GetDerivFactors(ptsKeys);
674
675 if (m_fields[0]
676 ->GetExp(el)
677 ->GetGeom()
678 ->GetMetricInfo()
679 ->GetGtype() == SpatialDomains::eDeformed)
680 {
681 for (int i = 0; i < n_points; i++)
682 {
683 stdVelocity[0][i] = gmat[0][i] * inarray[0][i + cnt] +
684 gmat[3][i] * inarray[1][i + cnt] +
685 gmat[6][i] * inarray[2][i + cnt];
686
687 stdVelocity[1][i] = gmat[1][i] * inarray[0][i + cnt] +
688 gmat[4][i] * inarray[1][i + cnt] +
689 gmat[7][i] * inarray[2][i + cnt];
690
691 stdVelocity[2][i] = gmat[2][i] * inarray[0][i + cnt] +
692 gmat[5][i] * inarray[1][i + cnt] +
693 gmat[8][i] * inarray[2][i + cnt];
694 }
695 }
696 else
697 {
698 for (int i = 0; i < n_points; i++)
699 {
700 stdVelocity[0][i] = gmat[0][0] * inarray[0][i + cnt] +
701 gmat[3][0] * inarray[1][i + cnt] +
702 gmat[6][0] * inarray[2][i + cnt];
703
704 stdVelocity[1][i] = gmat[1][0] * inarray[0][i + cnt] +
705 gmat[4][0] * inarray[1][i + cnt] +
706 gmat[7][0] * inarray[2][i + cnt];
707
708 stdVelocity[2][i] = gmat[2][0] * inarray[0][i + cnt] +
709 gmat[5][0] * inarray[1][i + cnt] +
710 gmat[8][0] * inarray[2][i + cnt];
711 }
712 }
713
714 cnt += n_points;
715
716 for (int i = 0; i < n_points; i++)
717 {
718 pntVelocity = stdVelocity[0][i] * stdVelocity[0][i] +
719 stdVelocity[1][i] * stdVelocity[1][i] +
720 stdVelocity[2][i] * stdVelocity[2][i];
721
722 if (pntVelocity > maxV[el])
723 {
724 maxV[el] = pntVelocity;
725 }
726 }
727
728 maxV[el] = sqrt(maxV[el]);
729 }
730 }
731
732 return maxV;
733}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:231
@ eDeformed
Geometry is curved or has non-constant factors.
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References ASSERTL0, Nektar::SpatialDomains::eDeformed, Nektar::SolverUtils::EquationSystem::m_fields, and tinysimd::sqrt().

Referenced by GetSubstepTimeStep().

◆ GetSubstepTimeStep()

NekDouble Nektar::UnsteadyAdvectionDiffusion::GetSubstepTimeStep ( )
protected

Definition at line 385 of file UnsteadyAdvectionDiffusion.cpp.

386{
387 int n_element = m_fields[0]->GetExpSize();
388
389 const Array<OneD, int> ExpOrder = m_fields[0]->EvalBasisNumModesMaxPerExp();
390 Array<OneD, int> ExpOrderList(n_element, ExpOrder);
391
392 const NekDouble cLambda = 0.2; // Spencer book pag. 317
393
394 Array<OneD, NekDouble> tstep(n_element, 0.0);
395 Array<OneD, NekDouble> stdVelocity(n_element, 0.0);
396
397 stdVelocity = GetMaxStdVelocity(m_velocity);
398
399 for (int el = 0; el < n_element; ++el)
400 {
401 tstep[el] =
402 m_cflSafetyFactor / (stdVelocity[el] * cLambda *
403 (ExpOrder[el] - 1) * (ExpOrder[el] - 1));
404 }
405
406 NekDouble TimeStep = Vmath::Vmin(n_element, tstep, 1);
407 m_session->GetComm()->AllReduce(TimeStep, LibUtilities::ReduceMin);
408
409 return TimeStep;
410}
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
NekDouble m_cflSafetyFactor
CFL safety factor (comprise between 0 to 1).
Array< OneD, NekDouble > GetMaxStdVelocity(const Array< OneD, Array< OneD, NekDouble > > inarray)
Array< OneD, Array< OneD, NekDouble > > m_velocity
Advection velocity.
T Vmin(int n, const T *x, const int incx)
Return the minimum element in x - called vmin to avoid conflict with min.
Definition: Vmath.hpp:725

References GetMaxStdVelocity(), Nektar::SolverUtils::UnsteadySystem::m_cflSafetyFactor, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_session, Nektar::UnsteadyAdvection::m_velocity, Nektar::LibUtilities::ReduceMin, and Vmath::Vmin().

Referenced by SubStepAdvance().

◆ SetUpSubSteppingTimeIntegration()

void Nektar::UnsteadyAdvectionDiffusion::SetUpSubSteppingTimeIntegration ( const LibUtilities::TimeIntegrationSchemeSharedPtr IntegrationScheme)
protected

Definition at line 412 of file UnsteadyAdvectionDiffusion.cpp.

414{
415 // Set to 1 for first step and it will then be increased in
416 // time advance routines
417 unsigned int order = IntegrationScheme->GetOrder();
418
419 // Set to 1 for first step and it will then be increased in
420 // time advance routines
421 if ((IntegrationScheme->GetName() == "Euler" &&
422 IntegrationScheme->GetVariant() == "Backward") ||
423 (IntegrationScheme->GetName() == "BDFImplicit" &&
424 (order == 1 || order == 2)))
425 {
426 // Note RK first order SSP is just Forward Euler.
429 "RungeKutta", "SSP", order, std::vector<NekDouble>());
430 }
431 else
432 {
434 "Integration method not suitable: "
435 "Options include BackwardEuler or BDFImplicitOrder1");
436 }
437
438 m_intSteps = IntegrationScheme->GetNumIntegrationPhases();
439
440 // set explicit time-integration class operators
445}
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
void SubStepAdvection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
void SubStepProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
LibUtilities::TimeIntegrationSchemeSharedPtr m_subStepIntegrationScheme
LibUtilities::TimeIntegrationSchemeOperators m_subStepIntegrationOps
TimeIntegrationSchemeFactory & GetTimeIntegrationSchemeFactory()

References Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineOdeRhs(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineProjection(), Nektar::ErrorUtil::efatal, Nektar::LibUtilities::GetTimeIntegrationSchemeFactory(), m_intSteps, m_subStepIntegrationOps, m_subStepIntegrationScheme, NEKERROR, SubStepAdvection(), and SubStepProjection().

Referenced by v_InitObject().

◆ SubStepAdvance()

void Nektar::UnsteadyAdvectionDiffusion::SubStepAdvance ( int  nstep,
NekDouble  time 
)
protected

Definition at line 329 of file UnsteadyAdvectionDiffusion.cpp.

330{
331 int nsubsteps;
332
333 NekDouble dt;
334
335 Array<OneD, Array<OneD, NekDouble>> fields, velfields;
336
337 static int ncalls = 1;
338 int nint = std::min(ncalls++, m_intSteps);
339
340 Array<OneD, NekDouble> CFL(m_fields[0]->GetExpSize(), m_cflSafetyFactor);
341
342 LibUtilities::CommSharedPtr comm = m_session->GetComm();
343
344 // Get the proper time step with CFL control
345 dt = GetSubstepTimeStep();
346
347 nsubsteps = (m_timestep > dt) ? ((int)(m_timestep / dt) + 1) : 1;
348 nsubsteps = std::max(m_minsubsteps, nsubsteps);
349
350 dt = m_timestep / nsubsteps;
351
352 if (m_infosteps && !((nstep + 1) % m_infosteps) && comm->GetRank() == 0)
353 {
354 std::cout << "Sub-integrating using " << nsubsteps
355 << " steps over Dt = " << m_timestep
356 << " (SubStep CFL=" << m_cflSafetyFactor << ")" << std::endl;
357 }
358
359 const TripleArray &solutionVector = m_intScheme->GetSolutionVector();
360
361 for (int m = 0; m < nint; ++m)
362 {
363 // We need to update the fields held by the m_intScheme
364 fields = solutionVector[m];
365
366 // Initialise NS solver which is set up to use a GLM method
367 // with calls to EvaluateAdvection_SetPressureBCs and
368 // SolveUnsteadyStokesSystem
369 m_subStepIntegrationScheme->InitializeScheme(dt, fields, time,
371
372 for (int n = 0; n < nsubsteps; ++n)
373 {
374 fields = m_subStepIntegrationScheme->TimeIntegrate(n, dt);
375 }
376
377 // Reset time integrated solution in m_intScheme
378 m_intScheme->SetSolutionVector(m, fields);
379 }
380}
NekDouble m_timestep
Time step size.
int m_infosteps
Number of time steps between outputting status information.
SOLVER_UTILS_EXPORT int GetExpSize()
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
Wrapper to the time integration scheme.
AT< AT< AT< NekDouble > > > TripleArray
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55

References Nektar::SolverUtils::EquationSystem::GetExpSize(), GetSubstepTimeStep(), Nektar::SolverUtils::UnsteadySystem::m_cflSafetyFactor, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_infosteps, Nektar::SolverUtils::UnsteadySystem::m_intScheme, m_intSteps, m_minsubsteps, Nektar::SolverUtils::EquationSystem::m_session, m_subStepIntegrationOps, m_subStepIntegrationScheme, and Nektar::SolverUtils::EquationSystem::m_timestep.

Referenced by v_PreIntegrate().

◆ SubStepAdvection()

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

Explicit Advection terms used by SubStepAdvance time integration

Get the number of coefficients

Define an auxiliary variable to compute the RHS

Operations to compute the RHS

Multiply the flux by the inverse of the mass matrix

Store in outarray the physical values of the RHS

Definition at line 450 of file UnsteadyAdvectionDiffusion.cpp.

453{
454 int i;
455 int nVariables = inarray.size();
456
457 /// Get the number of coefficients
458 int ncoeffs = m_fields[0]->GetNcoeffs();
459
460 /// Define an auxiliary variable to compute the RHS
461 Array<OneD, Array<OneD, NekDouble>> WeakAdv(nVariables);
462 WeakAdv[0] = Array<OneD, NekDouble>(ncoeffs * nVariables);
463 for (i = 1; i < nVariables; ++i)
464 {
465 WeakAdv[i] = WeakAdv[i - 1] + ncoeffs;
466 }
467
468 // Currently assume velocity field is time independent and does not
469 // therefore need extrapolating. RHS computation using the advection base
470 // class
471 m_advObject->Advect(nVariables, m_fields, m_velocity, inarray, outarray,
472 time);
473
474 for (i = 0; i < nVariables; ++i)
475 {
476 m_fields[i]->IProductWRTBase(outarray[i], WeakAdv[i]);
477 // negation requried due to sign of DoAdvection term to be consistent
478 Vmath::Neg(ncoeffs, WeakAdv[i], 1);
479 }
480
481 AddAdvectionPenaltyFlux(m_velocity, inarray, WeakAdv);
482
483 /// Operations to compute the RHS
484 for (i = 0; i < nVariables; ++i)
485 {
486 // Negate the RHS
487 Vmath::Neg(ncoeffs, WeakAdv[i], 1);
488
489 /// Multiply the flux by the inverse of the mass matrix
490 m_fields[i]->MultiplyByElmtInvMass(WeakAdv[i], WeakAdv[i]);
491
492 /// Store in outarray the physical values of the RHS
493 m_fields[i]->BwdTrans(WeakAdv[i], outarray[i]);
494 }
495}
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
void AddAdvectionPenaltyFlux(const Array< OneD, const Array< OneD, NekDouble > > &velfield, const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &Outarray)
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.hpp:292

References AddAdvectionPenaltyFlux(), Nektar::SolverUtils::AdvectionSystem::m_advObject, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::UnsteadyAdvection::m_velocity, and Vmath::Neg().

Referenced by SetUpSubSteppingTimeIntegration().

◆ SubStepProjection()

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

Projection used by SubStepAdvance time integration

Definition at line 500 of file UnsteadyAdvectionDiffusion.cpp.

504{
505 ASSERTL1(inarray.size() == outarray.size(),
506 "Inarray and outarray of different sizes ");
507
508 for (int i = 0; i < inarray.size(); ++i)
509 {
510 Vmath::Vcopy(inarray[i].size(), inarray[i], 1, outarray[i], 1);
511 }
512}
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

References ASSERTL1, and Vmath::Vcopy().

Referenced by SetUpSubSteppingTimeIntegration().

◆ v_ALEInitObject()

void Nektar::UnsteadyAdvectionDiffusion::v_ALEInitObject ( int  spaceDim,
Array< OneD, MultiRegions::ExpListSharedPtr > &  fields 
)
overridevirtual

Reimplemented from Nektar::UnsteadyAdvection.

Definition at line 735 of file UnsteadyAdvectionDiffusion.cpp.

737{
739 {
740 m_spaceDim = spaceDim;
741 m_fieldsALE = fields;
742
743 // Initialise grid velocities as 0s
744 m_gridVelocity = Array<OneD, Array<OneD, NekDouble>>(m_spaceDim);
745 m_gridVelocityTrace = Array<OneD, Array<OneD, NekDouble>>(m_spaceDim);
746 for (int i = 0; i < spaceDim; ++i)
747 {
748 m_gridVelocity[i] =
749 Array<OneD, NekDouble>(fields[0]->GetTotPoints(), 0.0);
750 m_gridVelocityTrace[i] = Array<OneD, NekDouble>(
751 fields[0]->GetTrace()->GetTotPoints(), 0.0);
752 }
753 }
754 ALEHelper::InitObject(spaceDim, fields);
755}
Array< OneD, MultiRegions::ExpListSharedPtr > m_fieldsALE
Definition: ALEHelper.h:89
SOLVER_UTILS_EXPORT void InitObject(int spaceDim, Array< OneD, MultiRegions::ExpListSharedPtr > &fields)
Definition: ALEHelper.cpp:48
Array< OneD, Array< OneD, NekDouble > > m_gridVelocityTrace
Definition: ALEHelper.h:91
Array< OneD, Array< OneD, NekDouble > > m_gridVelocity
Definition: ALEHelper.h:90
SOLVER_UTILS_EXPORT int GetTotPoints()

References Nektar::MultiRegions::eDiscontinuous, Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::SolverUtils::ALEHelper::InitObject(), Nektar::SolverUtils::ALEHelper::m_fieldsALE, Nektar::SolverUtils::ALEHelper::m_gridVelocity, Nektar::SolverUtils::ALEHelper::m_gridVelocityTrace, Nektar::SolverUtils::EquationSystem::m_projectionType, and Nektar::SolverUtils::ALEHelper::m_spaceDim.

◆ v_ExtraFldOutput()

void Nektar::UnsteadyAdvectionDiffusion::v_ExtraFldOutput ( std::vector< Array< OneD, NekDouble > > &  fieldcoeffs,
std::vector< std::string > &  variables 
)
overrideprotectedvirtual

Reimplemented from Nektar::UnsteadyAdvection.

Definition at line 300 of file UnsteadyAdvectionDiffusion.cpp.

303{
304 bool extraFields;
305 m_session->MatchSolverInfo("OutputExtraFields", "True", extraFields, true);
306
307 if (extraFields && m_ALESolver)
308 {
309 ExtraFldOutputGridVelocity(fieldcoeffs, variables);
310 }
311}
SOLVER_UTILS_EXPORT void ExtraFldOutputGridVelocity(std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
Definition: ALEHelper.cpp:392

References Nektar::SolverUtils::ALEHelper::ExtraFldOutputGridVelocity(), Nektar::SolverUtils::ALEHelper::m_ALESolver, and Nektar::SolverUtils::EquationSystem::m_session.

◆ v_GenerateSummary()

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

Print Summary.

Reimplemented from Nektar::UnsteadyAdvection.

Definition at line 288 of file UnsteadyAdvectionDiffusion.cpp.

289{
292 {
293 std::stringstream ss;
294 ss << "SVV (cut off = " << m_sVVCutoffRatio
295 << ", coeff = " << m_sVVDiffCoeff << ")";
296 SolverUtils::AddSummaryItem(s, "Smoothing", ss.str());
297 }
298}
void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print Summary.
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:47

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

◆ v_InitObject()

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

Initialise the object.

Initialisation object for the unsteady linear advection diffusion equation.

Reimplemented from Nektar::UnsteadyAdvection.

Definition at line 59 of file UnsteadyAdvectionDiffusion.cpp.

60{
62
63 m_session->LoadParameter("epsilon", m_epsilon, 1.0);
64
65 m_session->MatchSolverInfo("SpectralVanishingViscosity", "True",
66 m_useSpecVanVisc, false);
68 {
69 m_session->LoadParameter("SVVCutoffRatio", m_sVVCutoffRatio, 0.75);
70 m_session->LoadParameter("SVVDiffCoeff", m_sVVDiffCoeff, 0.1);
71 }
72
73 // turn on substepping
74 m_session->MatchSolverInfo("Extrapolation", "SubStepping",
75 m_subSteppingScheme, false);
76
77 m_session->LoadParameter("MinSubSteps", m_minsubsteps, 1);
78
79 // Type of advection and diffusion classes to be used
80 switch (m_projectionType)
81 {
82 // Discontinuous field
84 {
85 // Diffusion term
87 {
88 std::string diffName;
89 m_session->LoadSolverInfo("DiffusionType", diffName, "LDG");
91 diffName, diffName);
92 m_diffusion->SetFluxVector(
94 m_diffusion->InitObject(m_session, m_fields);
95 }
96
98 "SubSteppingScheme is not set up for DG projection");
99 break;
100 }
101 // Continuous field
104 {
105 // Advection term
106 std::string advName;
107 m_session->LoadSolverInfo("AdvectionType", advName,
108 "NonConservative");
109 if (advName.compare("WeakDG") == 0)
110 {
111 // Define the normal velocity fields
112 if (m_fields[0]->GetTrace())
113 {
114 m_traceVn = Array<OneD, NekDouble>(GetTraceNpoints());
115 }
116
117 std::string riemName;
118 m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
121 riemName, m_session);
122 m_riemannSolver->SetScalar(
124 m_advObject->SetRiemannSolver(m_riemannSolver);
125 m_advObject->InitObject(m_session, m_fields);
126 }
127
128 // In case of Galerkin explicit diffusion gives an error
130 {
131 ASSERTL0(false, "Explicit Galerkin diffusion not set up.");
132 }
133 // In case of Galerkin implicit diffusion: do nothing
134 break;
135 }
136 default:
137 {
138 ASSERTL0(false, "Unsupported projection type.");
139 break;
140 }
141 }
142
144 this);
146
147 if (m_subSteppingScheme) // Substepping
148 {
150 "Projection must be set to Mixed_CG_Discontinuous for "
151 "substepping");
153 }
154}
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
SOLVER_UTILS_EXPORT int GetTraceNpoints()
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the RHS.
void SetUpSubSteppingTimeIntegration(const LibUtilities::TimeIntegrationSchemeSharedPtr &IntegrationScheme)
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 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.
Array< OneD, NekDouble > m_traceVn
Array< OneD, NekDouble > & GetNormalVelocity()
Get the normal velocity.
void v_InitObject(bool DeclareFields=true) override
Initialise the object.
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:39
RiemannSolverFactory & GetRiemannSolverFactory()

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineImplicitSolve(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineOdeRhs(), DoImplicitSolve(), DoOdeRhs(), Nektar::MultiRegions::eDiscontinuous, Nektar::MultiRegions::eGalerkin, Nektar::MultiRegions::eMixed_CG_Discontinuous, Nektar::SolverUtils::GetDiffusionFactory(), GetFluxVectorDiff(), Nektar::UnsteadyAdvection::GetNormalVelocity(), Nektar::SolverUtils::GetRiemannSolverFactory(), Nektar::SolverUtils::EquationSystem::GetTraceNpoints(), Nektar::SolverUtils::AdvectionSystem::m_advObject, m_diffusion, m_epsilon, Nektar::SolverUtils::UnsteadySystem::m_explicitDiffusion, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::UnsteadySystem::m_intScheme, m_minsubsteps, Nektar::SolverUtils::UnsteadySystem::m_ode, Nektar::SolverUtils::EquationSystem::m_projectionType, Nektar::UnsteadyAdvection::m_riemannSolver, Nektar::SolverUtils::EquationSystem::m_session, m_subSteppingScheme, m_sVVCutoffRatio, m_sVVDiffCoeff, Nektar::UnsteadyAdvection::m_traceVn, m_useSpecVanVisc, SetUpSubSteppingTimeIntegration(), and Nektar::UnsteadyAdvection::v_InitObject().

◆ v_PreIntegrate()

bool Nektar::UnsteadyAdvectionDiffusion::v_PreIntegrate ( int  step)
overrideprotectedvirtual

PreIntegration step for substepping.

Perform the extrapolation.

Reimplemented from Nektar::UnsteadyAdvection.

Definition at line 316 of file UnsteadyAdvectionDiffusion.cpp.

317{
319 {
320 SubStepAdvance(step, m_time);
321 }
322
323 return false;
324}
NekDouble m_time
Current time of simulation.
void SubStepAdvance(int nstep, NekDouble time)

References m_subSteppingScheme, Nektar::SolverUtils::EquationSystem::m_time, and SubStepAdvance().

Friends And Related Function Documentation

◆ MemoryManager< UnsteadyAdvectionDiffusion >

friend class MemoryManager< UnsteadyAdvectionDiffusion >
friend

Definition at line 1 of file UnsteadyAdvectionDiffusion.h.

Member Data Documentation

◆ className

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

Name of class.

Definition at line 61 of file UnsteadyAdvectionDiffusion.h.

◆ m_diffusion

SolverUtils::DiffusionSharedPtr Nektar::UnsteadyAdvectionDiffusion::m_diffusion
protected

Definition at line 80 of file UnsteadyAdvectionDiffusion.h.

Referenced by DoOdeRhs(), and v_InitObject().

◆ m_epsilon

NekDouble Nektar::UnsteadyAdvectionDiffusion::m_epsilon
private

Definition at line 147 of file UnsteadyAdvectionDiffusion.h.

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

◆ m_intSteps

int Nektar::UnsteadyAdvectionDiffusion::m_intSteps
protected

◆ m_minsubsteps

int Nektar::UnsteadyAdvectionDiffusion::m_minsubsteps
protected

Definition at line 144 of file UnsteadyAdvectionDiffusion.h.

Referenced by SubStepAdvance(), and v_InitObject().

◆ m_subStepIntegrationOps

LibUtilities::TimeIntegrationSchemeOperators Nektar::UnsteadyAdvectionDiffusion::m_subStepIntegrationOps
protected

◆ m_subStepIntegrationScheme

LibUtilities::TimeIntegrationSchemeSharedPtr Nektar::UnsteadyAdvectionDiffusion::m_subStepIntegrationScheme
protected

◆ m_subSteppingScheme

bool Nektar::UnsteadyAdvectionDiffusion::m_subSteppingScheme
protected

Definition at line 71 of file UnsteadyAdvectionDiffusion.h.

Referenced by v_InitObject(), and v_PreIntegrate().

◆ m_sVVCutoffRatio

NekDouble Nektar::UnsteadyAdvectionDiffusion::m_sVVCutoffRatio
protected

Definition at line 76 of file UnsteadyAdvectionDiffusion.h.

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

◆ m_sVVDiffCoeff

NekDouble Nektar::UnsteadyAdvectionDiffusion::m_sVVDiffCoeff
protected

Definition at line 78 of file UnsteadyAdvectionDiffusion.h.

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

◆ m_useSpecVanVisc

bool Nektar::UnsteadyAdvectionDiffusion::m_useSpecVanVisc
protected

Definition at line 74 of file UnsteadyAdvectionDiffusion.h.

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