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

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

Static Public Member Functions

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

Static Public Attributes

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

Protected Member Functions

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

Protected Attributes

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

Private Attributes

NekDouble m_waveFreq
 
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 46 of file UnsteadyAdvectionDiffusion.h.

Constructor & Destructor Documentation

◆ ~UnsteadyAdvectionDiffusion()

Nektar::UnsteadyAdvectionDiffusion::~UnsteadyAdvectionDiffusion ( )
virtual

Destructor.

Unsteady linear advection diffusion equation destructor.

Definition at line 206 of file UnsteadyAdvectionDiffusion.cpp.

207{
208}

◆ UnsteadyAdvectionDiffusion()

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

Session reader.

Definition at line 50 of file UnsteadyAdvectionDiffusion.cpp.

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

References m_planeNumber.

Member Function Documentation

◆ 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 679 of file UnsteadyAdvectionDiffusion.cpp.

683{
684 ASSERTL1(physfield.size() == Outarray.size(),
685 "Physfield and outarray are of different dimensions");
686
687 int i;
688
689 /// Number of trace points
690 int nTracePts = m_fields[0]->GetTrace()->GetNpoints();
691
692 /// Forward state array
693 Array<OneD, NekDouble> Fwd(3 * nTracePts);
694
695 /// Backward state array
696 Array<OneD, NekDouble> Bwd = Fwd + nTracePts;
697
698 /// upwind numerical flux state array
699 Array<OneD, NekDouble> numflux = Bwd + nTracePts;
700
701 /// Normal velocity array
702 Array<OneD, NekDouble> Vn = GetNormalVel(velfield);
703
704 for (i = 0; i < physfield.size(); ++i)
705 {
706 /// Extract forwards/backwards trace spaces
707 /// Note: Needs to have correct i value to get boundary conditions
708 m_fields[i]->GetFwdBwdTracePhys(physfield[i], Fwd, Bwd);
709
710 /// Upwind between elements
711 m_fields[0]->GetTrace()->Upwind(Vn, Fwd, Bwd, numflux);
712
713 /// Construct difference between numflux and Fwd,Bwd
714 Vmath::Vsub(nTracePts, numflux, 1, Fwd, 1, Fwd, 1);
715 Vmath::Vsub(nTracePts, numflux, 1, Bwd, 1, Bwd, 1);
716
717 /// Calculate the numerical fluxes multipling Fwd, Bwd and
718 /// numflux by the normal advection velocity
719 Vmath::Vmul(nTracePts, Fwd, 1, Vn, 1, Fwd, 1);
720 Vmath::Vmul(nTracePts, Bwd, 1, Vn, 1, Bwd, 1);
721
722 m_fields[0]->AddFwdBwdTraceIntegral(Fwd, Bwd, Outarray[i]);
723 }
724}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
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.cpp:207
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.cpp:414

References ASSERTL1, 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 52 of file UnsteadyAdvectionDiffusion.h.

55 {
58 pSession, pGraph);
59 p->InitObject();
60 return p;
61 }
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 355 of file UnsteadyAdvectionDiffusion.cpp.

359{
360 int nvariables = inarray.size();
361 int nq = m_fields[0]->GetNpoints();
362
365
366 // Add factor is GJP turned on
368 {
370 }
371
373 {
376 }
378 {
380 }
381
382 Array<OneD, Array<OneD, NekDouble>> F(nvariables);
383 for (int n = 0; n < nvariables; ++n)
384 {
385 F[n] = Array<OneD, NekDouble>(nq);
386 }
387
388 // We solve ( \nabla^2 - HHlambda ) Y[i] = rhs [i]
389 // inarray = input: \hat{rhs} -> output: \hat{Y}
390 // outarray = output: nabla^2 \hat{Y}
391 // where \hat = modal coeffs
392 for (int i = 0; i < nvariables; ++i)
393 {
394 // Multiply 1.0/timestep/lambda
395 Vmath::Smul(nq, -factors[StdRegions::eFactorLambda], inarray[i], 1,
396 F[i], 1);
397 }
398
399 // Setting boundary conditions
401
402 for (int i = 0; i < nvariables; ++i)
403 {
404 // Solve a system of equations with Helmholtz solver
405 m_fields[i]->HelmSolve(F[i], m_fields[i]->UpdateCoeffs(), factors);
406
407 m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(), outarray[i]);
408 }
409}
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:408
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.cpp:245

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, m_GJPJumpScale, Nektar::SolverUtils::EquationSystem::m_projectionType, m_sVVCutoffRatio, m_sVVDiffCoeff, m_useGJPStabilisation, m_useSpecVanVisc, Nektar::SolverUtils::EquationSystem::SetBoundaryConditions(), and Vmath::Smul().

Referenced by v_InitObject().

◆ DoOdeProjection()

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

308{
309 int i;
310 int nvariables = inarray.size();
312 switch (m_projectionType)
313 {
315 {
316 // Just copy over array
317 if (inarray != outarray)
318 {
319 int npoints = GetNpoints();
320
321 for (i = 0; i < nvariables; ++i)
322 {
323 Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
324 }
325 }
326 break;
327 }
330 {
331 Array<OneD, NekDouble> coeffs(m_fields[0]->GetNcoeffs());
332
333 for (i = 0; i < nvariables; ++i)
334 {
335 m_fields[i]->FwdTrans(inarray[i], coeffs);
336 m_fields[i]->BwdTrans(coeffs, outarray[i]);
337 }
338 break;
339 }
340 default:
341 {
342 ASSERTL0(false, "Unknown projection scheme");
343 break;
344 }
345 }
346}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
SOLVER_UTILS_EXPORT int GetNpoints()
SOLVER_UTILS_EXPORT int GetNcoeffs()
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191

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

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 252 of file UnsteadyAdvectionDiffusion.cpp.

255{
256 // Number of fields (variables of the problem)
257 int nVariables = inarray.size();
258
259 // Number of solution points
260 int nSolutionPts = GetNpoints();
261
262 // RHS computation using the new advection base class
263 m_advObject->Advect(nVariables, m_fields, m_velocity, inarray, outarray,
264 time);
265
266 // Negate the RHS
267 for (int i = 0; i < nVariables; ++i)
268 {
269 Vmath::Neg(nSolutionPts, outarray[i], 1);
270 }
271
273 {
274 Array<OneD, Array<OneD, NekDouble>> outarrayDiff(nVariables);
275 for (int i = 0; i < nVariables; ++i)
276 {
277 outarrayDiff[i] = Array<OneD, NekDouble>(nSolutionPts, 0.0);
278 }
279
280 m_diffusion->Diffuse(nVariables, m_fields, inarray, outarrayDiff);
281
282 for (int i = 0; i < nVariables; ++i)
283 {
284 Vmath::Vadd(nSolutionPts, &outarrayDiff[i][0], 1, &outarray[i][0],
285 1, &outarray[i][0], 1);
286 }
287 }
288
289 // Add forcing terms
290 for (auto &x : m_forcing)
291 {
292 // set up non-linear terms
293 x->Apply(m_fields, inarray, outarray, time);
294 }
295}
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
bool m_explicitDiffusion
Indicates if explicit or implicit treatment of diffusion is used.
SolverUtils::DiffusionSharedPtr m_diffusion
Array< OneD, Array< OneD, NekDouble > > m_velocity
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Forcing terms.
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:513
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:354

References Nektar::SolverUtils::EquationSystem::GetNpoints(), Nektar::SolverUtils::AdvectionSystem::m_advObject, m_diffusion, Nektar::SolverUtils::UnsteadySystem::m_explicitDiffusion, Nektar::SolverUtils::EquationSystem::m_fields, m_forcing, m_velocity, Vmath::Neg(), and Vmath::Vadd().

Referenced by v_InitObject().

◆ GetFluxVectorAdv()

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

420{
421 ASSERTL1(flux[0].size() == m_velocity.size(),
422 "Dimension of flux array and velocity array do not match");
423
424 const int nq = m_fields[0]->GetNpoints();
425
426 for (int i = 0; i < flux.size(); ++i)
427 {
428 for (int j = 0; j < flux[0].size(); ++j)
429 {
430 Vmath::Vmul(nq, physfield[i], 1, m_velocity[j], 1, flux[i][j], 1);
431 }
432 }
433}

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

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.

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

Definition at line 444 of file UnsteadyAdvectionDiffusion.cpp.

448{
449 boost::ignore_unused(inarray);
450
451 unsigned int nDim = qfield.size();
452 unsigned int nConvectiveFields = qfield[0].size();
453 unsigned int nPts = qfield[0][0].size();
454 for (unsigned int j = 0; j < nDim; ++j)
455 {
456 for (unsigned int i = 0; i < nConvectiveFields; ++i)
457 {
458 Vmath::Smul(nPts, m_epsilon, qfield[j][i], 1, viscousTensor[j][i],
459 1);
460 }
461 }
462}

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 726 of file UnsteadyAdvectionDiffusion.cpp.

728{
729
730 int n_points_0 = m_fields[0]->GetExp(0)->GetTotPoints();
731 int n_element = m_fields[0]->GetExpSize();
732 int nvel = inarray.size();
733 int cnt;
734
735 ASSERTL0(nvel >= 2, "Method not implemented for 1D");
736
737 NekDouble pntVelocity;
738
739 // Getting the standard velocity vector on the 2D normal space
740 Array<OneD, Array<OneD, NekDouble>> stdVelocity(nvel);
741 Array<OneD, NekDouble> maxV(n_element, 0.0);
743
744 for (int i = 0; i < nvel; ++i)
745 {
746 stdVelocity[i] = Array<OneD, NekDouble>(n_points_0);
747 }
748
749 if (nvel == 2)
750 {
751 cnt = 0.0;
752 for (int el = 0; el < n_element; ++el)
753 {
754 int n_points = m_fields[0]->GetExp(el)->GetTotPoints();
755 ptsKeys = m_fields[0]->GetExp(el)->GetPointsKeys();
756
757 // reset local space if necessary
758 if (n_points != n_points_0)
759 {
760 for (int j = 0; j < nvel; ++j)
761 {
762 stdVelocity[j] = Array<OneD, NekDouble>(n_points);
763 }
764 n_points_0 = n_points;
765 }
766
767 Array<TwoD, const NekDouble> gmat = m_fields[0]
768 ->GetExp(el)
769 ->GetGeom()
770 ->GetMetricInfo()
771 ->GetDerivFactors(ptsKeys);
772
773 if (m_fields[0]
774 ->GetExp(el)
775 ->GetGeom()
776 ->GetMetricInfo()
777 ->GetGtype() == SpatialDomains::eDeformed)
778 {
779 for (int i = 0; i < n_points; i++)
780 {
781 stdVelocity[0][i] = gmat[0][i] * inarray[0][i + cnt] +
782 gmat[2][i] * inarray[1][i + cnt];
783
784 stdVelocity[1][i] = gmat[1][i] * inarray[0][i + cnt] +
785 gmat[3][i] * inarray[1][i + cnt];
786 }
787 }
788 else
789 {
790 for (int i = 0; i < n_points; i++)
791 {
792 stdVelocity[0][i] = gmat[0][0] * inarray[0][i + cnt] +
793 gmat[2][0] * inarray[1][i + cnt];
794
795 stdVelocity[1][i] = gmat[1][0] * inarray[0][i + cnt] +
796 gmat[3][0] * inarray[1][i + cnt];
797 }
798 }
799
800 cnt += n_points;
801
802 for (int i = 0; i < n_points; i++)
803 {
804 pntVelocity = stdVelocity[0][i] * stdVelocity[0][i] +
805 stdVelocity[1][i] * stdVelocity[1][i];
806
807 if (pntVelocity > maxV[el])
808 {
809 maxV[el] = pntVelocity;
810 }
811 }
812 maxV[el] = sqrt(maxV[el]);
813 }
814 }
815 else
816 {
817 cnt = 0;
818 for (int el = 0; el < n_element; ++el)
819 {
820
821 int n_points = m_fields[0]->GetExp(el)->GetTotPoints();
822 ptsKeys = m_fields[0]->GetExp(el)->GetPointsKeys();
823
824 // reset local space if necessary
825 if (n_points != n_points_0)
826 {
827 for (int j = 0; j < nvel; ++j)
828 {
829 stdVelocity[j] = Array<OneD, NekDouble>(n_points);
830 }
831 n_points_0 = n_points;
832 }
833
834 Array<TwoD, const NekDouble> gmat = m_fields[0]
835 ->GetExp(el)
836 ->GetGeom()
837 ->GetMetricInfo()
838 ->GetDerivFactors(ptsKeys);
839
840 if (m_fields[0]
841 ->GetExp(el)
842 ->GetGeom()
843 ->GetMetricInfo()
844 ->GetGtype() == SpatialDomains::eDeformed)
845 {
846 for (int i = 0; i < n_points; i++)
847 {
848 stdVelocity[0][i] = gmat[0][i] * inarray[0][i + cnt] +
849 gmat[3][i] * inarray[1][i + cnt] +
850 gmat[6][i] * inarray[2][i + cnt];
851
852 stdVelocity[1][i] = gmat[1][i] * inarray[0][i + cnt] +
853 gmat[4][i] * inarray[1][i + cnt] +
854 gmat[7][i] * inarray[2][i + cnt];
855
856 stdVelocity[2][i] = gmat[2][i] * inarray[0][i + cnt] +
857 gmat[5][i] * inarray[1][i + cnt] +
858 gmat[8][i] * inarray[2][i + cnt];
859 }
860 }
861 else
862 {
863 for (int i = 0; i < n_points; i++)
864 {
865 stdVelocity[0][i] = gmat[0][0] * inarray[0][i + cnt] +
866 gmat[3][0] * inarray[1][i + cnt] +
867 gmat[6][0] * inarray[2][i + cnt];
868
869 stdVelocity[1][i] = gmat[1][0] * inarray[0][i + cnt] +
870 gmat[4][0] * inarray[1][i + cnt] +
871 gmat[7][0] * inarray[2][i + cnt];
872
873 stdVelocity[2][i] = gmat[2][0] * inarray[0][i + cnt] +
874 gmat[5][0] * inarray[1][i + cnt] +
875 gmat[8][0] * inarray[2][i + cnt];
876 }
877 }
878
879 cnt += n_points;
880
881 for (int i = 0; i < n_points; i++)
882 {
883 pntVelocity = stdVelocity[0][i] * stdVelocity[0][i] +
884 stdVelocity[1][i] * stdVelocity[1][i] +
885 stdVelocity[2][i] * stdVelocity[2][i];
886
887 if (pntVelocity > maxV[el])
888 {
889 maxV[el] = pntVelocity;
890 }
891 }
892
893 maxV[el] = sqrt(maxV[el]);
894 // cout << maxV[el]*maxV[el] << endl;
895 }
896 }
897
898 return maxV;
899}
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:236
@ 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().

◆ GetNormalVel()

Array< OneD, NekDouble > & Nektar::UnsteadyAdvectionDiffusion::GetNormalVel ( const Array< OneD, const Array< OneD, NekDouble > > &  velfield)
protected

Get the normal velocity based on input velfield.

Definition at line 219 of file UnsteadyAdvectionDiffusion.cpp.

221{
222 // Number of trace (interface) points
223 int i;
224 int nTracePts = GetTraceNpoints();
225
226 // Auxiliary variable to compute the normal velocity
227 Array<OneD, NekDouble> tmp(nTracePts);
228 m_traceVn = Array<OneD, NekDouble>(nTracePts, 0.0);
229
230 // Reset the normal velocity
231 Vmath::Zero(nTracePts, m_traceVn, 1);
232
233 for (i = 0; i < velfield.size(); ++i)
234 {
235 m_fields[0]->ExtractTracePhys(velfield[i], tmp);
236
237 Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, tmp, 1, m_traceVn, 1,
238 m_traceVn, 1);
239 }
240
241 return m_traceVn;
242}
SOLVER_UTILS_EXPORT int GetTraceNpoints()
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:569
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:487

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

Referenced by AddAdvectionPenaltyFlux(), and GetNormalVelocity().

◆ GetNormalVelocity()

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

Get the normal velocity based on m_velocity.

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

Definition at line 214 of file UnsteadyAdvectionDiffusion.cpp.

215{
216 return GetNormalVel(m_velocity);
217}

References GetNormalVel(), and m_velocity.

Referenced by v_InitObject().

◆ GetSubstepTimeStep()

NekDouble Nektar::UnsteadyAdvectionDiffusion::GetSubstepTimeStep ( )
protected

Definition at line 549 of file UnsteadyAdvectionDiffusion.cpp.

550{
551 int n_element = m_fields[0]->GetExpSize();
552
553 const Array<OneD, int> ExpOrder = m_fields[0]->EvalBasisNumModesMaxPerExp();
554 Array<OneD, int> ExpOrderList(n_element, ExpOrder);
555
556 const NekDouble cLambda = 0.2; // Spencer book pag. 317
557
558 Array<OneD, NekDouble> tstep(n_element, 0.0);
559 Array<OneD, NekDouble> stdVelocity(n_element, 0.0);
560
561 stdVelocity = GetMaxStdVelocity(m_velocity);
562
563 for (int el = 0; el < n_element; ++el)
564 {
565 tstep[el] =
566 m_cflSafetyFactor / (stdVelocity[el] * cLambda *
567 (ExpOrder[el] - 1) * (ExpOrder[el] - 1));
568 }
569
570 NekDouble TimeStep = Vmath::Vmin(n_element, tstep, 1);
571 m_session->GetComm()->AllReduce(TimeStep, LibUtilities::ReduceMin);
572
573 return TimeStep;
574}
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
Array< OneD, NekDouble > GetMaxStdVelocity(const Array< OneD, Array< OneD, NekDouble > > inarray)
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.cpp:1045

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

Referenced by SubStepAdvance().

◆ SetUpSubSteppingTimeIntegration()

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

Definition at line 576 of file UnsteadyAdvectionDiffusion.cpp.

578{
579 // Set to 1 for first step and it will then be increased in
580 // time advance routines
581 unsigned int order = IntegrationScheme->GetOrder();
582
583 // Set to 1 for first step and it will then be increased in
584 // time advance routines
585 if ((IntegrationScheme->GetName() == "Euler" &&
586 IntegrationScheme->GetVariant() == "Backward") ||
587 (IntegrationScheme->GetName() == "BDFImplicit" &&
588 (order == 1 || order == 2)))
589 {
590 // Note RK first order SSP is just Forward Euler.
593 "RungeKutta", "SSP", order, std::vector<NekDouble>());
594 }
595 else
596 {
598 "Integration method not suitable: "
599 "Options include BackwardEuler or BDFImplicitOrder1");
600 }
601
602 m_intSteps = IntegrationScheme->GetNumIntegrationPhases();
603
604 // set explicit time-integration class operators
609}
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
void 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 492 of file UnsteadyAdvectionDiffusion.cpp.

493{
494 int n;
495 int nsubsteps;
496
497 NekDouble dt;
498
499 Array<OneD, Array<OneD, NekDouble>> fields, velfields;
500
501 static int ncalls = 1;
502 int nint = std::min(ncalls++, m_intSteps);
503
504 Array<OneD, NekDouble> CFL(m_fields[0]->GetExpSize(), m_cflSafetyFactor);
505
506 LibUtilities::CommSharedPtr comm = m_session->GetComm();
507
508 // Get the proper time step with CFL control
509 dt = GetSubstepTimeStep();
510
511 nsubsteps = (m_timestep > dt) ? ((int)(m_timestep / dt) + 1) : 1;
512 nsubsteps = std::max(m_minsubsteps, nsubsteps);
513
514 dt = m_timestep / nsubsteps;
515
516 if (m_infosteps && !((nstep + 1) % m_infosteps) && comm->GetRank() == 0)
517 {
518 std::cout << "Sub-integrating using " << nsubsteps
519 << " steps over Dt = " << m_timestep
520 << " (SubStep CFL=" << m_cflSafetyFactor << ")" << std::endl;
521 }
522
523 const TripleArray &solutionVector = m_intScheme->GetSolutionVector();
524
525 for (int m = 0; m < nint; ++m)
526 {
527 // We need to update the fields held by the m_intScheme
528 fields = solutionVector[m];
529
530 // Initialise NS solver which is set up to use a GLM method
531 // with calls to EvaluateAdvection_SetPressureBCs and
532 // SolveUnsteadyStokesSystem
533 m_subStepIntegrationScheme->InitializeScheme(dt, fields, time,
535
536 for (n = 0; n < nsubsteps; ++n)
537 {
538 fields = m_subStepIntegrationScheme->TimeIntegrate(n, dt);
539 }
540
541 // Reset time integrated solution in m_intScheme
542 m_intScheme->SetSolutionVector(m, fields);
543 }
544}
NekDouble m_timestep
Time step size.
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:57

References Nektar::SolverUtils::EquationSystem::GetExpSize(), GetSubstepTimeStep(), m_cflSafetyFactor, Nektar::SolverUtils::EquationSystem::m_fields, 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 614 of file UnsteadyAdvectionDiffusion.cpp.

617{
618 int i;
619 int nVariables = inarray.size();
620
621 /// Get the number of coefficients
622 int ncoeffs = m_fields[0]->GetNcoeffs();
623
624 /// Define an auxiliary variable to compute the RHS
625 Array<OneD, Array<OneD, NekDouble>> WeakAdv(nVariables);
626 WeakAdv[0] = Array<OneD, NekDouble>(ncoeffs * nVariables);
627 for (i = 1; i < nVariables; ++i)
628 {
629 WeakAdv[i] = WeakAdv[i - 1] + ncoeffs;
630 }
631
632 // Currently assume velocity field is time independent and does not
633 // therefore need extrapolating. RHS computation using the advection base
634 // class
635 m_advObject->Advect(nVariables, m_fields, m_velocity, inarray, outarray,
636 time);
637
638 for (i = 0; i < nVariables; ++i)
639 {
640 m_fields[i]->IProductWRTBase(outarray[i], WeakAdv[i]);
641 // negation requried due to sign of DoAdvection term to be consistent
642 Vmath::Neg(ncoeffs, WeakAdv[i], 1);
643 }
644
645 AddAdvectionPenaltyFlux(m_velocity, inarray, WeakAdv);
646
647 /// Operations to compute the RHS
648 for (i = 0; i < nVariables; ++i)
649 {
650 // Negate the RHS
651 Vmath::Neg(ncoeffs, WeakAdv[i], 1);
652
653 /// Multiply the flux by the inverse of the mass matrix
654 m_fields[i]->MultiplyByElmtInvMass(WeakAdv[i], WeakAdv[i]);
655
656 /// Store in outarray the physical values of the RHS
657 m_fields[i]->BwdTrans(WeakAdv[i], outarray[i]);
658 }
659}
void AddAdvectionPenaltyFlux(const Array< OneD, const Array< OneD, NekDouble > > &velfield, const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &Outarray)

References AddAdvectionPenaltyFlux(), Nektar::SolverUtils::AdvectionSystem::m_advObject, Nektar::SolverUtils::EquationSystem::m_fields, 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 664 of file UnsteadyAdvectionDiffusion.cpp.

667{
668 boost::ignore_unused(time);
669
670 ASSERTL1(inarray.size() == outarray.size(),
671 "Inarray and outarray of different sizes ");
672
673 for (int i = 0; i < inarray.size(); ++i)
674 {
675 Vmath::Vcopy(inarray[i].size(), inarray[i], 1, outarray[i], 1);
676 }
677}

References ASSERTL1, and Vmath::Vcopy().

Referenced by SetUpSubSteppingTimeIntegration().

◆ v_GenerateSummary()

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

Print Summary.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 464 of file UnsteadyAdvectionDiffusion.cpp.

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

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

◆ v_InitObject()

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

Initialise the object.

Initialisation object for the unsteady linear advection diffusion equation.

Reimplemented from Nektar::SolverUtils::AdvectionSystem.

Definition at line 62 of file UnsteadyAdvectionDiffusion.cpp.

63{
64 AdvectionSystem::v_InitObject(DeclareFields);
65
66 m_session->LoadParameter("wavefreq", m_waveFreq, 0.0);
67 m_session->LoadParameter("epsilon", m_epsilon, 1.0);
68
69 // turn on substepping
70 m_session->MatchSolverInfo("Extrapolation", "SubStepping",
71 m_subSteppingScheme, false);
72
73 // Define Velocity fields
74 m_velocity = Array<OneD, Array<OneD, NekDouble>>(m_spacedim);
75 std::vector<std::string> vel;
76 vel.push_back("Vx");
77 vel.push_back("Vy");
78 vel.push_back("Vz");
79 vel.resize(m_spacedim);
80
81 GetFunction("AdvectionVelocity")->Evaluate(vel, m_velocity);
82
83 m_session->MatchSolverInfo("SpectralVanishingViscosity", "True",
84 m_useSpecVanVisc, false);
85
86 // check to see if it is explicity turned off
87 m_session->MatchSolverInfo("GJPStabilisation", "False",
89
90 // if GJPStabilisation set to False bool will be true and
91 // if not false so negate/revese bool
93
94 m_session->LoadParameter("GJPJumpScale", m_GJPJumpScale, 1.0);
95
97 {
98 m_session->LoadParameter("SVVCutoffRatio", m_sVVCutoffRatio, 0.75);
99 m_session->LoadParameter("SVVDiffCoeff", m_sVVDiffCoeff, 0.1);
100 }
101
102 // Type of advection and diffusion classes to be used
103 switch (m_projectionType)
104 {
105 // Discontinuous field
107 {
108 // Do not forwards transform initial condition
109 m_homoInitialFwd = false;
110
111 // Advection term
112 std::string advName;
113 std::string riemName;
114 m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
116 advName, advName);
117 m_advObject->SetFluxVector(
119 m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
122 riemName, m_session);
123 m_riemannSolver->SetScalar(
125 m_advObject->SetRiemannSolver(m_riemannSolver);
126 m_advObject->InitObject(m_session, m_fields);
127
128 // Diffusion term
130 {
131 std::string diffName;
132 m_session->LoadSolverInfo("DiffusionType", diffName, "LDG");
134 diffName, diffName);
135 m_diffusion->SetFluxVector(
137 m_diffusion->InitObject(m_session, m_fields);
138 }
139
141 "SubSteppingScheme is not set up for DG projection");
142 break;
143 }
144 // Continuous field
147 {
148 // Advection term
149 std::string advName;
150 m_session->LoadSolverInfo("AdvectionType", advName,
151 "NonConservative");
153 advName, advName);
154 m_advObject->SetFluxVector(
156
157 if (advName.compare("WeakDG") == 0)
158 {
159 std::string riemName;
160 m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
163 riemName, m_session);
164 m_riemannSolver->SetScalar(
166 m_advObject->SetRiemannSolver(m_riemannSolver);
167 m_advObject->InitObject(m_session, m_fields);
168 }
169
170 // In case of Galerkin explicit diffusion gives an error
172 {
173 ASSERTL0(false, "Explicit Galerkin diffusion not set up.");
174 }
175 // In case of Galerkin implicit diffusion: do nothing
176 break;
177 }
178 default:
179 {
180 ASSERTL0(false, "Unsupported projection type.");
181 break;
182 }
183 }
184
185 // Forcing terms
186 m_forcing = SolverUtils::Forcing::Load(m_session, shared_from_this(),
187 m_fields, m_fields.size());
188
190 this);
193
194 if (m_subSteppingScheme) // Substepping
195 {
197 "Projection must be set to Mixed_CG_Discontinuous for "
198 "substepping");
200 }
201}
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
virtual SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareField=true) override
Init object for UnsteadySystem class.
int m_spacedim
Spatial dimension (>= expansion dim).
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
Get a SessionFunction by name.
static SOLVER_UTILS_EXPORT std::vector< ForcingSharedPtr > Load(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields=0)
Definition: Forcing.cpp:120
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
bool m_homoInitialFwd
Flag to determine if simulation should start in homogeneous forward transformed state.
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Perform the projection.
Array< OneD, NekDouble > & GetNormalVelocity()
Get the normal velocity based on m_velocity.
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 GetFluxVectorAdv(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
Evaluate the flux at each solution point for the advection part.
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
void GetFluxVectorDiff(const Array< OneD, Array< OneD, NekDouble > > &inarray, const Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &viscousTensor)
Evaluate the flux at each solution point for the diffusion part.
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:41
RiemannSolverFactory & GetRiemannSolverFactory()

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineImplicitSolve(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineOdeRhs(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineProjection(), DoImplicitSolve(), DoOdeProjection(), DoOdeRhs(), Nektar::MultiRegions::eDiscontinuous, Nektar::MultiRegions::eGalerkin, Nektar::MultiRegions::eMixed_CG_Discontinuous, Nektar::SolverUtils::GetAdvectionFactory(), Nektar::SolverUtils::GetDiffusionFactory(), GetFluxVectorAdv(), GetFluxVectorDiff(), Nektar::SolverUtils::EquationSystem::GetFunction(), GetNormalVelocity(), 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, m_GJPJumpScale, Nektar::SolverUtils::UnsteadySystem::m_homoInitialFwd, Nektar::SolverUtils::UnsteadySystem::m_intScheme, Nektar::SolverUtils::UnsteadySystem::m_ode, Nektar::SolverUtils::EquationSystem::m_projectionType, m_riemannSolver, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_spacedim, m_subSteppingScheme, m_sVVCutoffRatio, m_sVVDiffCoeff, m_useGJPStabilisation, m_useSpecVanVisc, m_velocity, m_waveFreq, SetUpSubSteppingTimeIntegration(), and Nektar::SolverUtils::AdvectionSystem::v_InitObject().

◆ v_PreIntegrate()

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

PreIntegration step for substepping.

Perform the extrapolation.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 479 of file UnsteadyAdvectionDiffusion.cpp.

480{
482 {
483 SubStepAdvance(step, m_time);
484 }
485
486 return false;
487}
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)
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
EquationSystemFactory & GetEquationSystemFactory()

Name of class.

Definition at line 63 of file UnsteadyAdvectionDiffusion.h.

◆ m_cflSafetyFactor

NekDouble Nektar::UnsteadyAdvectionDiffusion::m_cflSafetyFactor
protected

Definition at line 164 of file UnsteadyAdvectionDiffusion.h.

Referenced by GetSubstepTimeStep(), and SubStepAdvance().

◆ m_diffusion

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

Definition at line 78 of file UnsteadyAdvectionDiffusion.h.

Referenced by DoOdeRhs(), and v_InitObject().

◆ m_epsilon

NekDouble Nektar::UnsteadyAdvectionDiffusion::m_epsilon
private

Definition at line 170 of file UnsteadyAdvectionDiffusion.h.

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

◆ m_forcing

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

Forcing terms.

Definition at line 87 of file UnsteadyAdvectionDiffusion.h.

Referenced by DoOdeRhs(), and v_InitObject().

◆ m_GJPJumpScale

NekDouble Nektar::UnsteadyAdvectionDiffusion::m_GJPJumpScale
protected

Definition at line 73 of file UnsteadyAdvectionDiffusion.h.

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

◆ m_infosteps

int Nektar::UnsteadyAdvectionDiffusion::m_infosteps
protected

Definition at line 165 of file UnsteadyAdvectionDiffusion.h.

Referenced by SubStepAdvance().

◆ m_intSteps

int Nektar::UnsteadyAdvectionDiffusion::m_intSteps
protected

◆ m_minsubsteps

int Nektar::UnsteadyAdvectionDiffusion::m_minsubsteps
protected

Definition at line 166 of file UnsteadyAdvectionDiffusion.h.

Referenced by SubStepAdvance().

◆ m_planeNumber

int Nektar::UnsteadyAdvectionDiffusion::m_planeNumber
protected

Definition at line 84 of file UnsteadyAdvectionDiffusion.h.

Referenced by UnsteadyAdvectionDiffusion().

◆ m_riemannSolver

SolverUtils::RiemannSolverSharedPtr Nektar::UnsteadyAdvectionDiffusion::m_riemannSolver
protected

Definition at line 77 of file UnsteadyAdvectionDiffusion.h.

Referenced by 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 69 of file UnsteadyAdvectionDiffusion.h.

Referenced by v_InitObject(), and v_PreIntegrate().

◆ m_sVVCutoffRatio

NekDouble Nektar::UnsteadyAdvectionDiffusion::m_sVVCutoffRatio
protected

Definition at line 75 of file UnsteadyAdvectionDiffusion.h.

Referenced by DoImplicitSolve(), and v_InitObject().

◆ m_sVVDiffCoeff

NekDouble Nektar::UnsteadyAdvectionDiffusion::m_sVVDiffCoeff
protected

Definition at line 76 of file UnsteadyAdvectionDiffusion.h.

Referenced by DoImplicitSolve(), and v_InitObject().

◆ m_traceVn

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

Definition at line 80 of file UnsteadyAdvectionDiffusion.h.

Referenced by GetNormalVel().

◆ m_useGJPStabilisation

bool Nektar::UnsteadyAdvectionDiffusion::m_useGJPStabilisation
protected

Definition at line 71 of file UnsteadyAdvectionDiffusion.h.

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

◆ m_useSpecVanVisc

bool Nektar::UnsteadyAdvectionDiffusion::m_useSpecVanVisc
protected

Definition at line 70 of file UnsteadyAdvectionDiffusion.h.

Referenced by DoImplicitSolve(), and v_InitObject().

◆ m_velocity

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

◆ m_waveFreq

NekDouble Nektar::UnsteadyAdvectionDiffusion::m_waveFreq
private

Definition at line 169 of file UnsteadyAdvectionDiffusion.h.

Referenced by v_InitObject().