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

#include <UnsteadyAdvectionDiffusion.h>

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

Public Member Functions

virtual ~UnsteadyAdvectionDiffusion ()
 Destructor. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::AdvectionSystem
SOLVER_UTILS_EXPORT AdvectionSystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 
virtual SOLVER_UTILS_EXPORT ~AdvectionSystem ()
 
AdvectionSharedPtr GetAdvObject ()
 Returns the advection object held by this instance. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
virtual SOLVER_UTILS_EXPORT ~UnsteadySystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT NekDouble GetTimeStep (const Array< OneD, const Array< OneD, NekDouble > > &inarray)
 Calculate the larger time-step mantaining the problem stable. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT void SetUpTraceNormals (void)
 
SOLVER_UTILS_EXPORT void InitObject ()
 Initialises the members of this object. More...
 
SOLVER_UTILS_EXPORT void DoInitialise ()
 Perform any initialisation necessary before solving the problem. More...
 
SOLVER_UTILS_EXPORT void DoSolve ()
 Solve the problem. More...
 
SOLVER_UTILS_EXPORT void TransCoeffToPhys ()
 Transform from coefficient to physical space. More...
 
SOLVER_UTILS_EXPORT void TransPhysToCoeff ()
 Transform from physical to coefficient space. More...
 
SOLVER_UTILS_EXPORT void Output ()
 Perform output operations after solve. More...
 
SOLVER_UTILS_EXPORT NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation. More...
 
SOLVER_UTILS_EXPORT std::string GetSessionName ()
 Get Session name. More...
 
template<class T >
boost::shared_ptr< T > as ()
 
SOLVER_UTILS_EXPORT void ResetSessionName (std::string newname)
 Reset Session name. More...
 
SOLVER_UTILS_EXPORT
LibUtilities::SessionReaderSharedPtr 
GetSession ()
 Get Session name. More...
 
SOLVER_UTILS_EXPORT
MultiRegions::ExpListSharedPtr 
GetPressure ()
 Get pressure field if available. More...
 
SOLVER_UTILS_EXPORT void PrintSummary (std::ostream &out)
 Print a summary of parameters and solver characteristics. More...
 
SOLVER_UTILS_EXPORT void SetLambda (NekDouble lambda)
 Set parameter m_lambda. More...
 
SOLVER_UTILS_EXPORT void EvaluateFunction (Array< OneD, Array< OneD, NekDouble > > &pArray, std::string pFunctionName, const NekDouble pTime=0.0, const int domain=0)
 Evaluates a function as specified in the session file. More...
 
SOLVER_UTILS_EXPORT void EvaluateFunction (std::vector< std::string > pFieldNames, Array< OneD, Array< OneD, NekDouble > > &pFields, const std::string &pName, const NekDouble &pTime=0.0, const int domain=0)
 Populate given fields with the function from session. More...
 
SOLVER_UTILS_EXPORT void EvaluateFunction (std::vector< std::string > pFieldNames, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const std::string &pName, const NekDouble &pTime=0.0, const int domain=0)
 Populate given fields with the function from session. More...
 
SOLVER_UTILS_EXPORT void EvaluateFunction (std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, const NekDouble &pTime=0.0, const int domain=0)
 
SOLVER_UTILS_EXPORT std::string DescribeFunction (std::string pFieldName, const std::string &pFunctionName, const int domain)
 Provide a description of a function for a given field name. More...
 
SOLVER_UTILS_EXPORT void InitialiseBaseFlow (Array< OneD, Array< OneD, NekDouble > > &base)
 Perform initialisation of the base flow. More...
 
SOLVER_UTILS_EXPORT void SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 Initialise the data in the dependent fields. More...
 
SOLVER_UTILS_EXPORT void EvaluateExactSolution (int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 Evaluates an exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln, bool Normalised=false)
 Compute the L2 error between fields and a given exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, bool Normalised=false)
 Compute the L2 error of the fields. More...
 
SOLVER_UTILS_EXPORT Array
< OneD, NekDouble
ErrorExtraPoints (unsigned int field)
 Compute error (L2 and L_inf) over an larger set of quadrature points return [L2 Linf]. More...
 
SOLVER_UTILS_EXPORT void WeakAdvectionGreensDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\nabla \phi \cdot F) $. More...
 
SOLVER_UTILS_EXPORT void WeakAdvectionDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\phi, \nabla \cdot F) $. More...
 
SOLVER_UTILS_EXPORT void WeakAdvectionNonConservativeForm (const Array< OneD, Array< OneD, NekDouble > > &V, const Array< OneD, const NekDouble > &u, Array< OneD, NekDouble > &outarray, bool UseContCoeffs=false)
 Compute the inner product $ (\phi, V\cdot \nabla u) $. More...
 
f SOLVER_UTILS_EXPORT void AdvectionNonConservativeForm (const Array< OneD, Array< OneD, NekDouble > > &V, const Array< OneD, const NekDouble > &u, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wk=NullNekDouble1DArray)
 Compute the non-conservative advection. More...
 
SOLVER_UTILS_EXPORT void WeakDGAdvection (const Array< OneD, Array< OneD, NekDouble > > &InField, Array< OneD, Array< OneD, NekDouble > > &OutField, bool NumericalFluxIncludesNormal=true, bool InFieldIsInPhysSpace=false, int nvariables=0)
 Calculate the weak discontinuous Galerkin advection. More...
 
SOLVER_UTILS_EXPORT void WeakDGDiffusion (const Array< OneD, Array< OneD, NekDouble > > &InField, Array< OneD, Array< OneD, NekDouble > > &OutField, bool NumericalFluxIncludesNormal=true, bool InFieldIsInPhysSpace=false)
 Calculate weak DG Diffusion in the LDG form. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n)
 Write checkpoint file of m_fields. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 Write checkpoint file of custom data fields. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_BaseFlow (const int n)
 Write base flow file of m_fields. More...
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname)
 Write field data to the given filename. More...
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 Write input fields to the given filename. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 Input field data from the given file. More...
 
SOLVER_UTILS_EXPORT void ImportFldToMultiDomains (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int ndomains)
 Input field data from the given file to multiple domains. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, std::vector< std::string > &fieldStr, Array< OneD, Array< OneD, NekDouble > > &coeffs)
 Output a field. Input field data into array from the given file. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, MultiRegions::ExpListSharedPtr &pField, std::string &pFieldName)
 Output a field. Input field data into ExpList from the given file. More...
 
SOLVER_UTILS_EXPORT void ScanForHistoryPoints ()
 Builds map of which element holds each history point. More...
 
SOLVER_UTILS_EXPORT void WriteHistoryData (std::ostream &out)
 Probe each history point and write to file. More...
 
SOLVER_UTILS_EXPORT void SessionSummary (SummaryList &vSummary)
 Write out a session summary. More...
 
SOLVER_UTILS_EXPORT Array
< OneD,
MultiRegions::ExpListSharedPtr > & 
UpdateFields ()
 
SOLVER_UTILS_EXPORT
LibUtilities::FieldMetaDataMap
UpdateFieldMetaDataMap ()
 Get hold of FieldInfoMap so it can be updated. More...
 
SOLVER_UTILS_EXPORT NekDouble GetFinalTime ()
 Return final time. More...
 
SOLVER_UTILS_EXPORT int GetNcoeffs ()
 
SOLVER_UTILS_EXPORT int GetNcoeffs (const int eid)
 
SOLVER_UTILS_EXPORT int GetNumExpModes ()
 
SOLVER_UTILS_EXPORT const
Array< OneD, int > 
GetNumExpModesPerExp ()
 
SOLVER_UTILS_EXPORT int GetNvariables ()
 
SOLVER_UTILS_EXPORT const
std::string 
GetVariable (unsigned int i)
 
SOLVER_UTILS_EXPORT int GetTraceTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTraceNpoints ()
 
SOLVER_UTILS_EXPORT int GetExpSize ()
 
SOLVER_UTILS_EXPORT int GetPhys_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetCoeff_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTotPoints (int n)
 
SOLVER_UTILS_EXPORT int GetNpoints ()
 
SOLVER_UTILS_EXPORT int GetNumElmVelocity ()
 
SOLVER_UTILS_EXPORT int GetSteps ()
 
SOLVER_UTILS_EXPORT NekDouble GetTimeStep ()
 
SOLVER_UTILS_EXPORT void CopyFromPhysField (const int i, Array< OneD, NekDouble > &output)
 
SOLVER_UTILS_EXPORT void CopyToPhysField (const int i, Array< OneD, NekDouble > &output)
 
SOLVER_UTILS_EXPORT void SetSteps (const int steps)
 
SOLVER_UTILS_EXPORT void ZeroPhysFields ()
 
SOLVER_UTILS_EXPORT void FwdTransFields ()
 
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
 
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &fluxX, Array< OneD, Array< OneD, NekDouble > > &fluxY)
 
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, const int j, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
 
SOLVER_UTILS_EXPORT void NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numflux)
 
SOLVER_UTILS_EXPORT void NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxX, Array< OneD, Array< OneD, NekDouble > > &numfluxY)
 
SOLVER_UTILS_EXPORT void NumFluxforScalar (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &uflux)
 
SOLVER_UTILS_EXPORT void NumFluxforVector (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
 
SOLVER_UTILS_EXPORT void SetModifiedBasis (const bool modbasis)
 
SOLVER_UTILS_EXPORT int NoCaseStringCompare (const string &s1, const string &s2)
 Perform a case-insensitive string comparison. More...
 
SOLVER_UTILS_EXPORT int GetCheckpointNumber ()
 
SOLVER_UTILS_EXPORT void SetCheckpointNumber (int num)
 
SOLVER_UTILS_EXPORT int GetCheckpointSteps ()
 
SOLVER_UTILS_EXPORT void SetCheckpointSteps (int num)
 
SOLVER_UTILS_EXPORT void SetTime (const NekDouble time)
 
SOLVER_UTILS_EXPORT void SetInitialStep (const int step)
 
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time. More...
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp ()
 Virtual function to identify if operator is negated in DoSolve. More...
 

Static Public Member Functions

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

Static Public Attributes

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

Protected Member Functions

 UnsteadyAdvectionDiffusion (const LibUtilities::SessionReaderSharedPtr &pSession)
 Session reader. More...
 
void GetFluxVectorAdv (const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
 Evaluate the flux at each solution point for the advection part. More...
 
void GetFluxVectorDiff (const int i, const int j, const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &derivatives, Array< OneD, Array< OneD, NekDouble > > &flux)
 Evaluate the flux at each solution point for the diffusion part. More...
 
virtual void DoOdeRhs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 Compute the RHS. More...
 
void DoOdeProjection (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 Perform the projection. More...
 
virtual void DoImplicitSolve (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, NekDouble time, NekDouble lambda)
 Solve implicitly the diffusion term. More...
 
Array< OneD, NekDouble > & GetNormalVelocity ()
 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::UnsteadySystem
SOLVER_UTILS_EXPORT UnsteadySystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 Initialises UnsteadySystem class members. More...
 
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve ()
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise ()
 Sets up initial conditions. More...
 
virtual SOLVER_UTILS_EXPORT void v_AppendOutput1D (Array< OneD, Array< OneD, NekDouble > > &solution1D)
 Print the solution at each solution point in a txt file. More...
 
virtual SOLVER_UTILS_EXPORT void v_NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numflux)
 
virtual SOLVER_UTILS_EXPORT void v_NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxX, Array< OneD, Array< OneD, NekDouble > > &numfluxY)
 
virtual SOLVER_UTILS_EXPORT void v_NumFluxforScalar (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &uflux)
 
virtual SOLVER_UTILS_EXPORT void v_NumFluxforVector (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
 
virtual SOLVER_UTILS_EXPORT
NekDouble 
v_GetTimeStep (const Array< OneD, const Array< OneD, NekDouble > > &inarray)
 Return the timestep to be used for the next step in the time-marching loop. More...
 
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT bool v_SteadyStateCheck (int step)
 
SOLVER_UTILS_EXPORT void CheckForRestartTime (NekDouble &time)
 
SOLVER_UTILS_EXPORT void SVVVarDiffCoeff (const Array< OneD, Array< OneD, NekDouble > > vel, StdRegions::VarCoeffMap &varCoeffMap)
 Evaluate the SVV diffusion coefficient according to Moura's paper where it should proportional to h time velocity. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 Initialises EquationSystem class members. More...
 
int nocase_cmp (const string &s1, const string &s2)
 
virtual SOLVER_UTILS_EXPORT
NekDouble 
v_LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Virtual function for the L_inf error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT
NekDouble 
v_L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false)
 Virtual function for the L_2 error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransCoeffToPhys ()
 Virtual function for transformation to physical space. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransPhysToCoeff ()
 Virtual function for transformation to coefficient space. More...
 
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 
virtual SOLVER_UTILS_EXPORT void v_EvaluateExactSolution (unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 
SOLVER_UTILS_EXPORT void SetUpBaseFields (SpatialDomains::MeshGraphSharedPtr &mesh)
 
SOLVER_UTILS_EXPORT void ImportFldBase (std::string pInfile, SpatialDomains::MeshGraphSharedPtr pGraph)
 
virtual SOLVER_UTILS_EXPORT void v_Output (void)
 
virtual SOLVER_UTILS_EXPORT
MultiRegions::ExpListSharedPtr 
v_GetPressure (void)
 
virtual SOLVER_UTILS_EXPORT void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 

Protected Attributes

bool m_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...
 
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
LibUtilities::TimeIntegrationSolutionSharedPtr m_intSoln
 
NekDouble m_epsilon
 
bool m_explicitDiffusion
 Indicates if explicit or implicit treatment of diffusion is used. More...
 
bool m_explicitAdvection
 Indicates if explicit or implicit treatment of advection is used. More...
 
bool m_explicitReaction
 Indicates if explicit or implicit treatment of reaction is used. More...
 
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state. More...
 
std::vector< int > m_intVariables
 
std::vector< FilterSharedPtrm_filters
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output. More...
 
map< std::string, Array< OneD,
Array< OneD, float > > > 
m_interpWeights
 Map of the interpolation weights for a specific filename. More...
 
map< std::string, Array< OneD,
Array< OneD, unsigned int > > > 
m_interpInds
 Map of the interpolation indices for a specific filename. More...
 
Array< OneD,
MultiRegions::ExpListSharedPtr
m_fields
 Array holding all dependent variables. More...
 
Array< OneD,
MultiRegions::ExpListSharedPtr
m_base
 Base fields. More...
 
Array< OneD,
MultiRegions::ExpListSharedPtr
m_derivedfields
 Array holding all dependent variables. More...
 
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
 Pointer to boundary conditions object. More...
 
SpatialDomains::MeshGraphSharedPtr m_graph
 Pointer to graph defining mesh. More...
 
std::string m_sessionName
 Name of the session. More...
 
NekDouble m_time
 Current time of simulation. More...
 
int m_initialStep
 Number of the step where the simulation should begin. More...
 
NekDouble m_fintime
 Finish time of the simulation. More...
 
NekDouble m_timestep
 Time step size. More...
 
NekDouble m_lambda
 Lambda constant in real system if one required. More...
 
std::set< std::string > m_loadedFields
 
NekDouble m_checktime
 Time between checkpoints. More...
 
int m_nchk
 Number of checkpoints written so far. More...
 
int m_steps
 Number of steps to take. More...
 
int m_checksteps
 Number of steps between checkpoints. More...
 
int m_spacedim
 Spatial dimension (>= expansion dim). More...
 
int m_expdim
 Expansion dimension. More...
 
bool m_singleMode
 Flag to determine if single homogeneous mode is used. More...
 
bool m_halfMode
 Flag to determine if half homogeneous mode is used. More...
 
bool m_multipleModes
 Flag to determine if use multiple homogenenous modes are used. More...
 
bool m_useFFT
 Flag to determine if FFT is used for homogeneous transform. More...
 
bool m_homogen_dealiasing
 Flag to determine if dealiasing is used for homogeneous simulations. More...
 
bool m_specHP_dealiasing
 Flag to determine if dealisising is usde for the Spectral/hp element discretisation. More...
 
enum MultiRegions::ProjectionType m_projectionType
 Type of projection; e.g continuous or discontinuous. More...
 
Array< OneD, Array< OneD,
NekDouble > > 
m_traceNormals
 Array holding trace normals for DG simulations in the forwards direction. More...
 
Array< OneD, Array< OneD,
Array< OneD, NekDouble > > > 
m_gradtan
 1 x nvariable x nq More...
 
Array< OneD, Array< OneD,
Array< OneD, NekDouble > > > 
m_tanbasis
 2 x m_spacedim x nq More...
 
Array< OneD, bool > m_checkIfSystemSingular
 Flag to indicate if the fields should be checked for singularity. More...
 
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
 Map to identify relevant solver info to dump in output fields. More...
 
int m_NumQuadPointsError
 Number of Quadrature points used to work out the error. More...
 
enum HomogeneousType m_HomogeneousType
 
NekDouble m_LhomX
 physical length in X direction (if homogeneous) More...
 
NekDouble m_LhomY
 physical length in Y direction (if homogeneous) More...
 
NekDouble m_LhomZ
 physical length in Z direction (if homogeneous) More...
 
int m_npointsX
 number of points in X direction (if homogeneous) More...
 
int m_npointsY
 number of points in Y direction (if homogeneous) More...
 
int m_npointsZ
 number of points in Z direction (if homogeneous) More...
 
int m_HomoDirec
 number of homogenous directions More...
 
int m_NumMode
 Mode to use in case of single mode analysis. More...
 

Private Attributes

NekDouble m_waveFreq
 
NekDouble m_epsilon
 

Friends

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

Detailed Description

Definition at line 46 of file UnsteadyAdvectionDiffusion.h.

Constructor & Destructor Documentation

Nektar::UnsteadyAdvectionDiffusion::~UnsteadyAdvectionDiffusion ( )
virtual

Destructor.

Unsteady linear advection diffusion equation destructor.

Definition at line 194 of file UnsteadyAdvectionDiffusion.cpp.

195  {
196  }
Nektar::UnsteadyAdvectionDiffusion::UnsteadyAdvectionDiffusion ( const LibUtilities::SessionReaderSharedPtr pSession)
protected

Session reader.

Definition at line 47 of file UnsteadyAdvectionDiffusion.cpp.

References m_planeNumber.

49  : UnsteadySystem(pSession),
50  AdvectionSystem(pSession)
51  {
52  m_planeNumber = 0;
53  }
SOLVER_UTILS_EXPORT AdvectionSystem(const LibUtilities::SessionReaderSharedPtr &pSession)
SOLVER_UTILS_EXPORT UnsteadySystem(const LibUtilities::SessionReaderSharedPtr &pSession)
Initialises UnsteadySystem class members.

Member Function Documentation

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

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

Referenced by SubStepAdvection().

657  {
658  ASSERTL1(physfield.num_elements() == Outarray.num_elements(),
659  "Physfield and outarray are of different dimensions");
660 
661  int i;
662 
663  /// Number of trace points
664  int nTracePts = m_fields[0]->GetTrace()->GetNpoints();
665 
666  /// Forward state array
667  Array<OneD, NekDouble> Fwd(3*nTracePts);
668 
669  /// Backward state array
670  Array<OneD, NekDouble> Bwd = Fwd + nTracePts;
671 
672  /// upwind numerical flux state array
673  Array<OneD, NekDouble> numflux = Bwd + nTracePts;
674 
675  /// Normal velocity array
676  Array<OneD, NekDouble> Vn = GetNormalVel(velfield);
677 
678  for(i = 0; i < physfield.num_elements(); ++i)
679  {
680  /// Extract forwards/backwards trace spaces
681  /// Note: Needs to have correct i value to get boundary conditions
682  m_fields[i]->GetFwdBwdTracePhys(physfield[i], Fwd, Bwd);
683 
684  /// Upwind between elements
685  m_fields[0]->GetTrace()->Upwind(Vn, Fwd, Bwd, numflux);
686 
687  /// Construct difference between numflux and Fwd,Bwd
688  Vmath::Vsub(nTracePts, numflux, 1, Fwd, 1, Fwd, 1);
689  Vmath::Vsub(nTracePts, numflux, 1, Bwd, 1, Bwd, 1);
690 
691  /// Calculate the numerical fluxes multipling Fwd, Bwd and
692  /// numflux by the normal advection velocity
693  Vmath::Vmul(nTracePts, Fwd, 1, Vn, 1, Fwd, 1);
694  Vmath::Vmul(nTracePts, Bwd, 1, Vn, 1, Bwd, 1);
695 
696  m_fields[0]->AddFwdBwdTraceIntegral(Fwd,Bwd,Outarray[i]);
697  }
698  }
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:329
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:191
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:169
static SolverUtils::EquationSystemSharedPtr Nektar::UnsteadyAdvectionDiffusion::create ( const LibUtilities::SessionReaderSharedPtr pSession)
inlinestatic

Creates an instance of this class.

Definition at line 52 of file UnsteadyAdvectionDiffusion.h.

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

53  {
56  AllocateSharedPtr(pSession);
57  p->InitObject();
58  return p;
59  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
boost::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.
void Nektar::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 343 of file UnsteadyAdvectionDiffusion.cpp.

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

Referenced by v_InitObject().

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

298  {
299  int i;
300  int nvariables = inarray.num_elements();
301  SetBoundaryConditions(time);
302  switch(m_projectionType)
303  {
305  {
306  // Just copy over array
307  int npoints = GetNpoints();
308 
309  for(i = 0; i < nvariables; ++i)
310  {
311  Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
312  }
313  break;
314  }
317  {
319 
320  for(i = 0; i < nvariables; ++i)
321  {
322  m_fields[i]->FwdTrans(inarray[i], coeffs);
323  m_fields[i]->BwdTrans_IterPerExp(coeffs, outarray[i]);
324  }
325  break;
326  }
327  default:
328  {
329  ASSERTL0(false, "Unknown projection scheme");
330  break;
331  }
332  }
333  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
SOLVER_UTILS_EXPORT int GetNpoints()
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:1047
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 244 of file UnsteadyAdvectionDiffusion.cpp.

References Nektar::MultiRegions::eDiscontinuous, Nektar::SolverUtils::EquationSystem::GetNpoints(), Nektar::SolverUtils::AdvectionSystem::m_advObject, m_diffusion, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_projectionType, m_velocity, Vmath::Neg(), and Vmath::Vadd().

Referenced by v_InitObject().

248  {
249  // Number of fields (variables of the problem)
250  int nVariables = inarray.num_elements();
251 
252  // Number of solution points
253  int nSolutionPts = GetNpoints();
254 
255  Array<OneD, Array<OneD, NekDouble> > outarrayDiff(nVariables);
256 
257  for (int i = 0; i < nVariables; ++i)
258  {
259  outarrayDiff[i] = Array<OneD, NekDouble>(nSolutionPts, 0.0);
260  }
261 
262  // RHS computation using the new advection base class
263  m_advObject->Advect(nVariables, m_fields, m_velocity,
264  inarray, outarray, time);
265 
266  // Negate the RHS
267  for (int i = 0; i < nVariables; ++i)
268  {
269  Vmath::Neg(nSolutionPts, outarray[i], 1);
270  }
271 
272  // No explicit diffusion for CG
274  {
275  m_diffusion->Diffuse(nVariables, m_fields, inarray, outarrayDiff);
276 
277  for (int i = 0; i < nVariables; ++i)
278  {
279  Vmath::Vadd(nSolutionPts, &outarray[i][0], 1,
280  &outarrayDiff[i][0], 1, &outarray[i][0], 1);
281  }
282  }
283 
284  }
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
Array< OneD, Array< OneD, NekDouble > > m_velocity
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SolverUtils::DiffusionSharedPtr m_diffusion
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:382
SOLVER_UTILS_EXPORT int GetNpoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:285
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 399 of file UnsteadyAdvectionDiffusion.cpp.

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

Referenced by v_InitObject().

402  {
403  ASSERTL1(flux[0].num_elements() == m_velocity.num_elements(),
404  "Dimension of flux array and velocity array do not match");
405 
406  const int nq = m_fields[0]->GetNpoints();
407 
408  for (int i = 0; i < flux.num_elements(); ++i)
409  {
410  for (int j = 0; j < flux[0].num_elements(); ++j)
411  {
412  Vmath::Vmul(nq, physfield[i], 1, m_velocity[j], 1,
413  flux[i][j], 1);
414  }
415  }
416  }
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:191
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:169
void Nektar::UnsteadyAdvectionDiffusion::GetFluxVectorDiff ( const int  i,
const int  j,
const Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, NekDouble > > &  derivatives,
Array< OneD, Array< OneD, NekDouble > > &  flux 
)
protected

Evaluate the flux at each solution point for the diffusion part.

Return the flux vector for the diffusion part.

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

Definition at line 427 of file UnsteadyAdvectionDiffusion.cpp.

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

Referenced by v_InitObject().

433  {
434  for (int k = 0; k < flux.num_elements(); ++k)
435  {
436  Vmath::Zero(GetNpoints(), flux[k], 1);
437  }
438  Vmath::Vcopy(GetNpoints(), physfield[i], 1, flux[j], 1);
439  }
SOLVER_UTILS_EXPORT int GetNpoints()
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:359
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
Array< OneD, NekDouble > Nektar::UnsteadyAdvectionDiffusion::GetMaxStdVelocity ( const Array< OneD, Array< OneD, NekDouble > >  inarray)
protected

Definition at line 701 of file UnsteadyAdvectionDiffusion.cpp.

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

Referenced by GetSubstepTimeStep().

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

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

References GetNormalVel(), and m_velocity.

Referenced by v_InitObject().

203  {
204  return GetNormalVel(m_velocity);
205  }
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.
NekDouble Nektar::UnsteadyAdvectionDiffusion::GetSubstepTimeStep ( )
protected

Definition at line 527 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().

528  {
529  int n_element = m_fields[0]->GetExpSize();
530 
531  const Array<OneD, int> ExpOrder=m_fields[0]->EvalBasisNumModesMaxPerExp();
532  Array<OneD, int> ExpOrderList (n_element, ExpOrder);
533 
534  const NekDouble cLambda = 0.2; // Spencer book pag. 317
535 
536  Array<OneD, NekDouble> tstep (n_element, 0.0);
537  Array<OneD, NekDouble> stdVelocity(n_element, 0.0);
538 
539  stdVelocity = GetMaxStdVelocity(m_velocity);
540 
541  for(int el = 0; el < n_element; ++el)
542  {
543  tstep[el] = m_cflSafetyFactor /
544  (stdVelocity[el] * cLambda *
545  (ExpOrder[el]-1) * (ExpOrder[el]-1));
546  }
547 
548  NekDouble TimeStep = Vmath::Vmin(n_element, tstep, 1);
549  m_session->GetComm()->AllReduce(TimeStep,LibUtilities::ReduceMin);
550 
551  return TimeStep;
552  }
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:857
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.
void Nektar::UnsteadyAdvectionDiffusion::SetUpSubSteppingTimeIntegration ( int  intMethod,
const LibUtilities::TimeIntegrationWrapperSharedPtr IntegrationScheme 
)
protected

Definition at line 554 of file UnsteadyAdvectionDiffusion.cpp.

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, >::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().

557  {
558  // Set to 1 for first step and it will then be increased in
559  // time advance routines
560  switch(intMethod)
561  {
564  {
566 
567  }
568  break;
570  {
572  }
573  break;
574  default:
575  ASSERTL0(0,"Integration method not suitable: Options include BackwardEuler or BDFImplicitOrder1");
576  break;
577  }
578  m_intSteps = IntegrationScheme->GetIntegrationSteps();
579 
580  // set explicit time-integration class operators
583  }
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:161
BDF multi-step scheme of order 1 (implicit)
tBaseSharedPtr CreateInstance(tKey idKey BOOST_PP_COMMA_IF(MAX_PARAM) BOOST_PP_ENUM_BINARY_PARAMS(MAX_PARAM, tParam, x))
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:162
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)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
TimeIntegrationWrapperFactory & GetTimeIntegrationWrapperFactory()
BDF multi-step scheme of order 2 (implicit)
LibUtilities::TimeIntegrationSchemeOperators m_subStepIntegrationOps
void Nektar::UnsteadyAdvectionDiffusion::SubStepAdvance ( const LibUtilities::TimeIntegrationSolutionSharedPtr integrationSoln,
int  nstep,
NekDouble  time 
)
protected

Definition at line 465 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().

469  {
470  int n;
471  int nsubsteps;
472 
473  NekDouble dt;
474 
475  Array<OneD, Array<OneD, NekDouble> > fields, velfields;
476 
477  static int ncalls = 1;
478  int nint = min(ncalls++, m_intSteps);
479 
482 
483  LibUtilities::CommSharedPtr comm = m_session->GetComm();
484 
485  // Get the proper time step with CFL control
486  dt = GetSubstepTimeStep();
487 
488  nsubsteps = (m_timestep > dt)? ((int)(m_timestep/dt)+1):1;
489  nsubsteps = max(m_minsubsteps, nsubsteps);
490 
491  dt = m_timestep/nsubsteps;
492 
493  if (m_infosteps && !((nstep+1)%m_infosteps) && comm->GetRank() == 0)
494  {
495  cout << "Sub-integrating using "<< nsubsteps
496  << " steps over Dt = " << m_timestep
497  << " (SubStep CFL=" << m_cflSafetyFactor << ")"<< endl;
498  }
499 
500  for (int m = 0; m < nint; ++m)
501  {
502  // We need to update the fields held by the m_integrationSoln
503  fields = integrationSoln->UpdateSolutionVector()[m];
504 
505  // Initialise NS solver which is set up to use a GLM method
506  // with calls to EvaluateAdvection_SetPressureBCs and
507  // SolveUnsteadyStokesSystem
509  SubIntegrationSoln = m_subStepIntegrationScheme->
510  InitializeScheme(dt, fields, time, m_subStepIntegrationOps);
511 
512  for(n = 0; n < nsubsteps; ++n)
513  {
514  fields = m_subStepIntegrationScheme->TimeIntegrate(n, dt, SubIntegrationSoln,
516  }
517 
518  // Reset time integrated solution in m_integrationSoln
519  integrationSoln->SetSolVector(m,fields);
520  }
521  }
LibUtilities::TimeIntegrationWrapperSharedPtr m_subStepIntegrationScheme
NekDouble m_timestep
Time step size.
boost::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.
boost::shared_ptr< TimeIntegrationSolution > TimeIntegrationSolutionSharedPtr
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 588 of file UnsteadyAdvectionDiffusion.cpp.

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

Referenced by SetUpSubSteppingTimeIntegration().

592  {
593  int i;
594  int nVariables = inarray.num_elements();
595 
596  /// Get the number of coefficients
597  int ncoeffs = m_fields[0]->GetNcoeffs();
598 
599  /// Define an auxiliary variable to compute the RHS
600  Array<OneD, Array<OneD, NekDouble> > WeakAdv(nVariables);
601  WeakAdv[0] = Array<OneD, NekDouble> (ncoeffs*nVariables);
602  for(i = 1; i < nVariables; ++i)
603  {
604  WeakAdv[i] = WeakAdv[i-1] + ncoeffs;
605  }
606 
607  // Currently assume velocity field is time independent and does not therefore
608  // need extrapolating.
609  // RHS computation using the advection base class
610  m_advObject->Advect(nVariables, m_fields, m_velocity,
611  inarray, outarray, time);
612 
613  for(i = 0; i < nVariables; ++i)
614  {
615  m_fields[i]->IProductWRTBase(outarray[i],WeakAdv[i]);
616  // negation requried due to sign of DoAdvection term to be consistent
617  Vmath::Neg(ncoeffs, WeakAdv[i], 1);
618  }
619 
620  AddAdvectionPenaltyFlux(m_velocity, inarray, WeakAdv);
621 
622 
623  /// Operations to compute the RHS
624  for(i = 0; i < nVariables; ++i)
625  {
626  // Negate the RHS
627  Vmath::Neg(ncoeffs, WeakAdv[i], 1);
628 
629  /// Multiply the flux by the inverse of the mass matrix
630  m_fields[i]->MultiplyByElmtInvMass(WeakAdv[i], WeakAdv[i]);
631 
632  /// Store in outarray the physical values of the RHS
633  m_fields[i]->BwdTrans(WeakAdv[i], outarray[i]);
634  }
635  }
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:382
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)
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 640 of file UnsteadyAdvectionDiffusion.cpp.

References ASSERTL1, and Vmath::Vcopy().

Referenced by SetUpSubSteppingTimeIntegration().

644  {
645  ASSERTL1(inarray.num_elements() == outarray.num_elements(),"Inarray and outarray of different sizes ");
646 
647  for(int i = 0; i < inarray.num_elements(); ++i)
648  {
649  Vmath::Vcopy(inarray[i].num_elements(),inarray[i],1,outarray[i],1);
650  }
651  }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
void Nektar::UnsteadyAdvectionDiffusion::v_GenerateSummary ( SolverUtils::SummaryList s)
protectedvirtual

Print Summary.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 441 of file UnsteadyAdvectionDiffusion.cpp.

443  {
444  AdvectionSystem::v_GenerateSummary(s);
445  }
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 59 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::EquationSystem::EvaluateFunction(), Nektar::SolverUtils::GetAdvectionFactory(), Nektar::SolverUtils::GetDiffusionFactory(), GetFluxVectorAdv(), GetFluxVectorDiff(), 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, and SetUpSubSteppingTimeIntegration().

60  {
61  AdvectionSystem::v_InitObject();
62 
63  m_session->LoadParameter("wavefreq", m_waveFreq, 0.0);
64  m_session->LoadParameter("epsilon", m_epsilon, 0.0);
65 
66  // turn on substepping
67  m_session->MatchSolverInfo("Extrapolation", "SubStepping",
68  m_subSteppingScheme, false);
69 
70  // Define Velocity fields
72  std::vector<std::string> vel;
73  vel.push_back("Vx");
74  vel.push_back("Vy");
75  vel.push_back("Vz");
76  vel.resize(m_spacedim);
77 
78  EvaluateFunction(vel, m_velocity, "AdvectionVelocity");
79 
80  m_session->MatchSolverInfo(
81  "SpectralVanishingViscosity", "True", m_useSpecVanVisc, false);
82 
84  {
85  m_session->LoadParameter("SVVCutoffRatio",m_sVVCutoffRatio,0.75);
86  m_session->LoadParameter("SVVDiffCoeff",m_sVVDiffCoeff,0.1);
87  }
88 
89  // Type of advection and diffusion classes to be used
90  switch(m_projectionType)
91  {
92  // Discontinuous field
94  {
95  // Do not forwards transform initial condition
96  m_homoInitialFwd = false;
97 
98  // Advection term
99  string advName;
100  string riemName;
101  m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
103  CreateInstance(advName, advName);
104  m_advObject->SetFluxVector(&UnsteadyAdvectionDiffusion::
105  GetFluxVectorAdv, this);
106  m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
108  CreateInstance(riemName);
109  m_riemannSolver->SetScalar("Vn", &UnsteadyAdvectionDiffusion::
110  GetNormalVelocity, this);
111  m_advObject->SetRiemannSolver(m_riemannSolver);
112  m_advObject->InitObject (m_session, m_fields);
113 
114  // Diffusion term
115  std::string diffName;
116  m_session->LoadSolverInfo("DiffusionType", diffName, "LDG");
118  CreateInstance(diffName, diffName);
119  m_diffusion->SetFluxVector(&UnsteadyAdvectionDiffusion::
120  GetFluxVectorDiff, this);
121  m_diffusion->InitObject(m_session, m_fields);
122 
123  ASSERTL0(m_subSteppingScheme == false,"SubSteppingScheme is not set up for DG projection");
124  break;
125  }
126  // Continuous field
129  {
130  // Advection term
131  std::string advName;
132  m_session->LoadSolverInfo("AdvectionType", advName,
133  "NonConservative");
135  CreateInstance(advName, advName);
136  m_advObject->SetFluxVector(&UnsteadyAdvectionDiffusion::
137  GetFluxVectorAdv, this);
138 
139  if(advName.compare("WeakDG") == 0)
140  {
141  string riemName;
142  m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
144  CreateInstance(riemName);
145  m_riemannSolver->SetScalar("Vn",
147  GetNormalVelocity, this);
148  m_advObject->SetRiemannSolver(m_riemannSolver);
149  m_advObject->InitObject (m_session, m_fields);
150  }
151 
152  // In case of Galerkin explicit diffusion gives an error
154  {
155  ASSERTL0(false, "Explicit Galerkin diffusion not set up.");
156  }
157  // In case of Galerkin implicit diffusion: do nothing
158  break;
159  }
160  default:
161  {
162  ASSERTL0(false, "Unsupported projection type.");
163  break;
164  }
165  }
166 
169 
170  if(m_subSteppingScheme) // Substepping
171  {
173  "Projection must be set to Mixed_CG_Discontinuous for "
174  "substepping");
176  m_intScheme->GetIntegrationMethod(), m_intScheme);
177 
178  }
179  else // Standard velocity correction scheme
180  {
182  }
183 
185  m_explicitDiffusion == 1)
186  {
188  }
189  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
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.
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:42
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
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:46
SOLVER_UTILS_EXPORT void EvaluateFunction(Array< OneD, Array< OneD, NekDouble > > &pArray, std::string pFunctionName, const NekDouble pTime=0.0, const int domain=0)
Evaluates a function as specified in the session file.
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.
void GetFluxVectorDiff(const int i, const int j, const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &derivatives, Array< OneD, Array< OneD, NekDouble > > &flux)
Evaluate the flux at each solution point for the diffusion part.
bool Nektar::UnsteadyAdvectionDiffusion::v_PreIntegrate ( int  step)
protectedvirtual

PreIntegration step for substepping.

Perform the extrapolation.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 451 of file UnsteadyAdvectionDiffusion.cpp.

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

452  {
454  {
456  }
457 
458  return false;
459  }
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

friend class MemoryManager< UnsteadyAdvectionDiffusion >
friend

Definition at line 49 of file UnsteadyAdvectionDiffusion.h.

Member Data Documentation

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

Name of class.

Definition at line 61 of file UnsteadyAdvectionDiffusion.h.

NekDouble Nektar::UnsteadyAdvectionDiffusion::m_cflSafetyFactor
protected

Definition at line 159 of file UnsteadyAdvectionDiffusion.h.

Referenced by GetSubstepTimeStep(), and SubStepAdvance().

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

Definition at line 72 of file UnsteadyAdvectionDiffusion.h.

Referenced by DoOdeRhs(), and v_InitObject().

NekDouble Nektar::UnsteadyAdvectionDiffusion::m_epsilon
private

Definition at line 165 of file UnsteadyAdvectionDiffusion.h.

Referenced by DoImplicitSolve(), and v_InitObject().

int Nektar::UnsteadyAdvectionDiffusion::m_infosteps
protected

Definition at line 160 of file UnsteadyAdvectionDiffusion.h.

Referenced by SubStepAdvance().

int Nektar::UnsteadyAdvectionDiffusion::m_intSteps
protected
int Nektar::UnsteadyAdvectionDiffusion::m_minsubsteps
protected

Definition at line 161 of file UnsteadyAdvectionDiffusion.h.

Referenced by SubStepAdvance().

int Nektar::UnsteadyAdvectionDiffusion::m_planeNumber
protected

Definition at line 78 of file UnsteadyAdvectionDiffusion.h.

Referenced by UnsteadyAdvectionDiffusion().

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

Definition at line 71 of file UnsteadyAdvectionDiffusion.h.

Referenced by v_InitObject().

LibUtilities::TimeIntegrationSchemeOperators Nektar::UnsteadyAdvectionDiffusion::m_subStepIntegrationOps
protected
LibUtilities::TimeIntegrationWrapperSharedPtr Nektar::UnsteadyAdvectionDiffusion::m_subStepIntegrationScheme
protected
bool Nektar::UnsteadyAdvectionDiffusion::m_subSteppingScheme
protected

Definition at line 67 of file UnsteadyAdvectionDiffusion.h.

Referenced by v_InitObject(), and v_PreIntegrate().

NekDouble Nektar::UnsteadyAdvectionDiffusion::m_sVVCutoffRatio
protected

Definition at line 69 of file UnsteadyAdvectionDiffusion.h.

Referenced by DoImplicitSolve(), and v_InitObject().

NekDouble Nektar::UnsteadyAdvectionDiffusion::m_sVVDiffCoeff
protected

Definition at line 70 of file UnsteadyAdvectionDiffusion.h.

Referenced by DoImplicitSolve(), and v_InitObject().

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

Definition at line 74 of file UnsteadyAdvectionDiffusion.h.

Referenced by GetNormalVel().

bool Nektar::UnsteadyAdvectionDiffusion::m_useSpecVanVisc
protected

Definition at line 68 of file UnsteadyAdvectionDiffusion.h.

Referenced by DoImplicitSolve(), and v_InitObject().

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

Definition at line 164 of file UnsteadyAdvectionDiffusion.h.

Referenced by v_InitObject().