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

#include <RinglebFlow.h>

Inheritance diagram for Nektar::RinglebFlow:
[legend]

Public Member Functions

virtual ~RinglebFlow ()
 Destructor for EulerCFE class. More...
 
- Public Member Functions inherited from Nektar::EulerCFE
virtual ~EulerCFE ()
 Destructor for EulerCFE class. More...
 
- Public Member Functions inherited from Nektar::CompressibleFlowSystem
virtual ~CompressibleFlowSystem ()
 Destructor for CompressibleFlowSystem class. More...
 
NekDouble GetStabilityLimit (int n)
 Function to calculate the stability limit for DG/CG. More...
 
Array< OneD, NekDoubleGetStabilityLimitVector (const Array< OneD, int > &ExpOrder)
 Function to calculate the stability limit for DG/CG (a vector of them). More...
 
Array< OneD, NekDoubleGetElmtMinHP (void)
 Function to get estimate of min h/p factor per element. More...
 
virtual void GetPressure (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &pressure)
 Extract array with pressure from physfield. More...
 
virtual void GetDensity (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &density)
 Extract array with density from physfield. More...
 
virtual bool HasConstantDensity ()
 
virtual void GetVelocity (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &velocity)
 Extract array with velocity from physfield. 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 Member Functions inherited from Nektar::EulerCFE
static SolverUtils::EquationSystemSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Creates an instance of this class. More...
 

Static Public Attributes

static std::string className
 Name of class. More...
 
- Static Public Attributes inherited from Nektar::EulerCFE
static std::string className
 Name of class. More...
 
static std::string className2
 

Protected Member Functions

 RinglebFlow (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
virtual void v_GenerateSummary (SolverUtils::SummaryList &s)
 Print a summary of time stepping parameters. More...
 
virtual void v_EvaluateExactSolution (unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time=0.0)
 Get the exact solutions for isentropic vortex and Ringleb flow problems. More...
 
- Protected Member Functions inherited from Nektar::EulerCFE
 EulerCFE (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
- Protected Member Functions inherited from Nektar::CompressibleFlowSystem
 CompressibleFlowSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
virtual void v_InitObject ()
 Initialization object for CompressibleFlowSystem class. More...
 
void InitialiseParameters ()
 Load CFS parameters from the session file. More...
 
void InitAdvection ()
 Create advection and diffusion objects for CFS. More...
 
void DoOdeRhs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 Compute the right-hand side. More...
 
void DoOdeProjection (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 Compute the projection and call the method for imposing the boundary conditions in case of discontinuous projection. More...
 
void DoAdvection (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
 Compute the advection terms for the right-hand side. More...
 
void DoDiffusion (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
 Add the diffusions terms to the right-hand side. More...
 
void GetFluxVector (const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
 Return the flux vector for the compressible Euler equations. More...
 
void GetFluxVectorDeAlias (const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
 Return the flux vector for the compressible Euler equations by using the de-aliasing technique. More...
 
void SetBoundaryConditions (Array< OneD, Array< OneD, NekDouble > > &physarray, NekDouble time)
 
void GetElmtTimeStep (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &tstep)
 Calculate the maximum timestep on each element subject to CFL restrictions. More...
 
virtual NekDouble v_GetTimeStep (const Array< OneD, const Array< OneD, NekDouble > > &inarray)
 Calculate the maximum timestep subject to CFL restrictions. More...
 
virtual void v_SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 Set up logic for residual calculation. More...
 
NekDouble GetGamma ()
 
const Array< OneD, const Array< OneD, NekDouble > > & GetVecLocs ()
 
const Array< OneD, const Array< OneD, NekDouble > > & GetNormals ()
 
virtual void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 
virtual void v_DoDiffusion (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
 
virtual Array< OneD, NekDoublev_GetMaxStdVelocity ()
 Compute the advection velocity in the standard space for each element of the expansion. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::AdvectionSystem
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step)
 
- 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_PreIntegrate (int step)
 
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_Output (void)
 
virtual SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr v_GetPressure (void)
 

Private Member Functions

void GetExactRinglebFlow (int field, Array< OneD, NekDouble > &outarray)
 Ringleb Flow Test Case. More...
 

Friends

class MemoryManager< RinglebFlow >
 

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...
 
- Protected Attributes inherited from Nektar::CompressibleFlowSystem
SolverUtils::DiffusionSharedPtr m_diffusion
 
ArtificialDiffusionSharedPtr m_artificialDiffusion
 
Array< OneD, Array< OneD, NekDouble > > m_vecLocs
 
NekDouble m_gamma
 
std::string m_shockCaptureType
 
NekDouble m_filterAlpha
 
NekDouble m_filterExponent
 
NekDouble m_filterCutoff
 
bool m_useFiltering
 
bool m_useLocalTimeStep
 
VariableConverterSharedPtr m_varConv
 
std::vector< CFSBndCondSharedPtrm_bndConds
 
std::vector< SolverUtils::ForcingSharedPtrm_forcing
 
- 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...
 
- Static Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
static std::string equationSystemTypeLookupIds []
 

Detailed Description

Definition at line 43 of file RinglebFlow.h.

Constructor & Destructor Documentation

◆ ~RinglebFlow()

Nektar::RinglebFlow::~RinglebFlow ( )
virtual

Destructor for EulerCFE class.

Definition at line 59 of file RinglebFlow.cpp.

60  {
61  }

◆ RinglebFlow()

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

Definition at line 48 of file RinglebFlow.cpp.

51  : UnsteadySystem(pSession, pGraph),
52  EulerCFE(pSession, pGraph)
53  {
54  }
SOLVER_UTILS_EXPORT UnsteadySystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises UnsteadySystem class members.
EulerCFE(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Definition: EulerCFE.cpp:53

Member Function Documentation

◆ create()

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

Creates an instance of this class.

Definition at line 49 of file RinglebFlow.h.

References CellMLToNektar.cellml_metadata::p.

52  {
54  RinglebFlow>::AllocateSharedPtr(pSession, pGraph);
55  p->InitObject();
56  return p;
57  }
RinglebFlow(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Definition: RinglebFlow.cpp:48
std::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.

◆ GetExactRinglebFlow()

void Nektar::RinglebFlow::GetExactRinglebFlow ( int  field,
Array< OneD, NekDouble > &  outarray 
)
private

Ringleb Flow Test Case.

Compute the exact solution for the Ringleb flow problem.

Definition at line 88 of file RinglebFlow.cpp.

References ASSERTL0, Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::CompressibleFlowSystem::m_gamma, and class_topology::P.

Referenced by v_EvaluateExactSolution().

91  {
92  int nTotQuadPoints = GetTotPoints();
93 
94  Array<OneD, NekDouble> rho(nTotQuadPoints, 100.0);
95  Array<OneD, NekDouble> rhou(nTotQuadPoints);
96  Array<OneD, NekDouble> rhov(nTotQuadPoints);
97  Array<OneD, NekDouble> E(nTotQuadPoints);
98  Array<OneD, NekDouble> x(nTotQuadPoints);
99  Array<OneD, NekDouble> y(nTotQuadPoints);
100  Array<OneD, NekDouble> z(nTotQuadPoints);
101 
102  m_fields[0]->GetCoords(x, y, z);
103 
104  // Flow parameters
105  NekDouble c, k, phi, r, J, VV, pp, sint, P, ss;
106  NekDouble J11, J12, J21, J22, det;
107  NekDouble Fx, Fy;
108  NekDouble xi, yi;
109  NekDouble dV;
110  NekDouble dtheta;
111  NekDouble par1;
112  NekDouble theta = M_PI / 4.0;
113  NekDouble kExt = 0.7;
114  NekDouble V = kExt * sin(theta);
115  NekDouble toll = 1.0e-8;
116  NekDouble errV = 1.0;
117  NekDouble errTheta = 1.0;
118  NekDouble gamma = m_gamma;
119  NekDouble gamma_1_2 = (gamma - 1.0) / 2.0;
120 
121  for (int i = 0; i < nTotQuadPoints; ++i)
122  {
123  while ((abs(errV) > toll) || (abs(errTheta) > toll))
124  {
125  VV = V * V;
126  sint = sin(theta);
127  c = sqrt(1.0 - gamma_1_2 * VV);
128  k = V / sint;
129  phi = 1.0 / k;
130  pp = phi * phi;
131  J = 1.0 / c + 1.0 / (3.0 * c * c * c) +
132  1.0 / (5.0 * c * c * c * c * c) -
133  0.5 * log((1.0 + c) / (1.0 - c));
134 
135  r = pow(c, 1.0 / gamma_1_2);
136  xi = 1.0 / (2.0 * r) * (1.0 / VV - 2.0 * pp) + J / 2.0;
137  yi = phi / (r * V) * sqrt(1.0 - VV * pp);
138  par1 = 25.0 - 5.0 * VV;
139  ss = sint * sint;
140 
141  Fx = xi - x[i];
142  Fy = yi - y[i];
143 
144  J11 = 39062.5 / pow(par1, 3.5) * (1.0 / VV - 2.0 / VV * ss) *
145  V + 1562.5 / pow(par1, 2.5) * (-2.0 / (VV * V) + 4.0 /
146  (VV * V) * ss) + 12.5 / pow(par1, 1.5) * V + 312.5 /
147  pow(par1, 2.5) * V + 7812.5 / pow(par1, 3.5) * V -
148  0.25 * (-1.0 / pow(par1, 0.5) * V/(1.0 - 0.2 *
149  pow(par1, 0.5)) - (1.0 + 0.2 * pow(par1, 0.5)) /
150  pow((1.0 - 0.2 * pow(par1, 0.5)), 2.0) /
151  pow(par1, 0.5) * V) / (1.0 + 0.2 * pow(par1, 0.5)) *
152  (1.0 - 0.2 * pow(par1, 0.5));
153 
154  J12 = -6250.0 / pow(par1, 2.5) / VV * sint * cos(theta);
155  J21 = -6250.0 / (VV * V) * sint /
156  pow(par1, 2.5) * pow((1.0 - ss), 0.5) +
157  78125.0 / V * sint / pow(par1, 3.5) *
158  pow((1.0 - ss), 0.5);
159 
160  // the matrix is singular when theta = pi/2
161  if(abs(y[i])<toll && abs(cos(theta))<toll)
162  {
163  J22 = -39062.5 / pow(par1, 3.5) / V + 3125 /
164  pow(par1, 2.5) / (VV * V) + 12.5 / pow(par1, 1.5) *
165  V + 312.5 / pow(par1, 2.5) * V + 7812.5 /
166  pow(par1, 3.5) * V - 0.25 * (-1.0 / pow(par1, 0.5) *
167  V / (1.0 - 0.2 * pow(par1, 0.5)) - (1.0 + 0.2 *
168  pow(par1, 0.5)) / pow((1.0 - 0.2 *
169  pow(par1, 0.5)), 2.0) / pow(par1, 0.5) * V) /
170  (1.0 + 0.2 * pow(par1, 0.5)) * (1.0 - 0.2 *
171  pow(par1,0.5));
172 
173  // dV = -dV/dx * Fx
174  dV = -1.0 / J22 * Fx;
175  dtheta = 0.0;
176  theta = M_PI / 2.0;
177  }
178  else
179  {
180  J22 = 3125.0 / VV * cos(theta) / pow(par1, 2.5) *
181  pow((1.0 - ss), 0.5) - 3125.0 / VV * ss /
182  pow(par1, 2.5) / pow((1.0 - ss), 0.5) * cos(theta);
183 
184  det = -1.0 / (J11 * J22 - J12 * J21);
185 
186  // [dV dtheta]' = -[invJ]*[Fx Fy]'
187  dV = det * ( J22 * Fx - J12 * Fy);
188  dtheta = det * (-J21 * Fx + J11 * Fy);
189  }
190 
191  V = V + dV;
192  theta = theta + dtheta;
193 
194  errV = abs(dV);
195  errTheta = abs(dtheta);
196 
197  }
198 
199  c = sqrt(1.0 - gamma_1_2 * VV);
200  r = pow(c, 1.0 / gamma_1_2);
201 
202  rho[i] = r;
203  rhou[i] = rho[i] * V * cos(theta);
204  rhov[i] = rho[i] * V * sin(theta);
205  P = (c * c) * rho[i] / gamma;
206  E[i] = P / (gamma - 1.0) + 0.5 *
207  (rhou[i] * rhou[i] / rho[i] + rhov[i] * rhov[i] / rho[i]);
208 
209  // Resetting the guess value
210  errV = 1.0;
211  errTheta = 1.0;
212  theta = M_PI/4.0;
213  V = kExt*sin(theta);
214  }
215 
216  switch (field)
217  {
218  case 0:
219  outarray = rho;
220  break;
221  case 1:
222  outarray = rhou;
223  break;
224  case 2:
225  outarray = rhov;
226  break;
227  case 3:
228  outarray = E;
229  break;
230  default:
231  ASSERTL0(false, "Error in variable number!");
232  break;
233  }
234  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
SOLVER_UTILS_EXPORT int GetTotPoints()
double NekDouble
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.

◆ v_EvaluateExactSolution()

void Nektar::RinglebFlow::v_EvaluateExactSolution ( unsigned int  field,
Array< OneD, NekDouble > &  outfield,
const NekDouble  time = 0.0 
)
protectedvirtual

Get the exact solutions for isentropic vortex and Ringleb flow problems.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 76 of file RinglebFlow.cpp.

References GetExactRinglebFlow().

80  {
81  boost::ignore_unused(time);
82  GetExactRinglebFlow( field, outfield);
83  }
void GetExactRinglebFlow(int field, Array< OneD, NekDouble > &outarray)
Ringleb Flow Test Case.
Definition: RinglebFlow.cpp:88

◆ v_GenerateSummary()

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

Print a summary of time stepping parameters.

Print out a summary with some relevant information.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 66 of file RinglebFlow.cpp.

References Nektar::SolverUtils::AddSummaryItem(), and Nektar::SolverUtils::UnsteadySystem::v_GenerateSummary().

67  {
69  SolverUtils::AddSummaryItem(s, "Problem Type", "RinglebFlow");
70  }
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s)
Print a summary of time stepping parameters.
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:49

Friends And Related Function Documentation

◆ MemoryManager< RinglebFlow >

friend class MemoryManager< RinglebFlow >
friend

Definition at line 46 of file RinglebFlow.h.

Member Data Documentation

◆ className

string Nektar::RinglebFlow::className
static
Initial value:
=
"RinglebFlow", RinglebFlow::create,
"Euler equations for Ringleb flow.")

Name of class.

Definition at line 59 of file RinglebFlow.h.