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

#include <EulerCFE.h>

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

Public Member Functions

virtual ~EulerCFE ()
 problem type selector
- Public Member Functions inherited from Nektar::CompressibleFlowSystem
virtual ~CompressibleFlowSystem ()
 Destructor for CompressibleFlowSystem class.
NekDouble GetStabilityLimit (int n)
 Function to calculate the stability limit for DG/CG.
Array< OneD, NekDoubleGetStabilityLimitVector (const Array< OneD, int > &ExpOrder)
 Function to calculate the stability limit for DG/CG (a vector of them).
- Public Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
virtual SOLVER_UTILS_EXPORT ~UnsteadySystem ()
 Destructor.
SOLVER_UTILS_EXPORT NekDouble GetTimeStep (const Array< OneD, const Array< OneD, NekDouble > > &inarray)
 Calculate the larger time-step mantaining the problem stable.
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor.
SOLVER_UTILS_EXPORT void SetUpTraceNormals (void)
SOLVER_UTILS_EXPORT void InitObject ()
 Initialises the members of this object.
SOLVER_UTILS_EXPORT void DoInitialise ()
 Perform any initialisation necessary before solving the problem.
SOLVER_UTILS_EXPORT void DoSolve ()
 Solve the problem.
SOLVER_UTILS_EXPORT void TransCoeffToPhys ()
 Transform from coefficient to physical space.
SOLVER_UTILS_EXPORT void TransPhysToCoeff ()
 Transform from physical to coefficient space.
SOLVER_UTILS_EXPORT void Output ()
 Perform output operations after solve.
SOLVER_UTILS_EXPORT NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation.
SOLVER_UTILS_EXPORT std::string GetSessionName ()
 Get Session name.
template<class T >
boost::shared_ptr< T > as ()
SOLVER_UTILS_EXPORT void ResetSessionName (std::string newname)
 Reset Session name.
SOLVER_UTILS_EXPORT
LibUtilities::SessionReaderSharedPtr 
GetSession ()
 Get Session name.
SOLVER_UTILS_EXPORT
MultiRegions::ExpListSharedPtr 
GetPressure ()
 Get pressure field if available.
SOLVER_UTILS_EXPORT void PrintSummary (std::ostream &out)
 Print a summary of parameters and solver characteristics.
SOLVER_UTILS_EXPORT void SetLambda (NekDouble lambda)
 Set parameter m_lambda.
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.
SOLVER_UTILS_EXPORT void EvaluateFunction (std::vector< std::string > pFieldNames, Array< OneD, Array< OneD, NekDouble > > &pFields, const std::string &pName, const int domain=0)
 Populate given fields with the function from session.
SOLVER_UTILS_EXPORT void EvaluateFunction (std::vector< std::string > pFieldNames, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const std::string &pName, const int domain=0)
 Populate given fields with the function from session.
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.
SOLVER_UTILS_EXPORT void InitialiseBaseFlow (Array< OneD, Array< OneD, NekDouble > > &base)
 Perform initialisation of the base flow.
SOLVER_UTILS_EXPORT void SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 Initialise the data in the dependent fields.
SOLVER_UTILS_EXPORT void EvaluateExactSolution (int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 Evaluates an exact solution.
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.
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, bool Normalised=false)
 Compute the L2 error of the fields.
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].
SOLVER_UTILS_EXPORT void WeakAdvectionGreensDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\nabla \phi \cdot F) $.
SOLVER_UTILS_EXPORT void WeakAdvectionDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\phi, \nabla \cdot F) $.
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) $.
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.
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.
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.
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n)
 Write checkpoint file of m_fields.
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.
SOLVER_UTILS_EXPORT void Checkpoint_BaseFlow (const int n)
 Write base flow file of m_fields.
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname)
 Write field data to the given filename.
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.
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 Input field data from the given file.
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.
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.
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.
SOLVER_UTILS_EXPORT void ScanForHistoryPoints ()
 Builds map of which element holds each history point.
SOLVER_UTILS_EXPORT void WriteHistoryData (std::ostream &out)
 Probe each history point and write to file.
SOLVER_UTILS_EXPORT void SessionSummary (SummaryList &vSummary)
 Write out a session summary.
SOLVER_UTILS_EXPORT Array
< OneD,
MultiRegions::ExpListSharedPtr > & 
UpdateFields ()
SOLVER_UTILS_EXPORT
LibUtilities::FieldMetaDataMap
UpdateFieldMetaDataMap ()
 Get hold of FieldInfoMap so it can be updated.
SOLVER_UTILS_EXPORT NekDouble GetFinalTime ()
 Return final time.
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 SetStepsToOne ()
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.

Static Public Member Functions

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

Public Attributes

ProblemType m_problemType

Static Public Attributes

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

Protected Member Functions

 EulerCFE (const LibUtilities::SessionReaderSharedPtr &pSession)
virtual void v_InitObject ()
 Initialization object for CompressibleFlowSystem class.
virtual void v_GenerateSummary (SolverUtils::SummaryList &s)
 Print a summary of time stepping parameters.
void DoOdeRhs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 Compute the right-hand side.
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.
virtual void v_SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 Set the initial conditions.
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.
- Protected Member Functions inherited from Nektar::CompressibleFlowSystem
 CompressibleFlowSystem (const LibUtilities::SessionReaderSharedPtr &pSession)
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.
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.
void GetViscousFluxVector (const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &derivatives, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &viscousTensor)
 Return the flux vector for the LDG diffusion problem.
void GetFluxVectorPDESC (const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
void GetViscousFluxVectorDeAlias (const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &derivatives, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &viscousTensor)
 Return the flux vector for the LDG diffusion problem.
void WallBC (int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &physarray)
 Wall boundary conditions for compressible flow problems.
void WallViscousBC (int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &physarray)
 Wall boundary conditions for viscous compressible flow problems.
void SymmetryBC (int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &physarray)
 Symmetry boundary conditions for compressible flow problems.
void RiemannInvariantBC (int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &physarray)
 Outflow characteristic boundary conditions for compressible flow problems.
void ExtrapOrder0BC (int bcRegion, int cnt, Array< OneD, Array< OneD, NekDouble > > &physarray)
 Extrapolation of order 0 for all the variables such that, at the boundaries, a trivial Riemann problem is solved.
void GetVelocityVector (const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &velocity)
void GetSoundSpeed (const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &pressure, Array< OneD, NekDouble > &soundspeed)
void GetMach (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &soundspeed, Array< OneD, NekDouble > &mach)
void GetTemperature (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &pressure, Array< OneD, NekDouble > &temperature)
void GetPressure (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &pressure)
 Calculate the pressure field $ p = (\gamma-1)(E-\frac{1}{2}\rho\| \mathbf{v} \|^2) $ assuming an ideal gas law.
void GetPressure (const Array< OneD, const Array< OneD, NekDouble > > &physfield, const Array< OneD, const Array< OneD, NekDouble > > &velocity, Array< OneD, NekDouble > &pressure)
void GetEnthalpy (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &pressure, Array< OneD, NekDouble > &enthalpy)
void GetEntropy (const Array< OneD, const Array< OneD, NekDouble > > &physfield, const Array< OneD, const NekDouble > &pressure, const Array< OneD, const NekDouble > &temperature, Array< OneD, NekDouble > &entropy)
void GetSmoothArtificialViscosity (const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &eps_bar)
void GetDynamicViscosity (const Array< OneD, const NekDouble > &temperature, Array< OneD, NekDouble > &mu)
void GetStdVelocity (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &stdV)
void GetSensor (const Array< OneD, const Array< OneD, NekDouble > > &physarray, Array< OneD, NekDouble > &Sensor, Array< OneD, NekDouble > &SensorKappa)
void GetElementDimensions (Array< OneD, Array< OneD, NekDouble > > &outarray, Array< OneD, NekDouble > &hmin)
void GetAbsoluteVelocity (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &Vtot)
void GetArtificialDynamicViscosity (const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &mu_var)
void SetVarPOrderElmt (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &PolyOrder)
void GetForcingTerm (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > outarrayForcing)
virtual 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.
NekDouble GetGasConstant ()
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)
- Protected Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
SOLVER_UTILS_EXPORT UnsteadySystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 Initialises UnsteadySystem class members.
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control.
virtual SOLVER_UTILS_EXPORT void v_DoSolve ()
 Solves an unsteady problem.
virtual SOLVER_UTILS_EXPORT void v_DoInitialise ()
 Sets up initial conditions.
virtual SOLVER_UTILS_EXPORT void v_AppendOutput1D (Array< OneD, Array< OneD, NekDouble > > &solution1D)
 Print the solution at each solution point in a txt file.
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 bool v_PreIntegrate (int step)
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step)
SOLVER_UTILS_EXPORT void CheckForRestartTime (NekDouble &time)
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 Initialises EquationSystem class members.
int nocase_cmp (const string &s1, const string &s2)
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time.
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.
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.
virtual SOLVER_UTILS_EXPORT void v_TransCoeffToPhys ()
 Virtual function for transformation to physical space.
virtual SOLVER_UTILS_EXPORT void v_TransPhysToCoeff ()
 Virtual function for transformation to coefficient space.
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)

Private Member Functions

void SetBoundaryConditions (Array< OneD, Array< OneD, NekDouble > > &physarray, NekDouble time)
 Set boundary conditions which can be: a) Wall and Symmerty BCs implemented at CompressibleFlowSystem level since they are compressible solver specific; b) Isentropic vortex and Ringleb flow BCs implemented at EulerCFE level since they are Euler solver specific; c) Time dependent BCs.
void EvaluateIsentropicVortex (const Array< OneD, NekDouble > &x, const Array< OneD, NekDouble > &y, const Array< OneD, NekDouble > &z, Array< OneD, Array< OneD, NekDouble > > &u, NekDouble time, const int o=0)
 Isentropic Vortex Test Case.
void GetExactIsentropicVortex (int field, Array< OneD, NekDouble > &outarray, NekDouble time)
 Compute the exact solution for the isentropic vortex problem.
void SetInitialIsentropicVortex (NekDouble initialtime)
 Set the initial condition for the isentropic vortex problem.
void SetBoundaryIsentropicVortex (int bcRegion, NekDouble time, int cnt, Array< OneD, Array< OneD, NekDouble > > &physarray)
 Set the boundary conditions for the isentropic vortex problem.
void GetExactRinglebFlow (int field, Array< OneD, NekDouble > &outarray)
 Ringleb Flow Test Case.
void SetInitialRinglebFlow (void)
 Set the initial condition for the Ringleb flow problem.
void SetBoundaryRinglebFlow (int bcRegion, NekDouble time, int cnt, Array< OneD, Array< OneD, NekDouble > > &physarray)
 Set the boundary conditions for the Ringleb flow problem.

Friends

class MemoryManager< EulerCFE >

Additional Inherited Members

- 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::RiemannSolverSharedPtr m_riemannSolver
SolverUtils::RiemannSolverSharedPtr m_riemannSolverLDG
SolverUtils::AdvectionSharedPtr m_advection
SolverUtils::DiffusionSharedPtr m_diffusion
Array< OneD, Array< OneD,
NekDouble > > 
m_vecLocs
NekDouble m_gamma
NekDouble m_pInf
NekDouble m_rhoInf
NekDouble m_uInf
NekDouble m_vInf
NekDouble m_wInf
NekDouble m_gasConstant
NekDouble m_Twall
std::string m_ViscosityType
std::string m_shockCaptureType
std::string m_EqTypeStr
NekDouble m_mu
NekDouble m_Skappa
NekDouble m_Kappa
NekDouble m_mu0
NekDouble m_FacL
NekDouble m_FacH
NekDouble m_eps_max
NekDouble m_thermalConductivity
NekDouble m_Cp
NekDouble m_C1
NekDouble m_C2
NekDouble m_hFactor
NekDouble m_Prandtl
StdRegions::StdQuadExpSharedPtr m_OrthoQuadExp
StdRegions::StdHexExpSharedPtr m_OrthoHexExp
bool m_smoothDiffusion

Detailed Description

Definition at line 58 of file EulerCFE.h.

Constructor & Destructor Documentation

Nektar::EulerCFE::~EulerCFE ( )
virtual

problem type selector

Destructor for EulerCFE class.

Definition at line 95 of file EulerCFE.cpp.

{
}
Nektar::EulerCFE::EulerCFE ( const LibUtilities::SessionReaderSharedPtr pSession)
protected

Definition at line 52 of file EulerCFE.cpp.

{
}

Member Function Documentation

static SolverUtils::EquationSystemSharedPtr Nektar::EulerCFE::create ( const LibUtilities::SessionReaderSharedPtr pSession)
inlinestatic

Creates an instance of this class.

Reimplemented from Nektar::CompressibleFlowSystem.

Definition at line 64 of file EulerCFE.h.

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

Compute the projection and call the method for imposing the boundary conditions in case of discontinuous projection.

Definition at line 182 of file EulerCFE.cpp.

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

Referenced by v_InitObject().

{
int i;
int nvariables = inarray.num_elements();
{
{
// Just copy over array
int npoints = GetNpoints();
for(i = 0; i < nvariables; ++i)
{
Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
}
SetBoundaryConditions(outarray, time);
break;
}
{
ASSERTL0(false, "No Continuous Galerkin for Euler equations");
break;
}
default:
ASSERTL0(false, "Unknown projection scheme");
break;
}
}
void Nektar::EulerCFE::DoOdeRhs ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time 
)
protected

Compute the right-hand side.

Definition at line 158 of file EulerCFE.cpp.

References Nektar::SolverUtils::EquationSystem::GetNpoints(), Nektar::CompressibleFlowSystem::m_advection, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_spacedim, and Vmath::Neg().

Referenced by v_InitObject().

{
int i;
int nvariables = inarray.num_elements();
int npoints = GetNpoints();
Array<OneD, Array<OneD, NekDouble> > advVel(m_spacedim);
m_advection->Advect(nvariables, m_fields, advVel, inarray,
outarray, time);
for (i = 0; i < nvariables; ++i)
{
Vmath::Neg(npoints, outarray[i], 1);
}
}
void Nektar::EulerCFE::EvaluateIsentropicVortex ( const Array< OneD, NekDouble > &  x,
const Array< OneD, NekDouble > &  y,
const Array< OneD, NekDouble > &  z,
Array< OneD, Array< OneD, NekDouble > > &  u,
NekDouble  time,
const int  o = 0 
)
private

Isentropic Vortex Test Case.

Definition at line 331 of file EulerCFE.cpp.

References Nektar::CompressibleFlowSystem::m_gamma, Nektar::SolverUtils::EquationSystem::m_spacedim, and Vmath::Zero().

Referenced by GetExactIsentropicVortex(), SetBoundaryIsentropicVortex(), and SetInitialIsentropicVortex().

{
int nq = x.num_elements();
// Flow parameters
const NekDouble x0 = 5.0;
const NekDouble y0 = 0.0;
const NekDouble beta = 5.0;
const NekDouble u0 = 1.0;
const NekDouble v0 = 0.5;
const NekDouble gamma = m_gamma;
NekDouble r, xbar, ybar, tmp;
NekDouble fac = 1.0/(16.0*gamma*M_PI*M_PI);
// In 3D zero rhow field.
if (m_spacedim == 3)
{
Vmath::Zero(nq, &u[3][o], 1);
}
// Fill storage
for (int i = 0; i < nq; ++i)
{
xbar = x[i] - u0*time - x0;
ybar = y[i] - v0*time - y0;
r = sqrt(xbar*xbar + ybar*ybar);
tmp = beta*exp(1-r*r);
u[0][i+o] = pow(1.0 - (gamma-1.0)*tmp*tmp*fac, 1.0/(gamma-1.0));
u[1][i+o] = u[0][i+o]*(u0 - tmp*ybar/(2*M_PI));
u[2][i+o] = u[0][i+o]*(v0 + tmp*xbar/(2*M_PI));
u[m_spacedim+1][i+o] = pow(u[0][i+o], gamma)/(gamma-1.0) +
0.5*(u[1][i+o]*u[1][i+o] + u[2][i+o]*u[2][i+o]) / u[0][i+o];
}
}
void Nektar::EulerCFE::GetExactIsentropicVortex ( int  field,
Array< OneD, NekDouble > &  outarray,
NekDouble  time 
)
private

Compute the exact solution for the isentropic vortex problem.

Definition at line 375 of file EulerCFE.cpp.

References EvaluateIsentropicVortex(), Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_spacedim, and Vmath::Vcopy().

Referenced by v_EvaluateExactSolution().

{
int nTotQuadPoints = GetTotPoints();
Array<OneD, NekDouble> x(nTotQuadPoints);
Array<OneD, NekDouble> y(nTotQuadPoints);
Array<OneD, NekDouble> z(nTotQuadPoints);
Array<OneD, Array<OneD, NekDouble> > u(m_spacedim+2);
m_fields[0]->GetCoords(x, y, z);
for (int i = 0; i < m_spacedim + 2; ++i)
{
u[i] = Array<OneD, NekDouble>(nTotQuadPoints);
}
EvaluateIsentropicVortex(x, y, z, u, time);
Vmath::Vcopy(nTotQuadPoints, u[field], 1, outarray, 1);
}
void Nektar::EulerCFE::GetExactRinglebFlow ( int  field,
Array< OneD, NekDouble > &  outarray 
)
private

Ringleb Flow Test Case.

Compute the exact solution for the Ringleb flow problem.

Definition at line 481 of file EulerCFE.cpp.

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

Referenced by v_EvaluateExactSolution().

{
int nTotQuadPoints = GetTotPoints();
Array<OneD, NekDouble> rho(nTotQuadPoints, 100.0);
Array<OneD, NekDouble> rhou(nTotQuadPoints);
Array<OneD, NekDouble> rhov(nTotQuadPoints);
Array<OneD, NekDouble> E(nTotQuadPoints);
Array<OneD, NekDouble> x(nTotQuadPoints);
Array<OneD, NekDouble> y(nTotQuadPoints);
Array<OneD, NekDouble> z(nTotQuadPoints);
m_fields[0]->GetCoords(x, y, z);
// Flow parameters
NekDouble c, k, phi, r, J, VV, pp, sint, P, ss;
NekDouble J11, J12, J21, J22, det;
NekDouble Fx, Fy;
NekDouble xi, yi;
NekDouble dtheta;
NekDouble par1;
NekDouble theta = M_PI / 4.0;
NekDouble kExt = 0.7;
NekDouble V = kExt * sin(theta);
NekDouble toll = 1.0e-8;
NekDouble errV = 1.0;
NekDouble errTheta = 1.0;
NekDouble gamma = m_gamma;
NekDouble gamma_1_2 = (gamma - 1.0) / 2.0;
for (int i = 0; i < nTotQuadPoints; ++i)
{
while ((abs(errV) > toll) || (abs(errTheta) > toll))
{
VV = V * V;
sint = sin(theta);
c = sqrt(1.0 - gamma_1_2 * VV);
k = V / sint;
phi = 1.0 / k;
pp = phi * phi;
J = 1.0 / c + 1.0 / (3.0 * c * c * c) +
1.0 / (5.0 * c * c * c * c * c) -
0.5 * log((1.0 + c) / (1.0 - c));
r = pow(c, 1.0 / gamma_1_2);
xi = 1.0 / (2.0 * r) * (1.0 / VV - 2.0 * pp) + J / 2.0;
yi = phi / (r * V) * sqrt(1.0 - VV * pp);
par1 = 25.0 - 5.0 * VV;
ss = sint * sint;
Fx = xi - x[i];
Fy = yi - y[i];
J11 = 39062.5 / pow(par1, 3.5) * (1.0 / VV - 2.0 / VV * ss) *
V + 1562.5 / pow(par1, 2.5) * (-2.0 / (VV * V) + 4.0 /
(VV * V) * ss) + 12.5 / pow(par1, 1.5) * V + 312.5 /
pow(par1, 2.5) * V + 7812.5 / pow(par1, 3.5) * V -
0.25 * (-1.0 / pow(par1, 0.5) * V/(1.0 - 0.2 *
pow(par1, 0.5)) - (1.0 + 0.2 * pow(par1, 0.5)) /
pow((1.0 - 0.2 * pow(par1, 0.5)), 2.0) /
pow(par1, 0.5) * V) / (1.0 + 0.2 * pow(par1, 0.5)) *
(1.0 - 0.2 * pow(par1, 0.5));
J12 = -6250.0 / pow(par1, 2.5) / VV * sint * cos(theta);
J21 = -6250.0 / (VV * V) * sint /
pow(par1, 2.5) * pow((1.0 - ss), 0.5) +
78125.0 / V * sint / pow(par1, 3.5) *
pow((1.0 - ss), 0.5);
// the matrix is singular when theta = pi/2
if(abs(y[i])<toll && abs(cos(theta))<toll)
{
J22 = -39062.5 / pow(par1, 3.5) / V + 3125 /
pow(par1, 2.5) / (VV * V) + 12.5 / pow(par1, 1.5) *
V + 312.5 / pow(par1, 2.5) * V + 7812.5 /
pow(par1, 3.5) * V - 0.25 * (-1.0 / pow(par1, 0.5) *
V / (1.0 - 0.2 * pow(par1, 0.5)) - (1.0 + 0.2 *
pow(par1, 0.5)) / pow((1.0 - 0.2 *
pow(par1, 0.5)), 2.0) / pow(par1, 0.5) * V) /
(1.0 + 0.2 * pow(par1, 0.5)) * (1.0 - 0.2 *
pow(par1,0.5));
// dV = -dV/dx * Fx
dV = -1.0 / J22 * Fx;
dtheta = 0.0;
theta = M_PI / 2.0;
}
else
{
J22 = 3125.0 / VV * cos(theta) / pow(par1, 2.5) *
pow((1.0 - ss), 0.5) - 3125.0 / VV * ss /
pow(par1, 2.5) / pow((1.0 - ss), 0.5) * cos(theta);
det = -1.0 / (J11 * J22 - J12 * J21);
// [dV dtheta]' = -[invJ]*[Fx Fy]'
dV = det * ( J22 * Fx - J12 * Fy);
dtheta = det * (-J21 * Fx + J11 * Fy);
}
V = V + dV;
theta = theta + dtheta;
errV = abs(dV);
errTheta = abs(dtheta);
}
c = sqrt(1.0 - gamma_1_2 * VV);
r = pow(c, 1.0 / gamma_1_2);
rho[i] = r;
rhou[i] = rho[i] * V * cos(theta);
rhov[i] = rho[i] * V * sin(theta);
P = (c * c) * rho[i] / gamma;
E[i] = P / (gamma - 1.0) + 0.5 *
(rhou[i] * rhou[i] / rho[i] + rhov[i] * rhov[i] / rho[i]);
// Resetting the guess value
errV = 1.0;
errTheta = 1.0;
theta = M_PI/4.0;
V = kExt*sin(theta);
}
switch (field)
{
case 0:
outarray = rho;
break;
case 1:
outarray = rhou;
break;
case 2:
outarray = rhov;
break;
case 3:
outarray = E;
break;
default:
ASSERTL0(false, "Error in variable number!");
break;
}
}
void Nektar::EulerCFE::SetBoundaryConditions ( Array< OneD, Array< OneD, NekDouble > > &  inarray,
NekDouble  time 
)
private

Set boundary conditions which can be: a) Wall and Symmerty BCs implemented at CompressibleFlowSystem level since they are compressible solver specific; b) Isentropic vortex and Ringleb flow BCs implemented at EulerCFE level since they are Euler solver specific; c) Time dependent BCs.

Parameters
inarray,:fields array.
time,:time.

Definition at line 227 of file EulerCFE.cpp.

References ASSERTL0, Nektar::SpatialDomains::eExtrapOrder0, Nektar::SpatialDomains::eIsentropicVortex, Nektar::SpatialDomains::eRiemannInvariant, Nektar::SpatialDomains::eRinglebFlow, Nektar::SpatialDomains::eSymmetry, Nektar::SpatialDomains::eTimeDependent, Nektar::SpatialDomains::eWall, Nektar::SpatialDomains::eWallViscous, Nektar::CompressibleFlowSystem::ExtrapOrder0BC(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_session, Nektar::CompressibleFlowSystem::RiemannInvariantBC(), SetBoundaryIsentropicVortex(), SetBoundaryRinglebFlow(), Nektar::CompressibleFlowSystem::SymmetryBC(), and Nektar::CompressibleFlowSystem::WallBC().

Referenced by DoOdeProjection().

{
std::string varName;
int nvariables = m_fields.num_elements();
int cnt = 0;
// Loop over Boundary Regions
for (int n = 0; n < m_fields[0]->GetBndConditions().num_elements(); ++n)
{
// Wall Boundary Condition
if (m_fields[0]->GetBndConditions()[n]->GetUserDefined() ==
{
WallBC(n, cnt, inarray);
}
// Wall Boundary Condition
if (m_fields[0]->GetBndConditions()[n]->GetUserDefined() ==
{
ASSERTL0(false, "WallViscous is a wrong bc for the "
"Euler equations");
}
// Symmetric Boundary Condition
if (m_fields[0]->GetBndConditions()[n]->GetUserDefined() ==
{
SymmetryBC(n, cnt, inarray);
}
// Riemann invariant characteristic Boundary Condition (CBC)
if (m_fields[0]->GetBndConditions()[n]->GetUserDefined() ==
{
RiemannInvariantBC(n, cnt, inarray);
}
// Extrapolation of the data at the boundaries
if (m_fields[0]->GetBndConditions()[n]->GetUserDefined() ==
{
ExtrapOrder0BC(n, cnt, inarray);
}
// Time Dependent Boundary Condition (specified in meshfile)
if (m_fields[0]->GetBndConditions()[n]->GetUserDefined() ==
{
for (int i = 0; i < nvariables; ++i)
{
varName = m_session->GetVariable(i);
m_fields[i]->EvaluateBoundaryConditions(time, varName);
}
}
// Isentropic Vortex Boundary Condition
if (m_fields[0]->GetBndConditions()[n]->GetUserDefined() ==
{
SetBoundaryIsentropicVortex(n, time, cnt, inarray);
}
// Ringleb Flow Inflow and outflow Boundary Condition
if (m_fields[0]->GetBndConditions()[n]->GetUserDefined() ==
{
SetBoundaryRinglebFlow(n, time, cnt, inarray);
}
cnt += m_fields[0]->GetBndCondExpansions()[n]->GetExpSize();
}
}
void Nektar::EulerCFE::SetBoundaryIsentropicVortex ( int  bcRegion,
NekDouble  time,
int  cnt,
Array< OneD, Array< OneD, NekDouble > > &  physarray 
)
private

Set the boundary conditions for the isentropic vortex problem.

Definition at line 431 of file EulerCFE.cpp.

References EvaluateIsentropicVortex(), Nektar::SolverUtils::EquationSystem::GetTraceTotPoints(), Nektar::SolverUtils::EquationSystem::m_fields, and Vmath::Vcopy().

Referenced by SetBoundaryConditions().

{
int nvariables = physarray.num_elements();
int nTraceNumPoints = GetTraceTotPoints();
Array<OneD, Array<OneD, NekDouble> > Fwd(nvariables);
const Array<OneD, const int> &bndTraceMap = m_fields[0]->GetTraceBndMap();
// Get physical values of the forward trace (from exp to phys)
for (int i = 0; i < nvariables; ++i)
{
Fwd[i] = Array<OneD, NekDouble>(nTraceNumPoints);
m_fields[i]->ExtractTracePhys(physarray[i], Fwd[i]);
}
int id2, e_max;
e_max = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
for(int e = 0; e < e_max; ++e)
{
int npoints = m_fields[0]->
GetBndCondExpansions()[bcRegion]->GetExp(e)->GetTotPoints();
int id1 = m_fields[0]->
GetBndCondExpansions()[bcRegion]->GetPhys_Offset(e);
id2 = m_fields[0]->GetTrace()->GetPhys_Offset(bndTraceMap[cnt++]);
Array<OneD,NekDouble> x(npoints, 0.0);
Array<OneD,NekDouble> y(npoints, 0.0);
Array<OneD,NekDouble> z(npoints, 0.0);
m_fields[0]->GetBndCondExpansions()[bcRegion]->
GetExp(e)->GetCoords(x, y, z);
EvaluateIsentropicVortex(x, y, z, Fwd, time, id2);
for (int i = 0; i < nvariables; ++i)
{
Vmath::Vcopy(npoints, &Fwd[i][id2], 1,
&(m_fields[i]->GetBndCondExpansions()[bcRegion]->
UpdatePhys())[id1], 1);
}
}
}
void Nektar::EulerCFE::SetBoundaryRinglebFlow ( int  bcRegion,
NekDouble  time,
int  cnt,
Array< OneD, Array< OneD, NekDouble > > &  physarray 
)
private

Set the boundary conditions for the Ringleb flow problem.

Definition at line 876 of file EulerCFE.cpp.

References Nektar::LibUtilities::eFunctionTypeFile, Nektar::SolverUtils::EquationSystem::eHomogeneous1D, Nektar::SolverUtils::EquationSystem::GetTraceTotPoints(), Nektar::SolverUtils::EquationSystem::m_expdim, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::CompressibleFlowSystem::m_gamma, Nektar::SolverUtils::EquationSystem::m_HomogeneousType, Nektar::SolverUtils::EquationSystem::m_session, and Vmath::Vcopy().

Referenced by SetBoundaryConditions().

{
int nvariables = physarray.num_elements();
int nTraceNumPoints = GetTraceTotPoints();
Array<OneD, Array<OneD, NekDouble> > Fwd(nvariables);
// For 3DHomogenoeus1D
int n_planes = 1;
{
int nPointsTot = m_fields[0]->GetTotPoints();
int nPointsTot_plane = m_fields[0]->GetPlane(0)->GetTotPoints();
n_planes = nPointsTot/nPointsTot_plane;
nTraceNumPoints = nTraceNumPoints * n_planes;
}
// Get physical values of the forward trace (from exp to phys)
for (int i = 0; i < nvariables; ++i)
{
Fwd[i] = Array<OneD, NekDouble>(nTraceNumPoints);
m_fields[i]->ExtractTracePhys(physarray[i], Fwd[i]);
}
int id2, id2_plane, e_max;
e_max = m_fields[0]->GetBndCondExpansions()[bcRegion]->GetExpSize();
for(int e = 0; e < e_max; ++e)
{
int npoints = m_fields[0]->
GetBndCondExpansions()[bcRegion]->GetExp(e)->GetTotPoints();
int id1 = m_fields[0]->
GetBndCondExpansions()[bcRegion]->GetPhys_Offset(e);
// For 3DHomogenoeus1D
{
int cnt_plane = cnt/n_planes;
int e_plane;
int e_max_plane = e_max/n_planes;
int nTracePts_plane = GetTraceTotPoints();
int planeID = floor((e + 0.5 )/ e_max_plane );
e_plane = e - e_max_plane*planeID;
id2_plane = m_fields[0]->GetTrace()->GetPhys_Offset(
m_fields[0]->GetTraceMap()->
GetBndCondCoeffsToGlobalCoeffsMap(
cnt_plane + e_plane));
id2 = id2_plane + planeID*nTracePts_plane;
}
else // For general case
{
id2 = m_fields[0]->
GetTrace()->GetPhys_Offset(m_fields[0]->GetTraceMap()->
GetBndCondTraceToGlobalTraceMap(cnt++));
}
Array<OneD,NekDouble> x0(npoints, 0.0);
Array<OneD,NekDouble> x1(npoints, 0.0);
Array<OneD,NekDouble> x2(npoints, 0.0);
m_fields[0]->GetBndCondExpansions()[bcRegion]->
GetExp(e)->GetCoords(x0, x1, x2);
// Flow parameters
NekDouble c, k, phi, r, J, VV, pp, sint, P, ss;
NekDouble J11, J12, J21, J22, det;
NekDouble Fx, Fy;
NekDouble xi, yi;
NekDouble dtheta;
NekDouble par1;
NekDouble theta = M_PI / 4.0;
NekDouble kExt = 0.7;
NekDouble V = kExt * sin(theta);
NekDouble toll = 1.0e-8;
NekDouble errV = 1.0;
NekDouble errTheta = 1.0;
NekDouble gamma = m_gamma;
NekDouble gamma_1_2 = (gamma - 1.0) / 2.0;
// Loop on all the points of that edge
for (int j = 0; j < npoints; j++)
{
while ((abs(errV) > toll) || (abs(errTheta) > toll))
{
VV = V * V;
sint = sin(theta);
c = sqrt(1.0 - gamma_1_2 * VV);
k = V / sint;
phi = 1.0 / k;
pp = phi * phi;
J = 1.0 / c + 1.0 / (3.0 * c * c * c) +
1.0 / (5.0 * c * c * c * c * c) -
0.5 * log((1.0 + c) / (1.0 - c));
r = pow(c, 1.0 / gamma_1_2);
xi = 1.0 / (2.0 * r) * (1.0 / VV - 2.0 * pp) + J / 2.0;
yi = phi / (r * V) * sqrt(1.0 - VV * pp);
par1 = 25.0 - 5.0 * VV;
ss = sint * sint;
Fx = xi - x0[j];
Fy = yi - x1[j];
J11 = 39062.5 / pow(par1, 3.5) *
(1.0 / VV - 2.0 / VV * ss) * V + 1562.5 /
pow(par1, 2.5) * (-2.0 / (VV * V) + 4.0 /
(VV * V) * ss) + 12.5 / pow(par1, 1.5) * V +
312.5 / pow(par1, 2.5) * V + 7812.5 /
pow(par1, 3.5) * V - 0.25 *
(-1.0 / pow(par1, 0.5) * V / (1.0 - 0.2 *
pow(par1, 0.5)) - (1.0 + 0.2 * pow(par1, 0.5)) /
pow((1.0 - 0.2 * pow(par1, 0.5)), 2.0) /
pow(par1, 0.5) * V) / (1.0 + 0.2 * pow(par1, 0.5)) *
(1.0 - 0.2 * pow(par1, 0.5));
J12 = -6250.0 / pow(par1, 2.5) / VV * sint * cos(theta);
J21 = -6250.0 / (VV * V) * sint / pow(par1, 2.5) *
pow((1.0 - ss), 0.5) + 78125.0 / V * sint /
pow(par1, 3.5) * pow((1.0 - ss), 0.5);
// the matrix is singular when theta = pi/2
if (abs(x1[j]) < toll && abs(cos(theta)) < toll)
{
J22 = -39062.5 / pow(par1, 3.5) / V + 3125 /
pow(par1, 2.5) / (VV * V) + 12.5 /
pow(par1, 1.5) * V + 312.5 / pow(par1, 2.5) *
V + 7812.5 / pow(par1, 3.5) * V - 0.25 *
(-1.0 / pow(par1, 0.5) * V / (1.0 - 0.2 *
pow(par1, 0.5)) - (1.0 + 0.2 * pow(par1, 0.5)) /
pow((1.0 - 0.2 * pow(par1, 0.5)), 2.0) /
pow(par1, 0.5) * V) / (1.0 + 0.2 *
pow(par1, 0.5)) * (1.0 - 0.2 * pow(par1, 0.5));
// dV = -dV/dx * Fx
dV = -1.0 / J22 * Fx;
dtheta = 0.0;
theta = M_PI / 2.0;
}
else
{
J22 = 3125.0 / VV * cos(theta) / pow(par1, 2.5) *
pow((1.0 - ss), 0.5) - 3125.0 / VV * ss /
pow(par1, 2.5) / pow((1.0 - ss), 0.5) *
cos(theta);
det = -1.0 / (J11 * J22 - J12 * J21);
// [dV dtheta]' = -[invJ]*[Fx Fy]'
dV = det * ( J22 * Fx - J12 * Fy);
dtheta = det * (-J21 * Fx + J11 * Fy);
}
V = V + dV;
theta = theta + dtheta;
errV = abs(dV);
errTheta = abs(dtheta);
}
c = sqrt(1.0 - gamma_1_2 * VV);
int kk = id2 + j;
NekDouble timeramp = 200.0;
std::string restartstr = "RESTART";
if (time<timeramp &&
!(m_session->DefinesFunction("InitialConditions") &&
m_session->GetFunctionType("InitialConditions", 0) ==
{
Fwd[0][kk] = pow(c, 1.0 / gamma_1_2) *
exp(-1.0 + time /timeramp);
Fwd[1][kk] = Fwd[0][kk] * V * cos(theta) *
exp(-1 + time / timeramp);
Fwd[2][kk] = Fwd[0][kk] * V * sin(theta) *
exp(-1 + time / timeramp);
}
else
{
Fwd[0][kk] = pow(c, 1.0 / gamma_1_2);
Fwd[1][kk] = Fwd[0][kk] * V * cos(theta);
Fwd[2][kk] = Fwd[0][kk] * V * sin(theta);
}
P = (c * c) * Fwd[0][kk] / gamma;
Fwd[3][kk] = P / (gamma - 1.0) + 0.5 *
(Fwd[1][kk] * Fwd[1][kk] / Fwd[0][kk] +
Fwd[2][kk] * Fwd[2][kk] / Fwd[0][kk]);
errV = 1.0;
errTheta = 1.0;
theta = M_PI / 4.0;
V = kExt * sin(theta);
}
for (int i = 0; i < nvariables; ++i)
{
Vmath::Vcopy(npoints, &Fwd[i][id2], 1,
&(m_fields[i]->GetBndCondExpansions()[bcRegion]->
UpdatePhys())[id1],1);
}
}
}
void Nektar::EulerCFE::SetInitialIsentropicVortex ( NekDouble  initialtime)
private

Set the initial condition for the isentropic vortex problem.

Definition at line 401 of file EulerCFE.cpp.

References EvaluateIsentropicVortex(), Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_spacedim, and Vmath::Vcopy().

Referenced by v_SetInitialConditions().

{
int nTotQuadPoints = GetTotPoints();
Array<OneD, NekDouble> x(nTotQuadPoints);
Array<OneD, NekDouble> y(nTotQuadPoints);
Array<OneD, NekDouble> z(nTotQuadPoints);
Array<OneD, Array<OneD, NekDouble> > u(m_spacedim+2);
m_fields[0]->GetCoords(x, y, z);
for (int i = 0; i < m_spacedim + 2; ++i)
{
u[i] = Array<OneD, NekDouble>(nTotQuadPoints);
}
EvaluateIsentropicVortex(x, y, z, u, initialtime);
// Forward transform to fill the coefficient space
for(int i = 0; i < m_fields.num_elements(); ++i)
{
Vmath::Vcopy(nTotQuadPoints, u[i], 1, m_fields[i]->UpdatePhys(), 1);
m_fields[i]->SetPhysState(true);
m_fields[i]->FwdTrans(m_fields[i]->GetPhys(),
m_fields[i]->UpdateCoeffs());
}
}
void Nektar::EulerCFE::SetInitialRinglebFlow ( void  )
private

Set the initial condition for the Ringleb flow problem.

Definition at line 632 of file EulerCFE.cpp.

References Nektar::SolverUtils::EquationSystem::m_fields, Nektar::CompressibleFlowSystem::m_gamma, Nektar::SolverUtils::EquationSystem::m_sessionName, and Nektar::SolverUtils::EquationSystem::WriteFld().

{
// Get number of different boundaries in the input file
int nbnd = m_fields[0]->GetBndConditions().num_elements();
// Loop on all the edges of the input file
for(int bcRegion=0; bcRegion < nbnd; ++bcRegion)
{
int npoints = m_fields[0]->
GetBndCondExpansions()[bcRegion]->GetNpoints();
Array<OneD,NekDouble> x0(npoints, 0.0);
Array<OneD,NekDouble> x1(npoints, 0.0);
Array<OneD,NekDouble> x2(npoints, 0.0);
Array<OneD, NekDouble> rho(npoints, 0.0);
Array<OneD, NekDouble> rhou(npoints, 0.0);
Array<OneD, NekDouble> rhov(npoints, 0.0);
Array<OneD, NekDouble> E(npoints, 0.0);
m_fields[0]->GetBndCondExpansions()[bcRegion]->
GetCoords(x0, x1, x2);
// Flow parameters
NekDouble c, k, phi, r, J, VV, pp, sint, P, ss;
NekDouble J11, J12, J21, J22, det;
NekDouble Fx, Fy;
NekDouble xi, yi;
NekDouble dtheta;
NekDouble par1;
NekDouble theta = M_PI / 4.0;
NekDouble kExt = 0.7;
NekDouble V = kExt * sin(theta);
NekDouble toll = 1.0e-8;
NekDouble errV = 1.0;
NekDouble errTheta = 1.0;
NekDouble gamma = m_gamma;
NekDouble gamma_1_2 = (gamma - 1.0) / 2.0;
// Loop on all the points of that edge
for (int j = 0; j < npoints; j++)
{
while ((abs(errV) > toll) || (abs(errTheta) > toll))
{
VV = V * V;
sint = sin(theta);
c = sqrt(1.0 - gamma_1_2 * VV);
k = V / sint;
phi = 1.0 / k;
pp = phi * phi;
J = 1.0 / c + 1.0 / (3.0 * c * c * c) +
1.0 / (5.0 * c * c * c * c * c) -
0.5 * log((1.0 + c) / (1.0 - c));
r = pow(c, 1.0 / gamma_1_2);
xi = 1.0 / (2.0 * r) * (1.0 / VV - 2.0 * pp) + J / 2.0;
yi = phi / (r * V) * sqrt(1.0 - VV * pp);
par1 = 25.0 - 5.0 * VV;
ss = sint * sint;
Fx = xi - x0[j];
Fy = yi - x1[j];
J11 = 39062.5 / pow(par1, 3.5) *
(1.0 / VV - 2.0 / VV * ss) * V +
1562.5 / pow(par1, 2.5) * (-2.0 /
(VV * V) + 4.0 / (VV * V) * ss) +
12.5 / pow(par1, 1.5) * V +
312.5 / pow(par1, 2.5) * V +
7812.5 / pow(par1, 3.5) * V -
0.25 * (-1.0 / pow(par1, 0.5) * V /
(1.0 - 0.2 * pow(par1, 0.5)) - (1.0 + 0.2 *
pow(par1, 0.5)) / pow((1.0 - 0.2 *
pow(par1, 0.5)), 2.0) /
pow(par1, 0.5) * V) /
(1.0 + 0.2 * pow(par1, 0.5)) *
(1.0 - 0.2 * pow(par1, 0.5));
J12 = -6250.0 / pow(par1, 2.5) / VV * sint * cos(theta);
J21 = -6250.0 / (VV * V) * sint / pow(par1, 2.5) *
pow((1.0 - ss), 0.5) + 78125.0 / V * sint /
pow(par1, 3.5) * pow((1.0 - ss), 0.5);
// the matrix is singular when theta = pi/2
if (abs(x1[j]) < toll && abs(cos(theta)) < toll)
{
J22 = -39062.5 / pow(par1, 3.5) / V +
3125 / pow(par1, 2.5) / (VV * V) + 12.5 /
pow(par1, 1.5) * V + 312.5 / pow(par1, 2.5) * V +
7812.5 / pow(par1, 3.5) * V -
0.25 * (-1.0 / pow(par1, 0.5) * V /
(1.0 - 0.2 * pow(par1, 0.5)) -
(1.0 + 0.2 * pow(par1, 0.5)) /
pow((1.0 - 0.2* pow(par1, 0.5)), 2.0) /
pow(par1, 0.5) *V) /
(1.0 + 0.2 * pow(par1, 0.5)) *
(1.0 - 0.2 * pow(par1, 0.5));
// dV = -dV/dx * Fx
dV = -1.0 / J22 * Fx;
dtheta = 0.0;
theta = M_PI / 2.0;
}
else
{
J22 = 3125.0 / VV * cos(theta) / pow(par1, 2.5) *
pow((1.0 - ss), 0.5) - 3125.0 / VV * ss /
pow(par1, 2.5) / pow((1.0 - ss), 0.5) *
cos(theta);
det = -1.0 / (J11 * J22 - J12 * J21);
// [dV dtheta]' = -[invJ]*[Fx Fy]'
dV = det * ( J22 * Fx - J12 * Fy);
dtheta = det * (-J21 * Fx + J11 * Fy);
}
V = V + dV;
theta = theta + dtheta;
errV = abs(dV);
errTheta = abs(dtheta);
}
c = sqrt(1.0 - gamma_1_2 * VV);
rho[j] = pow(c, 1.0 / gamma_1_2) * exp(-1.0);
rhou[j] = rho[j] * V * cos(theta) * exp(-1.0);
rhov[j] = rho[j] * V * sin(theta) * exp(-1.0);
P = (c * c) * rho[j] / gamma;
E[j] = P / (gamma - 1.0) + 0.5 *
(rhou[j] * rhou[j] /
rho[j] + rhov[j] * rhov[j] / rho[j]);
errV = 1.0;
errTheta = 1.0;
theta = M_PI / 4.0;
V = kExt * sin(theta);
}
// Fill the physical space
m_fields[0]->GetBndCondExpansions()[bcRegion]->SetPhys(rho);
m_fields[1]->GetBndCondExpansions()[bcRegion]->SetPhys(rhou);
m_fields[2]->GetBndCondExpansions()[bcRegion]->SetPhys(rhov);
m_fields[3]->GetBndCondExpansions()[bcRegion]->SetPhys(E);
// Forward transform to fill the coefficients space
for(int i = 0; i < m_fields.num_elements(); ++i)
{
m_fields[i]->GetBndCondExpansions()[bcRegion]->
FwdTrans_BndConstrained(
m_fields[i]->GetBndCondExpansions()[bcRegion]->
GetPhys(),
m_fields[i]->GetBndCondExpansions()[bcRegion]->
UpdateCoeffs());
}
}
// Calculation of the initial internal values as a weighted average
// over the distance from the boundaries
int nq = m_fields[0]->GetNpoints();
Array<OneD,NekDouble> x0(nq);
Array<OneD,NekDouble> x1(nq);
Array<OneD,NekDouble> x2(nq);
// Get the coordinates
// (assuming all fields have the same discretisation)
m_fields[0]->GetCoords(x0, x1, x2);
for(int j = 0; j < nq; j++)
{
NekDouble Dist = 0.0;
NekDouble rho = 0.0;
NekDouble rhou = 0.0;
NekDouble rhov = 0.0;
NekDouble E = 0.0;
NekDouble SumDist = 0.0;
// Calculation of all the distances
// Loop on all the edges of the input file
for (int bcRegion = 0; bcRegion < nbnd; ++bcRegion)
{
// Get quadrature points on the edge
int npoints = m_fields[0]->
GetBndCondExpansions()[bcRegion]->GetNpoints();
Array<OneD,NekDouble> xb0(npoints, 0.0);
Array<OneD,NekDouble> xb1(npoints, 0.0);
Array<OneD,NekDouble> xb2(npoints, 0.0);
m_fields[0]->GetBndCondExpansions()[bcRegion]->
GetCoords(xb0, xb1, xb2);
for (int k = 0; k < npoints; k++)
{
Dist = sqrt((xb0[k] - x0[j]) * (xb0[k] - x0[j]) +
(xb1[k] - x1[j]) * (xb1[k] - x1[j]));
SumDist += Dist;
rho += Dist * (m_fields[0]->
GetBndCondExpansions()[bcRegion]->GetPhys())[k];
rhou += Dist * (m_fields[1]->
GetBndCondExpansions()[bcRegion]->GetPhys())[k];
rhov += Dist * (m_fields[2]->
GetBndCondExpansions()[bcRegion]->GetPhys())[k];
E += Dist * (m_fields[3]->
GetBndCondExpansions()[bcRegion]->GetPhys())[k];
}
}
rho = rho / SumDist;
rhou = rhou / SumDist;
rhov = rhov / SumDist;
E = E / SumDist;
(m_fields[0]->UpdatePhys())[j] = rho;
(m_fields[1]->UpdatePhys())[j] = rhou;
(m_fields[2]->UpdatePhys())[j] = rhov;
(m_fields[3]->UpdatePhys())[j] = E;
}
for (int i = 0 ; i < m_fields.num_elements(); i++)
{
m_fields[i]->SetPhysState(true);
m_fields[i]->FwdTrans(m_fields[i]->GetPhys(),
m_fields[i]->UpdateCoeffs());
}
// Dump initial conditions to file
std::string outname = m_sessionName + "_initialRingleb.chk";
WriteFld(outname);
}
void Nektar::EulerCFE::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 307 of file EulerCFE.cpp.

References Nektar::eIsentropicVortex, Nektar::eRinglebFlow, GetExactIsentropicVortex(), GetExactRinglebFlow(), and m_problemType.

{
switch(m_problemType)
{
{
GetExactIsentropicVortex(field, outfield, time);
break;
}
{
GetExactRinglebFlow(field, outfield);
break;
}
default:
{
break;
}
}
}
void Nektar::EulerCFE::v_GenerateSummary ( SolverUtils::SummaryList s)
protectedvirtual

Print a summary of time stepping parameters.

Print out a summary with some relevant information.

Reimplemented from Nektar::CompressibleFlowSystem.

Definition at line 102 of file EulerCFE.cpp.

References Nektar::SolverUtils::AddSummaryItem(), m_problemType, and Nektar::ProblemTypeMap.

void Nektar::EulerCFE::v_InitObject ( )
protectedvirtual
void Nektar::EulerCFE::v_SetInitialConditions ( NekDouble  initialtime = 0.0,
bool  dumpInitialConditions = true,
const int  domain = 0 
)
protectedvirtual

Set the initial conditions.

Reimplemented from Nektar::CompressibleFlowSystem.

Definition at line 111 of file EulerCFE.cpp.

References Nektar::SolverUtils::EquationSystem::Checkpoint_Output(), Nektar::eIsentropicVortex, Vmath::FillWhiteNoise(), Nektar::SolverUtils::EquationSystem::m_comm, Nektar::SolverUtils::EquationSystem::m_fields, m_problemType, Nektar::SolverUtils::EquationSystem::m_session, SetInitialIsentropicVortex(), and Vmath::Vadd().

{
switch (m_problemType)
{
{
break;
}
default:
{
break;
}
}
//insert white noise in initial condition
NekDouble Noise;
int phystot = m_fields[0]->GetTotPoints();
Array<OneD, NekDouble> noise(phystot);
m_session->LoadParameter("Noise", Noise,0.0);
int m_nConvectiveFields = m_fields.num_elements();
if(Noise > 0.0)
{
for(int i = 0; i < m_nConvectiveFields; i++)
{
Vmath::FillWhiteNoise(phystot,Noise,noise,1,m_comm->GetColumnComm()->GetRank()+1);
Vmath::Vadd(phystot,m_fields[i]->GetPhys(),1,noise,1,m_fields[i]->UpdatePhys(),1);
m_fields[i]->FwdTrans_IterPerExp(m_fields[i]->GetPhys(),m_fields[i]->UpdateCoeffs());
}
}
if (dumpInitialConditions)
{
// Dump initial conditions to file
}
}

Friends And Related Function Documentation

friend class MemoryManager< EulerCFE >
friend

Definition at line 61 of file EulerCFE.h.

Member Data Documentation

string Nektar::EulerCFE::className
static
Initial value:
"EulerCFE", EulerCFE::create,
"Euler equations in conservative variables without "
"artificial diffusion.")

Name of class.

Definition at line 72 of file EulerCFE.h.

ProblemType Nektar::EulerCFE::m_problemType