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

Static Public Member Functions

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

Static Public Attributes

static std::string className
 Name of class. More...
 

Protected Member Functions

 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...
 
virtual void DoOdeRhs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 Compute the RHS. More...
 
void DoOdeProjection (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 Perform the projection. More...
 
virtual void DoImplicitSolve (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, NekDouble time, NekDouble lambda)
 Solve implicitly the diffusion term. More...
 
Array< OneD, NekDouble > & GetNormalVelocity ()
 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 ()
 Initialise the object. More...
 
virtual void v_GenerateSummary (SolverUtils::SummaryList &s)
 Print Summary. More...
 
virtual bool v_PreIntegrate (int step)
 PreIntegration step for substepping. More...
 
void SubStepAdvance (const LibUtilities::TimeIntegrationSolutionSharedPtr &integrationSoln, int nstep, NekDouble time)
 
NekDouble GetSubstepTimeStep ()
 
void SetUpSubSteppingTimeIntegration (int intMethod, const LibUtilities::TimeIntegrationWrapperSharedPtr &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)
 
virtual SOLVER_UTILS_EXPORT Array< OneD, NekDoublev_GetMaxStdVelocity ()
 
- Protected Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
SOLVER_UTILS_EXPORT UnsteadySystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises UnsteadySystem class members. More...
 
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve ()
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise ()
 Sets up initial conditions. More...
 
virtual SOLVER_UTILS_EXPORT void v_AppendOutput1D (Array< OneD, Array< OneD, NekDouble >> &solution1D)
 Print the solution at each solution point in a txt file. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_GetTimeStep (const Array< OneD, const Array< OneD, NekDouble >> &inarray)
 Return the timestep to be used for the next step in the time-marching loop. More...
 
virtual SOLVER_UTILS_EXPORT bool v_RequireFwdTrans ()
 
SOLVER_UTILS_EXPORT void CheckForRestartTime (NekDouble &time, int &nchk)
 
SOLVER_UTILS_EXPORT void SVVVarDiffCoeff (const Array< OneD, Array< OneD, NekDouble >> vel, StdRegions::VarCoeffMap &varCoeffMap)
 Evaluate the SVV diffusion coefficient according to Moura's paper where it should proportional to h time velocity. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises EquationSystem class members. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Virtual function for the L_inf error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false)
 Virtual function for the L_2 error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransCoeffToPhys ()
 Virtual function for transformation to physical space. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransPhysToCoeff ()
 Virtual function for transformation to coefficient space. More...
 
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 
virtual SOLVER_UTILS_EXPORT void v_EvaluateExactSolution (unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 
virtual SOLVER_UTILS_EXPORT void v_Output (void)
 
virtual SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr v_GetPressure (void)
 
virtual SOLVER_UTILS_EXPORT void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 

Protected Attributes

bool m_subSteppingScheme
 
bool m_useSpecVanVisc
 
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
 
LibUtilities::TimeIntegrationWrapperSharedPtr 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
int m_infosteps
 Number of time steps between outputting status information. More...
 
int m_abortSteps
 Number of steps between checks for abort conditions. More...
 
int m_filtersInfosteps
 Number of time steps between outputting filters information. More...
 
int m_nanSteps
 
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
LibUtilities::TimeIntegrationSolutionSharedPtr m_intSoln
 
NekDouble m_epsilon
 
bool m_explicitDiffusion
 Indicates if explicit or implicit treatment of diffusion is used. More...
 
bool m_explicitAdvection
 Indicates if explicit or implicit treatment of advection is used. More...
 
bool m_explicitReaction
 Indicates if explicit or implicit treatment of reaction is used. More...
 
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state. More...
 
NekDouble m_steadyStateTol
 Tolerance to which steady state should be evaluated at. More...
 
int m_steadyStateSteps
 Check for steady state at step interval. More...
 
Array< OneD, Array< OneD, NekDouble > > m_previousSolution
 Storage for previous solution for steady-state check. More...
 
std::ofstream m_errFile
 
std::vector< int > m_intVariables
 
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
 
NekDouble m_filterTimeWarning
 Number of time steps between outputting status information. More...
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader. More...
 
std::map< std::string, SolverUtils::SessionFunctionSharedPtrm_sessionFunctions
 Map of known SessionFunctions. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 Array holding all dependent variables. More...
 
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
 Pointer to boundary conditions object. More...
 
SpatialDomains::MeshGraphSharedPtr m_graph
 Pointer to graph defining mesh. More...
 
std::string m_sessionName
 Name of the session. More...
 
NekDouble m_time
 Current time of simulation. More...
 
int m_initialStep
 Number of the step where the simulation should begin. More...
 
NekDouble m_fintime
 Finish time of the simulation. More...
 
NekDouble m_timestep
 Time step size. More...
 
NekDouble m_lambda
 Lambda constant in real system if one required. More...
 
NekDouble m_checktime
 Time between checkpoints. More...
 
int m_nchk
 Number of checkpoints written so far. More...
 
int m_steps
 Number of steps to take. More...
 
int m_checksteps
 Number of steps between checkpoints. More...
 
int m_spacedim
 Spatial dimension (>= expansion dim). More...
 
int m_expdim
 Expansion dimension. More...
 
bool m_singleMode
 Flag to determine if single homogeneous mode is used. More...
 
bool m_halfMode
 Flag to determine if half homogeneous mode is used. More...
 
bool m_multipleModes
 Flag to determine if use multiple homogenenous modes are used. More...
 
bool m_useFFT
 Flag to determine if FFT is used for homogeneous transform. More...
 
bool m_homogen_dealiasing
 Flag to determine if dealiasing is used for homogeneous simulations. More...
 
bool m_specHP_dealiasing
 Flag to determine if dealisising is usde for the Spectral/hp element discretisation. More...
 
enum MultiRegions::ProjectionType m_projectionType
 Type of projection; e.g continuous or discontinuous. More...
 
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
 Array holding trace normals for DG simulations in the forwards direction. More...
 
Array< OneD, bool > m_checkIfSystemSingular
 Flag to indicate if the fields should be checked for singularity. More...
 
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
 Map to identify relevant solver info to dump in output fields. More...
 
int m_NumQuadPointsError
 Number of Quadrature points used to work out the error. More...
 
enum HomogeneousType m_HomogeneousType
 
NekDouble m_LhomX
 physical length in X direction (if homogeneous) More...
 
NekDouble m_LhomY
 physical length in Y direction (if homogeneous) More...
 
NekDouble m_LhomZ
 physical length in Z direction (if homogeneous) More...
 
int m_npointsX
 number of points in X direction (if homogeneous) More...
 
int m_npointsY
 number of points in Y direction (if homogeneous) More...
 
int m_npointsZ
 number of points in Z direction (if homogeneous) More...
 
int m_HomoDirec
 number of homogenous directions More...
 

Private Attributes

NekDouble m_waveFreq
 
NekDouble m_epsilon
 

Friends

class MemoryManager< UnsteadyAdvectionDiffusion >
 

Additional Inherited Members

- Public Attributes inherited from Nektar::SolverUtils::UnsteadySystem
NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1). More...
 
- Protected Types inherited from Nektar::SolverUtils::EquationSystem
enum  HomogeneousType { eHomogeneous1D, eHomogeneous2D, eHomogeneous3D, eNotHomogeneous }
 Parameter for homogeneous expansions. More...
 
- Static Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
static std::string equationSystemTypeLookupIds []
 

Detailed Description

Definition at line 45 of file UnsteadyAdvectionDiffusion.h.

Constructor & Destructor Documentation

◆ ~UnsteadyAdvectionDiffusion()

Nektar::UnsteadyAdvectionDiffusion::~UnsteadyAdvectionDiffusion ( )
virtual

Destructor.

Unsteady linear advection diffusion equation destructor.

Definition at line 193 of file UnsteadyAdvectionDiffusion.cpp.

194  {
195  }

◆ UnsteadyAdvectionDiffusion()

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

Session reader.

Definition at line 50 of file UnsteadyAdvectionDiffusion.cpp.

References m_planeNumber.

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

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

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

Referenced by SubStepAdvection().

664  {
665  ASSERTL1(physfield.num_elements() == Outarray.num_elements(),
666  "Physfield and outarray are of different dimensions");
667 
668  int i;
669 
670  /// Number of trace points
671  int nTracePts = m_fields[0]->GetTrace()->GetNpoints();
672 
673  /// Forward state array
674  Array<OneD, NekDouble> Fwd(3*nTracePts);
675 
676  /// Backward state array
677  Array<OneD, NekDouble> Bwd = Fwd + nTracePts;
678 
679  /// upwind numerical flux state array
680  Array<OneD, NekDouble> numflux = Bwd + nTracePts;
681 
682  /// Normal velocity array
683  Array<OneD, NekDouble> Vn = GetNormalVel(velfield);
684 
685  for(i = 0; i < physfield.num_elements(); ++i)
686  {
687  /// Extract forwards/backwards trace spaces
688  /// Note: Needs to have correct i value to get boundary conditions
689  m_fields[i]->GetFwdBwdTracePhys(physfield[i], Fwd, Bwd);
690 
691  /// Upwind between elements
692  m_fields[0]->GetTrace()->Upwind(Vn, Fwd, Bwd, numflux);
693 
694  /// Construct difference between numflux and Fwd,Bwd
695  Vmath::Vsub(nTracePts, numflux, 1, Fwd, 1, Fwd, 1);
696  Vmath::Vsub(nTracePts, numflux, 1, Bwd, 1, Bwd, 1);
697 
698  /// Calculate the numerical fluxes multipling Fwd, Bwd and
699  /// numflux by the normal advection velocity
700  Vmath::Vmul(nTracePts, Fwd, 1, Vn, 1, Fwd, 1);
701  Vmath::Vmul(nTracePts, Bwd, 1, Vn, 1, Bwd, 1);
702 
703  m_fields[0]->AddFwdBwdTraceIntegral(Fwd,Bwd,Outarray[i]);
704  }
705  }
Array< OneD, NekDouble > & GetNormalVel(const Array< OneD, const Array< OneD, NekDouble > > &velfield)
Get the normal velocity based on input velfield.
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:346
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:186

◆ 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 51 of file UnsteadyAdvectionDiffusion.h.

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

54  {
57  AllocateSharedPtr(pSession, pGraph);
58  p->InitObject();
59  return p;
60  }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.

◆ DoImplicitSolve()

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

Solve implicitly the diffusion term.

Definition at line 340 of file UnsteadyAdvectionDiffusion.cpp.

References Nektar::MultiRegions::eDiscontinuous, Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::eFactorSVVCutoffRatio, Nektar::StdRegions::eFactorSVVDiffCoeff, Nektar::StdRegions::eFactorTau, m_epsilon, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_projectionType, m_sVVCutoffRatio, m_sVVDiffCoeff, m_useSpecVanVisc, Nektar::NullFlagList, Nektar::SolverUtils::EquationSystem::SetBoundaryConditions(), and Vmath::Smul().

Referenced by v_InitObject().

345  {
346  int nvariables = inarray.num_elements();
347  int nq = m_fields[0]->GetNpoints();
348 
350  factors[StdRegions::eFactorLambda] = 1.0/lambda/m_epsilon;
351 
352  if(m_useSpecVanVisc)
353  {
356  }
358  {
359  factors[StdRegions::eFactorTau] = 1.0;
360  }
361 
362  Array<OneD, Array< OneD, NekDouble> > F(nvariables);
363  for (int n = 0; n < nvariables; ++n)
364  {
365  F[n] = Array<OneD, NekDouble> (nq);
366  }
367 
368  // We solve ( \nabla^2 - HHlambda ) Y[i] = rhs [i]
369  // inarray = input: \hat{rhs} -> output: \hat{Y}
370  // outarray = output: nabla^2 \hat{Y}
371  // where \hat = modal coeffs
372  for (int i = 0; i < nvariables; ++i)
373  {
374  // Multiply 1.0/timestep/lambda
376  inarray[i], 1, F[i], 1);
377  }
378 
379  //Setting boundary conditions
380  SetBoundaryConditions(time);
381 
382  for (int i = 0; i < nvariables; ++i)
383  {
384  // Solve a system of equations with Helmholtz solver
385  m_fields[i]->HelmSolve(F[i], m_fields[i]->UpdateCoeffs(),
386  NullFlagList, factors);
387 
388  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(), outarray[i]);
389  }
390  }
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:294
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:216
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
static FlagList NullFlagList
An empty flag list.

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

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

295  {
296  int i;
297  int nvariables = inarray.num_elements();
298  SetBoundaryConditions(time);
299  switch(m_projectionType)
300  {
302  {
303  // Just copy over array
304  int npoints = GetNpoints();
305 
306  for(i = 0; i < nvariables; ++i)
307  {
308  Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
309  }
310  break;
311  }
314  {
315  Array<OneD, NekDouble> coeffs(m_fields[0]->GetNcoeffs());
316 
317  for(i = 0; i < nvariables; ++i)
318  {
319  m_fields[i]->FwdTrans(inarray[i], coeffs);
320  m_fields[i]->BwdTrans_IterPerExp(coeffs, outarray[i]);
321  }
322  break;
323  }
324  default:
325  {
326  ASSERTL0(false, "Unknown projection scheme");
327  break;
328  }
329  }
330  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
SOLVER_UTILS_EXPORT int GetNpoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetNcoeffs()
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064

◆ DoOdeRhs()

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

Compute the RHS.

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

Parameters
inarrayGiven fields.
outarrayCalculated solution.
timeTime.

Definition at line 243 of file UnsteadyAdvectionDiffusion.cpp.

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

Referenced by v_InitObject().

247  {
248  // Number of fields (variables of the problem)
249  int nVariables = inarray.num_elements();
250 
251  // Number of solution points
252  int nSolutionPts = GetNpoints();
253 
254  // RHS computation using the new advection base class
255  m_advObject->Advect(nVariables, m_fields, m_velocity,
256  inarray, outarray, time);
257 
258  // Negate the RHS
259  for (int i = 0; i < nVariables; ++i)
260  {
261  Vmath::Neg(nSolutionPts, outarray[i], 1);
262  }
263 
265  {
266  Array<OneD, Array<OneD, NekDouble> > outarrayDiff(nVariables);
267  for (int i = 0; i < nVariables; ++i)
268  {
269  outarrayDiff[i] = Array<OneD, NekDouble>(nSolutionPts, 0.0);
270  }
271 
272  m_diffusion->Diffuse(nVariables, m_fields, inarray, outarrayDiff);
273 
274  for (int i = 0; i < nVariables; ++i)
275  {
276  Vmath::Vadd(nSolutionPts, &outarrayDiff[i][0], 1,
277  &outarray[i][0], 1, &outarray[i][0], 1);
278  }
279  }
280 
281  }
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
bool m_explicitDiffusion
Indicates if explicit or implicit treatment of diffusion is used.
Array< OneD, Array< OneD, NekDouble > > m_velocity
SolverUtils::DiffusionSharedPtr m_diffusion
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:399
SOLVER_UTILS_EXPORT int GetNpoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:302

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

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

Referenced by v_InitObject().

401  {
402  ASSERTL1(flux[0].num_elements() == m_velocity.num_elements(),
403  "Dimension of flux array and velocity array do not match");
404 
405  const int nq = m_fields[0]->GetNpoints();
406 
407  for (int i = 0; i < flux.num_elements(); ++i)
408  {
409  for (int j = 0; j < flux[0].num_elements(); ++j)
410  {
411  Vmath::Vmul(nq, physfield[i], 1, m_velocity[j], 1,
412  flux[i][j], 1);
413  }
414  }
415  }
Array< OneD, Array< OneD, NekDouble > > m_velocity
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:186

◆ GetFluxVectorDiff()

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

References m_epsilon, and Vmath::Smul().

Referenced by v_InitObject().

430  {
431  boost::ignore_unused(inarray);
432 
433  unsigned int nDim = qfield.num_elements();
434  unsigned int nConvectiveFields = qfield[0].num_elements();
435  unsigned int nPts = qfield[0][0].num_elements();
436  for (unsigned int j = 0; j < nDim; ++j)
437  {
438  for (unsigned int i = 0; i < nConvectiveFields; ++i)
439  {
440  Vmath::Smul(nPts, m_epsilon, qfield[j][i], 1,
441  viscousTensor[j][i], 1 );
442  }
443  }
444  }
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:216

◆ GetMaxStdVelocity()

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

Definition at line 708 of file UnsteadyAdvectionDiffusion.cpp.

References ASSERTL0, Nektar::SpatialDomains::eDeformed, and Nektar::SolverUtils::EquationSystem::m_fields.

Referenced by GetSubstepTimeStep().

710  {
711 
712  int n_points_0 = m_fields[0]->GetExp(0)->GetTotPoints();
713  int n_element = m_fields[0]->GetExpSize();
714  int nvel = inarray.num_elements();
715  int cnt;
716 
717  ASSERTL0(nvel >= 2, "Method not implemented for 1D");
718 
719  NekDouble pntVelocity;
720 
721  // Getting the standard velocity vector on the 2D normal space
722  Array<OneD, Array<OneD, NekDouble> > stdVelocity(nvel);
723  Array<OneD, NekDouble> maxV(n_element, 0.0);
725 
726  for (int i = 0; i < nvel; ++i)
727  {
728  stdVelocity[i] = Array<OneD, NekDouble>(n_points_0);
729  }
730 
731  if (nvel == 2)
732  {
733  cnt = 0.0;
734  for (int el = 0; el < n_element; ++el)
735  {
736  int n_points = m_fields[0]->GetExp(el)->GetTotPoints();
737  ptsKeys = m_fields[0]->GetExp(el)->GetPointsKeys();
738 
739  // reset local space if necessary
740  if(n_points != n_points_0)
741  {
742  for (int j = 0; j < nvel; ++j)
743  {
744  stdVelocity[j] = Array<OneD, NekDouble>(n_points);
745  }
746  n_points_0 = n_points;
747  }
748 
749  Array<TwoD, const NekDouble> gmat =
750  m_fields[0]->GetExp(el)->GetGeom()->GetMetricInfo()->GetDerivFactors(ptsKeys);
751 
752  if (m_fields[0]->GetExp(el)->GetGeom()->GetMetricInfo()->GetGtype()
754  {
755  for (int i = 0; i < n_points; i++)
756  {
757  stdVelocity[0][i] = gmat[0][i]*inarray[0][i+cnt]
758  + gmat[2][i]*inarray[1][i+cnt];
759 
760  stdVelocity[1][i] = gmat[1][i]*inarray[0][i+cnt]
761  + gmat[3][i]*inarray[1][i+cnt];
762  }
763  }
764  else
765  {
766  for (int i = 0; i < n_points; i++)
767  {
768  stdVelocity[0][i] = gmat[0][0]*inarray[0][i+cnt]
769  + gmat[2][0]*inarray[1][i+cnt];
770 
771  stdVelocity[1][i] = gmat[1][0]*inarray[0][i+cnt]
772  + gmat[3][0]*inarray[1][i+cnt];
773  }
774  }
775 
776  cnt += n_points;
777 
778 
779  for (int i = 0; i < n_points; i++)
780  {
781  pntVelocity = stdVelocity[0][i]*stdVelocity[0][i]
782  + stdVelocity[1][i]*stdVelocity[1][i];
783 
784  if (pntVelocity>maxV[el])
785  {
786  maxV[el] = pntVelocity;
787  }
788  }
789  maxV[el] = sqrt(maxV[el]);
790  }
791  }
792  else
793  {
794  cnt = 0;
795  for (int el = 0; el < n_element; ++el)
796  {
797 
798  int n_points = m_fields[0]->GetExp(el)->GetTotPoints();
799  ptsKeys = m_fields[0]->GetExp(el)->GetPointsKeys();
800 
801  // reset local space if necessary
802  if(n_points != n_points_0)
803  {
804  for (int j = 0; j < nvel; ++j)
805  {
806  stdVelocity[j] = Array<OneD, NekDouble>(n_points);
807  }
808  n_points_0 = n_points;
809  }
810 
811  Array<TwoD, const NekDouble> gmat =
812  m_fields[0]->GetExp(el)->GetGeom()->GetMetricInfo()->GetDerivFactors(ptsKeys);
813 
814  if (m_fields[0]->GetExp(el)->GetGeom()->GetMetricInfo()->GetGtype()
816  {
817  for (int i = 0; i < n_points; i++)
818  {
819  stdVelocity[0][i] = gmat[0][i]*inarray[0][i+cnt]
820  + gmat[3][i]*inarray[1][i+cnt]
821  + gmat[6][i]*inarray[2][i+cnt];
822 
823  stdVelocity[1][i] = gmat[1][i]*inarray[0][i+cnt]
824  + gmat[4][i]*inarray[1][i+cnt]
825  + gmat[7][i]*inarray[2][i+cnt];
826 
827  stdVelocity[2][i] = gmat[2][i]*inarray[0][i+cnt]
828  + gmat[5][i]*inarray[1][i+cnt]
829  + gmat[8][i]*inarray[2][i+cnt];
830  }
831  }
832  else
833  {
834  for (int i = 0; i < n_points; i++)
835  {
836  stdVelocity[0][i] = gmat[0][0]*inarray[0][i+cnt]
837  + gmat[3][0]*inarray[1][i+cnt]
838  + gmat[6][0]*inarray[2][i+cnt];
839 
840  stdVelocity[1][i] = gmat[1][0]*inarray[0][i+cnt]
841  + gmat[4][0]*inarray[1][i+cnt]
842  + gmat[7][0]*inarray[2][i+cnt];
843 
844  stdVelocity[2][i] = gmat[2][0]*inarray[0][i+cnt]
845  + gmat[5][0]*inarray[1][i+cnt]
846  + gmat[8][0]*inarray[2][i+cnt];
847  }
848  }
849 
850  cnt += n_points;
851 
852  for (int i = 0; i < n_points; i++)
853  {
854  pntVelocity = stdVelocity[0][i]*stdVelocity[0][i]
855  + stdVelocity[1][i]*stdVelocity[1][i]
856  + stdVelocity[2][i]*stdVelocity[2][i];
857 
858  if (pntVelocity > maxV[el])
859  {
860  maxV[el] = pntVelocity;
861  }
862  }
863 
864  maxV[el] = sqrt(maxV[el]);
865  //cout << maxV[el]*maxV[el] << endl;
866  }
867  }
868 
869  return maxV;
870  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:246
double NekDouble
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
Geometry is curved or has non-constant factors.

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

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

209  {
210  // Number of trace (interface) points
211  int i;
212  int nTracePts = GetTraceNpoints();
213 
214  // Auxiliary variable to compute the normal velocity
215  Array<OneD, NekDouble> tmp(nTracePts);
216  m_traceVn = Array<OneD, NekDouble>(nTracePts, 0.0);
217 
218  // Reset the normal velocity
219  Vmath::Zero(nTracePts, m_traceVn, 1);
220 
221  for (i = 0; i < velfield.num_elements(); ++i)
222  {
223  m_fields[0]->ExtractTracePhys(velfield[i], tmp);
224 
225  Vmath::Vvtvp(nTracePts,
226  m_traceNormals[i], 1,
227  tmp, 1,
228  m_traceVn, 1,
229  m_traceVn, 1);
230  }
231 
232  return m_traceVn;
233  }
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:445
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetTraceNpoints()
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:376

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

References GetNormalVel(), and m_velocity.

Referenced by v_InitObject().

202  {
203  return GetNormalVel(m_velocity);
204  }
Array< OneD, Array< OneD, NekDouble > > m_velocity
Array< OneD, NekDouble > & GetNormalVel(const Array< OneD, const Array< OneD, NekDouble > > &velfield)
Get the normal velocity based on input velfield.

◆ GetSubstepTimeStep()

NekDouble Nektar::UnsteadyAdvectionDiffusion::GetSubstepTimeStep ( )
protected

Definition at line 532 of file UnsteadyAdvectionDiffusion.cpp.

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

533  {
534  int n_element = m_fields[0]->GetExpSize();
535 
536  const Array<OneD, int> ExpOrder=m_fields[0]->EvalBasisNumModesMaxPerExp();
537  Array<OneD, int> ExpOrderList (n_element, ExpOrder);
538 
539  const NekDouble cLambda = 0.2; // Spencer book pag. 317
540 
541  Array<OneD, NekDouble> tstep (n_element, 0.0);
542  Array<OneD, NekDouble> stdVelocity(n_element, 0.0);
543 
544  stdVelocity = GetMaxStdVelocity(m_velocity);
545 
546  for(int el = 0; el < n_element; ++el)
547  {
548  tstep[el] = m_cflSafetyFactor /
549  (stdVelocity[el] * cLambda *
550  (ExpOrder[el]-1) * (ExpOrder[el]-1));
551  }
552 
553  NekDouble TimeStep = Vmath::Vmin(n_element, tstep, 1);
554  m_session->GetComm()->AllReduce(TimeStep,LibUtilities::ReduceMin);
555 
556  return TimeStep;
557  }
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:874
Array< OneD, Array< OneD, NekDouble > > m_velocity
double NekDouble
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.

◆ SetUpSubSteppingTimeIntegration()

void Nektar::UnsteadyAdvectionDiffusion::SetUpSubSteppingTimeIntegration ( int  intMethod,
const LibUtilities::TimeIntegrationWrapperSharedPtr IntegrationScheme 
)
protected

Definition at line 559 of file UnsteadyAdvectionDiffusion.cpp.

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineOdeRhs(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineProjection(), Nektar::LibUtilities::eBackwardEuler, Nektar::LibUtilities::eBDFImplicitOrder1, Nektar::LibUtilities::eBDFImplicitOrder2, Nektar::LibUtilities::GetTimeIntegrationWrapperFactory(), m_intSteps, m_subStepIntegrationOps, m_subStepIntegrationScheme, SubStepAdvection(), and SubStepProjection().

Referenced by v_InitObject().

562  {
563  // Set to 1 for first step and it will then be increased in
564  // time advance routines
565  switch(intMethod)
566  {
569  {
571 
572  }
573  break;
575  {
577  }
578  break;
579  default:
580  ASSERTL0(0,"Integration method not suitable: Options include BackwardEuler or BDFImplicitOrder1");
581  break;
582  }
583  m_intSteps = IntegrationScheme->GetIntegrationSteps();
584 
585  // set explicit time-integration class operators
588  }
void SubStepProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
BDF multi-step scheme of order 1 (implicit)
LibUtilities::TimeIntegrationWrapperSharedPtr m_subStepIntegrationScheme
void SubStepAdvection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
TimeIntegrationWrapperFactory & GetTimeIntegrationWrapperFactory()
BDF multi-step scheme of order 2 (implicit)
LibUtilities::TimeIntegrationSchemeOperators m_subStepIntegrationOps

◆ SubStepAdvance()

void Nektar::UnsteadyAdvectionDiffusion::SubStepAdvance ( const LibUtilities::TimeIntegrationSolutionSharedPtr integrationSoln,
int  nstep,
NekDouble  time 
)
protected

Definition at line 470 of file UnsteadyAdvectionDiffusion.cpp.

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

Referenced by v_PreIntegrate().

474  {
475  int n;
476  int nsubsteps;
477 
478  NekDouble dt;
479 
480  Array<OneD, Array<OneD, NekDouble> > fields, velfields;
481 
482  static int ncalls = 1;
483  int nint = min(ncalls++, m_intSteps);
484 
485  Array<OneD, NekDouble> CFL(m_fields[0]->GetExpSize(),
487 
488  LibUtilities::CommSharedPtr comm = m_session->GetComm();
489 
490  // Get the proper time step with CFL control
491  dt = GetSubstepTimeStep();
492 
493  nsubsteps = (m_timestep > dt)? ((int)(m_timestep/dt)+1):1;
494  nsubsteps = max(m_minsubsteps, nsubsteps);
495 
496  dt = m_timestep/nsubsteps;
497 
498  if (m_infosteps && !((nstep+1)%m_infosteps) && comm->GetRank() == 0)
499  {
500  cout << "Sub-integrating using "<< nsubsteps
501  << " steps over Dt = " << m_timestep
502  << " (SubStep CFL=" << m_cflSafetyFactor << ")"<< endl;
503  }
504 
505  for (int m = 0; m < nint; ++m)
506  {
507  // We need to update the fields held by the m_integrationSoln
508  fields = integrationSoln->UpdateSolutionVector()[m];
509 
510  // Initialise NS solver which is set up to use a GLM method
511  // with calls to EvaluateAdvection_SetPressureBCs and
512  // SolveUnsteadyStokesSystem
514  SubIntegrationSoln = m_subStepIntegrationScheme->
515  InitializeScheme(dt, fields, time, m_subStepIntegrationOps);
516 
517  for(n = 0; n < nsubsteps; ++n)
518  {
519  fields = m_subStepIntegrationScheme->TimeIntegrate(n, dt, SubIntegrationSoln,
521  }
522 
523  // Reset time integrated solution in m_integrationSoln
524  integrationSoln->SetSolVector(m,fields);
525  }
526  }
LibUtilities::TimeIntegrationWrapperSharedPtr m_subStepIntegrationScheme
NekDouble m_timestep
Time step size.
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:53
double NekDouble
LibUtilities::TimeIntegrationSchemeOperators m_subStepIntegrationOps
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetExpSize()
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
std::shared_ptr< TimeIntegrationSolution > TimeIntegrationSolutionSharedPtr

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

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

Referenced by SetUpSubSteppingTimeIntegration().

597  {
598  int i;
599  int nVariables = inarray.num_elements();
600 
601  /// Get the number of coefficients
602  int ncoeffs = m_fields[0]->GetNcoeffs();
603 
604  /// Define an auxiliary variable to compute the RHS
605  Array<OneD, Array<OneD, NekDouble> > WeakAdv(nVariables);
606  WeakAdv[0] = Array<OneD, NekDouble> (ncoeffs*nVariables);
607  for(i = 1; i < nVariables; ++i)
608  {
609  WeakAdv[i] = WeakAdv[i-1] + ncoeffs;
610  }
611 
612  // Currently assume velocity field is time independent and does not therefore
613  // need extrapolating.
614  // RHS computation using the advection base class
615  m_advObject->Advect(nVariables, m_fields, m_velocity,
616  inarray, outarray, time);
617 
618  for(i = 0; i < nVariables; ++i)
619  {
620  m_fields[i]->IProductWRTBase(outarray[i],WeakAdv[i]);
621  // negation requried due to sign of DoAdvection term to be consistent
622  Vmath::Neg(ncoeffs, WeakAdv[i], 1);
623  }
624 
625  AddAdvectionPenaltyFlux(m_velocity, inarray, WeakAdv);
626 
627 
628  /// Operations to compute the RHS
629  for(i = 0; i < nVariables; ++i)
630  {
631  // Negate the RHS
632  Vmath::Neg(ncoeffs, WeakAdv[i], 1);
633 
634  /// Multiply the flux by the inverse of the mass matrix
635  m_fields[i]->MultiplyByElmtInvMass(WeakAdv[i], WeakAdv[i]);
636 
637  /// Store in outarray the physical values of the RHS
638  m_fields[i]->BwdTrans(WeakAdv[i], outarray[i]);
639  }
640  }
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
Array< OneD, Array< OneD, NekDouble > > m_velocity
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:399
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
void AddAdvectionPenaltyFlux(const Array< OneD, const Array< OneD, NekDouble > > &velfield, const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &Outarray)

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

References ASSERTL1, and Vmath::Vcopy().

Referenced by SetUpSubSteppingTimeIntegration().

649  {
650  boost::ignore_unused(time);
651 
652  ASSERTL1(inarray.num_elements() == outarray.num_elements(),"Inarray and outarray of different sizes ");
653 
654  for(int i = 0; i < inarray.num_elements(); ++i)
655  {
656  Vmath::Vcopy(inarray[i].num_elements(),inarray[i],1,outarray[i],1);
657  }
658  }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064

◆ v_GenerateSummary()

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

Print Summary.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 446 of file UnsteadyAdvectionDiffusion.cpp.

References Nektar::SolverUtils::UnsteadySystem::v_GenerateSummary().

448  {
450  }
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s)
Print a summary of time stepping parameters.

◆ v_InitObject()

void Nektar::UnsteadyAdvectionDiffusion::v_InitObject ( )
protectedvirtual

Initialise the object.

Initialisation object for the unsteady linear advection diffusion equation.

Reimplemented from Nektar::SolverUtils::AdvectionSystem.

Definition at line 63 of file UnsteadyAdvectionDiffusion.cpp.

References ASSERTL0, Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineImplicitSolve(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineOdeRhs(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineProjection(), DoImplicitSolve(), DoOdeProjection(), DoOdeRhs(), Nektar::MultiRegions::eDiscontinuous, Nektar::MultiRegions::eGalerkin, Nektar::MultiRegions::eMixed_CG_Discontinuous, Nektar::SolverUtils::GetAdvectionFactory(), Nektar::SolverUtils::GetDiffusionFactory(), GetFluxVectorAdv(), GetFluxVectorDiff(), Nektar::SolverUtils::EquationSystem::GetFunction(), GetNormalVelocity(), Nektar::SolverUtils::GetRiemannSolverFactory(), Nektar::SolverUtils::AdvectionSystem::m_advObject, m_diffusion, m_epsilon, Nektar::SolverUtils::UnsteadySystem::m_explicitDiffusion, Nektar::SolverUtils::EquationSystem::m_fields, 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_useSpecVanVisc, m_velocity, m_waveFreq, SetUpSubSteppingTimeIntegration(), and Nektar::SolverUtils::AdvectionSystem::v_InitObject().

64  {
66 
67  m_session->LoadParameter("wavefreq", m_waveFreq, 0.0);
68  m_session->LoadParameter("epsilon", m_epsilon, 0.0);
69 
70  // turn on substepping
71  m_session->MatchSolverInfo("Extrapolation", "SubStepping",
72  m_subSteppingScheme, false);
73 
74  // Define Velocity fields
75  m_velocity = Array<OneD, Array<OneD, NekDouble> >(m_spacedim);
76  std::vector<std::string> vel;
77  vel.push_back("Vx");
78  vel.push_back("Vy");
79  vel.push_back("Vz");
80  vel.resize(m_spacedim);
81 
82  GetFunction( "AdvectionVelocity")->Evaluate(vel, m_velocity);
83 
84  m_session->MatchSolverInfo(
85  "SpectralVanishingViscosity", "True", m_useSpecVanVisc, false);
86 
88  {
89  m_session->LoadParameter("SVVCutoffRatio",m_sVVCutoffRatio,0.75);
90  m_session->LoadParameter("SVVDiffCoeff",m_sVVDiffCoeff,0.1);
91  }
92 
93  // Type of advection and diffusion classes to be used
94  switch(m_projectionType)
95  {
96  // Discontinuous field
98  {
99  // Do not forwards transform initial condition
100  m_homoInitialFwd = false;
101 
102  // Advection term
103  string advName;
104  string riemName;
105  m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
107  CreateInstance(advName, advName);
108  m_advObject->SetFluxVector(&UnsteadyAdvectionDiffusion::
109  GetFluxVectorAdv, this);
110  m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
112  CreateInstance(riemName, m_session);
113  m_riemannSolver->SetScalar("Vn", &UnsteadyAdvectionDiffusion::
114  GetNormalVelocity, this);
115  m_advObject->SetRiemannSolver(m_riemannSolver);
116  m_advObject->InitObject (m_session, m_fields);
117 
118  // Diffusion term
120  {
121  std::string diffName;
122  m_session->LoadSolverInfo("DiffusionType", diffName, "LDG");
124  CreateInstance(diffName, diffName);
125  m_diffusion->SetFluxVector(&UnsteadyAdvectionDiffusion::
126  GetFluxVectorDiff, this);
127  m_diffusion->InitObject(m_session, m_fields);
128  }
129 
130  ASSERTL0(m_subSteppingScheme == false,"SubSteppingScheme is not set up for DG projection");
131  break;
132  }
133  // Continuous field
136  {
137  // Advection term
138  std::string advName;
139  m_session->LoadSolverInfo("AdvectionType", advName,
140  "NonConservative");
142  CreateInstance(advName, advName);
143  m_advObject->SetFluxVector(&UnsteadyAdvectionDiffusion::
144  GetFluxVectorAdv, this);
145 
146  if(advName.compare("WeakDG") == 0)
147  {
148  string riemName;
149  m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
151  CreateInstance(riemName, m_session);
152  m_riemannSolver->SetScalar("Vn",
154  GetNormalVelocity, this);
155  m_advObject->SetRiemannSolver(m_riemannSolver);
156  m_advObject->InitObject (m_session, m_fields);
157  }
158 
159  // In case of Galerkin explicit diffusion gives an error
161  {
162  ASSERTL0(false, "Explicit Galerkin diffusion not set up.");
163  }
164  // In case of Galerkin implicit diffusion: do nothing
165  break;
166  }
167  default:
168  {
169  ASSERTL0(false, "Unsupported projection type.");
170  break;
171  }
172  }
173 
179 
180  if(m_subSteppingScheme) // Substepping
181  {
183  "Projection must be set to Mixed_CG_Discontinuous for "
184  "substepping");
186  m_intScheme->GetIntegrationMethod(), m_intScheme);
187  }
188  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
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.
bool m_homoInitialFwd
Flag to determine if simulation should start in homogeneous forward transformed state.
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
bool m_explicitDiffusion
Indicates if explicit or implicit treatment of diffusion is used.
UnsteadyAdvectionDiffusion(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Session reader.
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
Array< OneD, Array< OneD, NekDouble > > m_velocity
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:41
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
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 > & GetNormalVelocity()
Get the normal velocity based on m_velocity.
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Perform the projection.
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
virtual void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the RHS.
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
RiemannSolverFactory & GetRiemannSolverFactory()
int m_spacedim
Spatial dimension (>= expansion dim).
void SetUpSubSteppingTimeIntegration(int intMethod, const LibUtilities::TimeIntegrationWrapperSharedPtr &IntegrationScheme)
SolverUtils::DiffusionSharedPtr m_diffusion
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
Get a SessionFunction by name.
virtual void DoImplicitSolve(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, NekDouble time, NekDouble lambda)
Solve implicitly the diffusion term.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
Wrapper to the time integration scheme.
virtual SOLVER_UTILS_EXPORT void v_InitObject()
Init object for UnsteadySystem class.

◆ v_PreIntegrate()

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

PreIntegration step for substepping.

Perform the extrapolation.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 456 of file UnsteadyAdvectionDiffusion.cpp.

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

457  {
459  {
461  }
462 
463  return false;
464  }
NekDouble m_time
Current time of simulation.
void SubStepAdvance(const LibUtilities::TimeIntegrationSolutionSharedPtr &integrationSoln, int nstep, NekDouble time)
LibUtilities::TimeIntegrationSolutionSharedPtr m_intSoln

Friends And Related Function Documentation

◆ MemoryManager< UnsteadyAdvectionDiffusion >

friend class MemoryManager< UnsteadyAdvectionDiffusion >
friend

Definition at line 48 of file UnsteadyAdvectionDiffusion.h.

Member Data Documentation

◆ className

string Nektar::UnsteadyAdvectionDiffusion::className
static
Initial value:

Name of class.

Definition at line 62 of file UnsteadyAdvectionDiffusion.h.

◆ m_cflSafetyFactor

NekDouble Nektar::UnsteadyAdvectionDiffusion::m_cflSafetyFactor
protected

Definition at line 159 of file UnsteadyAdvectionDiffusion.h.

Referenced by GetSubstepTimeStep(), and SubStepAdvance().

◆ m_diffusion

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

Definition at line 73 of file UnsteadyAdvectionDiffusion.h.

Referenced by DoOdeRhs(), and v_InitObject().

◆ m_epsilon

NekDouble Nektar::UnsteadyAdvectionDiffusion::m_epsilon
private

Definition at line 165 of file UnsteadyAdvectionDiffusion.h.

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

◆ m_infosteps

int Nektar::UnsteadyAdvectionDiffusion::m_infosteps
protected

Definition at line 160 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 161 of file UnsteadyAdvectionDiffusion.h.

Referenced by SubStepAdvance().

◆ m_planeNumber

int Nektar::UnsteadyAdvectionDiffusion::m_planeNumber
protected

Definition at line 79 of file UnsteadyAdvectionDiffusion.h.

Referenced by UnsteadyAdvectionDiffusion().

◆ m_riemannSolver

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

Definition at line 72 of file UnsteadyAdvectionDiffusion.h.

Referenced by v_InitObject().

◆ m_subStepIntegrationOps

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

◆ m_subStepIntegrationScheme

LibUtilities::TimeIntegrationWrapperSharedPtr Nektar::UnsteadyAdvectionDiffusion::m_subStepIntegrationScheme
protected

◆ m_subSteppingScheme

bool Nektar::UnsteadyAdvectionDiffusion::m_subSteppingScheme
protected

Definition at line 68 of file UnsteadyAdvectionDiffusion.h.

Referenced by v_InitObject(), and v_PreIntegrate().

◆ m_sVVCutoffRatio

NekDouble Nektar::UnsteadyAdvectionDiffusion::m_sVVCutoffRatio
protected

Definition at line 70 of file UnsteadyAdvectionDiffusion.h.

Referenced by DoImplicitSolve(), and v_InitObject().

◆ m_sVVDiffCoeff

NekDouble Nektar::UnsteadyAdvectionDiffusion::m_sVVDiffCoeff
protected

Definition at line 71 of file UnsteadyAdvectionDiffusion.h.

Referenced by DoImplicitSolve(), and v_InitObject().

◆ m_traceVn

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

Definition at line 75 of file UnsteadyAdvectionDiffusion.h.

Referenced by GetNormalVel().

◆ m_useSpecVanVisc

bool Nektar::UnsteadyAdvectionDiffusion::m_useSpecVanVisc
protected

Definition at line 69 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 164 of file UnsteadyAdvectionDiffusion.h.

Referenced by v_InitObject().