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

#include <CompressibleFlowSystem.h>

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

Public Member Functions

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...
 
- 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 std::string &s1, const std::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...
 

Protected Member Functions

 CompressibleFlowSystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 
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 InitializeSteadyState ()
 
void SetBoundaryConditions (Array< OneD, Array< OneD, NekDouble > > &physarray, NekDouble time)
 
void GetStdVelocity (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &stdV)
 Compute the advection velocity in the standard space for each element of the expansion. More...
 
virtual bool v_PostIntegrate (int step)
 Perform post-integration checks, presently just to check steady state behaviour. More...
 
bool CalcSteadyState (bool output)
 Calculate whether the system has reached a steady state by observing residuals to a user-defined tolerance. 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)
 
- 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_GenerateSummary (SummaryList &s)
 Print a summary of time stepping parameters. 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 bool v_PreIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT bool v_SteadyStateCheck (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)
 Initialises EquationSystem class members. More...
 
int nocase_cmp (const std::string &s1, const std::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_EvaluateExactSolution (unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 
SOLVER_UTILS_EXPORT void EvaluateFunctionExp (std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, const NekDouble &pTime=0.0, const int domain=0)
 
SOLVER_UTILS_EXPORT void EvaluateFunctionFld (std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, const NekDouble &pTime=0.0, const int domain=0)
 
SOLVER_UTILS_EXPORT void EvaluateFunctionPts (std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, const NekDouble &pTime=0.0, const int domain=0)
 
SOLVER_UTILS_EXPORT void LoadPts (std::string funcFilename, std::string filename, Nektar::LibUtilities::PtsFieldSharedPtr &outPts)
 
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)
 

Protected Attributes

SolverUtils::AdvectionSharedPtr m_advection
 
SolverUtils::DiffusionSharedPtr m_diffusion
 
ArtificialDiffusionSharedPtr m_artificialDiffusion
 
Array< OneD, Array< OneD,
NekDouble > > 
m_vecLocs
 
NekDouble m_gamma
 
NekDouble m_pInf
 
NekDouble m_rhoInf
 
NekDouble m_UInf
 
std::string m_ViscosityType
 
std::string m_shockCaptureType
 
NekDouble m_mu
 
NekDouble m_thermalConductivity
 
NekDouble m_Cp
 
NekDouble m_Prandtl
 
VariableConverterSharedPtr m_varConv
 
std::vector< CFSBndCondSharedPtrm_bndConds
 
std::ofstream m_errFile
 
NekDouble m_steadyStateTol
 
std::vector
< SolverUtils::ForcingSharedPtr
m_forcing
 
Array< OneD, Array< OneD,
NekDouble > > 
m_un
 
- Protected Attributes inherited from Nektar::SolverUtils::UnsteadySystem
int m_infosteps
 Number of time steps between outputting status 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...
 
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...
 
std::map< std::string,
FieldUtils::Interpolator
m_interpolators
 Map of interpolator objects. More...
 
std::map< std::string,
std::pair< std::string,
LibUtilities::PtsFieldSharedPtr > > 
m_loadedPtsFields
 pts fields we already read from disk: {funcFilename: (filename, ptsfield)} More...
 
std::map< std::string,
std::pair< std::string,
loadedFldField > > 
m_loadedFldFields
 
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...
 
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...
 

Friends

class MemoryManager< CompressibleFlowSystem >
 

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 53 of file CompressibleFlowSystem.h.

Constructor & Destructor Documentation

Nektar::CompressibleFlowSystem::~CompressibleFlowSystem ( )
virtual

Destructor for CompressibleFlowSystem class.

Definition at line 129 of file CompressibleFlowSystem.cpp.

130  {
131 
132  }
Nektar::CompressibleFlowSystem::CompressibleFlowSystem ( const LibUtilities::SessionReaderSharedPtr pSession)
protected

Definition at line 45 of file CompressibleFlowSystem.cpp.

47  : UnsteadySystem(pSession)
48  {
49  }
SOLVER_UTILS_EXPORT UnsteadySystem(const LibUtilities::SessionReaderSharedPtr &pSession)
Initialises UnsteadySystem class members.

Member Function Documentation

bool Nektar::CompressibleFlowSystem::CalcSteadyState ( bool  output)
protected

Calculate whether the system has reached a steady state by observing residuals to a user-defined tolerance.

Definition at line 579 of file CompressibleFlowSystem.cpp.

References Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::SolverUtils::EquationSystem::m_comm, m_errFile, Nektar::SolverUtils::EquationSystem::m_fields, m_gamma, m_pInf, m_rhoInf, Nektar::SolverUtils::EquationSystem::m_session, m_steadyStateTol, Nektar::SolverUtils::EquationSystem::m_time, m_UInf, m_un, Nektar::LibUtilities::ReduceSum, Vmath::Vcopy(), Vmath::Vmax(), Vmath::Vmul(), Vmath::Vsub(), and Vmath::Vsum().

Referenced by v_PostIntegrate().

580  {
581  const int nPoints = GetTotPoints();
582  const int nFields = m_fields.num_elements();
583 
584  // Holds L2 errors.
585  Array<OneD, NekDouble> L2 (nFields);
586  Array<OneD, NekDouble> residual(nFields);
587 
588  for (int i = 0; i < nFields; ++i)
589  {
590  Array<OneD, NekDouble> diff(nPoints);
591 
592  Vmath::Vsub(nPoints, m_fields[i]->GetPhys(), 1, m_un[i], 1, diff, 1);
593  Vmath::Vmul(nPoints, diff, 1, diff, 1, diff, 1);
594  residual[i] = Vmath::Vsum(nPoints, diff, 1);
595  }
596 
597  m_comm->AllReduce(residual, LibUtilities::ReduceSum);
598 
599  // L2 error
600  L2[0] = sqrt(residual[0]) / m_rhoInf;
601 
602  for (int i = 1; i < nFields-1; ++i)
603  {
604  L2[i] = sqrt(residual[i]) / m_UInf / m_rhoInf;
605  }
606 
607  NekDouble Einf = m_pInf / (m_gamma-1.0) + 0.5 * m_rhoInf * m_UInf;
608  L2[nFields-1] = sqrt(residual[nFields-1]) / Einf;
609 
610  if (m_comm->GetRank() == 0 && output)
611  {
612  // Output time
613  m_errFile << setprecision(8) << setw(17) << scientific << m_time;
614 
615  // Output residuals
616  for (int i = 0; i < nFields; ++i)
617  {
618  m_errFile << setprecision(11) << setw(22) << scientific
619  << L2[i];
620  }
621 
622  m_errFile << endl;
623  }
624 
625  // Calculate maximum L2 error
626  NekDouble maxL2 = Vmath::Vmax(nFields, L2, 1);
627 
628  if (m_session->DefinesCmdLineArgument("verbose") &&
629  m_comm->GetRank() == 0 && output)
630  {
631  cout << "-- Maximum L^2 residual: " << maxL2 << endl;
632  }
633 
634  if (maxL2 <= m_steadyStateTol)
635  {
636  return true;
637  }
638 
639  for (int i = 0; i < m_fields.num_elements(); ++i)
640  {
641  Vmath::Vcopy(nPoints, m_fields[i]->GetPhys(), 1, m_un[i], 1);
642  }
643 
644  return false;
645  }
NekDouble m_time
Current time of simulation.
T Vmax(int n, const T *x, const int incx)
Return the maximum element in x – called vmax to avoid conflict with max.
Definition: Vmath.cpp:779
Array< OneD, Array< OneD, NekDouble > > m_un
LibUtilities::CommSharedPtr m_comm
Communicator.
SOLVER_UTILS_EXPORT int GetTotPoints()
double NekDouble
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:343
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.cpp:737
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
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:183
void Nektar::CompressibleFlowSystem::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 
)
protected

Compute the advection terms for the right-hand side.

Definition at line 343 of file CompressibleFlowSystem.cpp.

References m_advection, Nektar::SolverUtils::EquationSystem::m_fields, and Nektar::SolverUtils::EquationSystem::m_spacedim.

Referenced by DoOdeRhs().

349  {
350  int nvariables = inarray.num_elements();
351  Array<OneD, Array<OneD, NekDouble> > advVel(m_spacedim);
352 
353  m_advection->Advect(nvariables, m_fields, advVel, inarray,
354  outarray, time, pFwd, pBwd);
355  }
int m_spacedim
Spatial dimension (>= expansion dim).
SolverUtils::AdvectionSharedPtr m_advection
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
void Nektar::CompressibleFlowSystem::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 
)
protected

Add the diffusions terms to the right-hand side.

Definition at line 360 of file CompressibleFlowSystem.cpp.

References m_artificialDiffusion, m_shockCaptureType, and v_DoDiffusion().

Referenced by DoOdeRhs().

365  {
366  v_DoDiffusion(inarray, outarray, pFwd, pBwd);
367 
368  if (m_shockCaptureType != "Off")
369  {
370  m_artificialDiffusion->DoArtificialDiffusion(inarray, outarray);
371  }
372  }
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)
ArtificialDiffusionSharedPtr m_artificialDiffusion
void Nektar::CompressibleFlowSystem::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 305 of file CompressibleFlowSystem.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().

309  {
310  int i;
311  int nvariables = inarray.num_elements();
312 
313  switch(m_projectionType)
314  {
316  {
317  // Just copy over array
318  int npoints = GetNpoints();
319 
320  for(i = 0; i < nvariables; ++i)
321  {
322  Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
323  }
324  SetBoundaryConditions(outarray, time);
325  break;
326  }
329  {
330  ASSERTL0(false, "No Continuous Galerkin for full compressible "
331  "Navier-Stokes equations");
332  break;
333  }
334  default:
335  ASSERTL0(false, "Unknown projection scheme");
336  break;
337  }
338  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
void SetBoundaryConditions(Array< OneD, Array< OneD, NekDouble > > &physarray, NekDouble time)
SOLVER_UTILS_EXPORT int GetNpoints()
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
void Nektar::CompressibleFlowSystem::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 252 of file CompressibleFlowSystem.cpp.

References DoAdvection(), DoDiffusion(), Nektar::SolverUtils::EquationSystem::eHomogeneous1D, Nektar::SolverUtils::EquationSystem::GetNpoints(), Nektar::SolverUtils::EquationSystem::GetTraceTotPoints(), Nektar::SolverUtils::EquationSystem::m_fields, m_forcing, Nektar::SolverUtils::EquationSystem::m_HomogeneousType, Vmath::Neg(), and Nektar::NullNekDoubleArrayofArray.

Referenced by v_InitObject().

256  {
257  int i;
258  int nvariables = inarray.num_elements();
259  int npoints = GetNpoints();
260  int nTracePts = GetTraceTotPoints();
261 
262  // Store forwards/backwards space along trace space
263  Array<OneD, Array<OneD, NekDouble> > Fwd (nvariables);
264  Array<OneD, Array<OneD, NekDouble> > Bwd (nvariables);
265 
267  {
270  }
271  else
272  {
273  for(i = 0; i < nvariables; ++i)
274  {
275  Fwd[i] = Array<OneD, NekDouble>(nTracePts, 0.0);
276  Bwd[i] = Array<OneD, NekDouble>(nTracePts, 0.0);
277  m_fields[i]->GetFwdBwdTracePhys(inarray[i], Fwd[i], Bwd[i]);
278  }
279  }
280 
281  // Calculate advection
282  DoAdvection(inarray, outarray, time, Fwd, Bwd);
283 
284  // Negate results
285  for (i = 0; i < nvariables; ++i)
286  {
287  Vmath::Neg(npoints, outarray[i], 1);
288  }
289 
290  // Add diffusion terms
291  DoDiffusion(inarray, outarray, Fwd, Bwd);
292 
293  // Add forcing terms
294  std::vector<SolverUtils::ForcingSharedPtr>::const_iterator x;
295  for (x = m_forcing.begin(); x != m_forcing.end(); ++x)
296  {
297  (*x)->Apply(m_fields, inarray, outarray, time);
298  }
299  }
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.
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:396
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.
SOLVER_UTILS_EXPORT int GetTraceTotPoints()
SOLVER_UTILS_EXPORT int GetNpoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayofArray
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
enum HomogeneousType m_HomogeneousType
void Nektar::CompressibleFlowSystem::GetFluxVector ( const Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &  flux 
)
protected

Return the flux vector for the compressible Euler equations.

Parameters
physfieldFields.
fluxResulting flux.

Definition at line 405 of file CompressibleFlowSystem.cpp.

References Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_spacedim, m_varConv, Vmath::Vadd(), Vmath::Vcopy(), Vmath::Vmul(), and Vmath::Zero().

Referenced by InitAdvection().

408  {
409  int i, j;
410  int nq = physfield[0].num_elements();
411  int nVariables = m_fields.num_elements();
412 
413  Array<OneD, NekDouble> pressure(nq);
414  Array<OneD, Array<OneD, NekDouble> > velocity(m_spacedim);
415 
416  // Flux vector for the rho equation
417  for (i = 0; i < m_spacedim; ++i)
418  {
419  velocity[i] = Array<OneD, NekDouble>(nq);
420  Vmath::Vcopy(nq, physfield[i+1], 1, flux[0][i], 1);
421  }
422 
423  m_varConv->GetVelocityVector(physfield, velocity);
424  m_varConv->GetPressure(physfield, velocity, pressure);
425 
426  // Flux vector for the velocity fields
427  for (i = 0; i < m_spacedim; ++i)
428  {
429  for (j = 0; j < m_spacedim; ++j)
430  {
431  Vmath::Vmul(nq, velocity[j], 1, physfield[i+1], 1,
432  flux[i+1][j], 1);
433  }
434 
435  // Add pressure to appropriate field
436  Vmath::Vadd(nq, flux[i+1][i], 1, pressure, 1, flux[i+1][i], 1);
437  }
438 
439  // Flux vector for energy.
440  Vmath::Vadd(nq, physfield[m_spacedim+1], 1, pressure, 1,
441  pressure, 1);
442 
443  for (j = 0; j < m_spacedim; ++j)
444  {
445  Vmath::Vmul(nq, velocity[j], 1, pressure, 1,
446  flux[m_spacedim+1][j], 1);
447  }
448 
449  // For the smooth viscosity model
450  if (nVariables == m_spacedim+3)
451  {
452  // Add a zero row for the advective fluxes
453  for (j = 0; j < m_spacedim; ++j)
454  {
455  Vmath::Zero(nq, flux[m_spacedim+2][j], 1);
456  }
457  }
458  }
VariableConverterSharedPtr m_varConv
int m_spacedim
Spatial dimension (>= expansion dim).
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:373
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
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:299
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:183
void Nektar::CompressibleFlowSystem::GetFluxVectorDeAlias ( const Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &  flux 
)
protected

Return the flux vector for the compressible Euler equations by using the de-aliasing technique.

Parameters
physfieldFields.
fluxResulting flux.

Definition at line 467 of file CompressibleFlowSystem.cpp.

References Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_spacedim, m_varConv, Vmath::Vadd(), and Vmath::Vmul().

Referenced by InitAdvection().

470  {
471  int i, j;
472  int nq = physfield[0].num_elements();
473  int nVariables = m_fields.num_elements();
474 
475  // Factor to rescale 1d points in dealiasing
476  NekDouble OneDptscale = 2;
477  nq = m_fields[0]->Get1DScaledTotPoints(OneDptscale);
478 
479  Array<OneD, NekDouble> pressure(nq);
480  Array<OneD, Array<OneD, NekDouble> > velocity(m_spacedim);
481 
482  Array<OneD, Array<OneD, NekDouble> > physfield_interp(nVariables);
483  Array<OneD, Array<OneD, Array<OneD, NekDouble> > > flux_interp(
484  nVariables);
485 
486  for (i = 0; i < nVariables; ++ i)
487  {
488  physfield_interp[i] = Array<OneD, NekDouble>(nq);
489  flux_interp[i] = Array<OneD, Array<OneD, NekDouble> >(m_spacedim);
490  m_fields[0]->PhysInterp1DScaled(
491  OneDptscale, physfield[i], physfield_interp[i]);
492 
493  for (j = 0; j < m_spacedim; ++j)
494  {
495  flux_interp[i][j] = Array<OneD, NekDouble>(nq);
496  }
497  }
498 
499  // Flux vector for the rho equation
500  for (i = 0; i < m_spacedim; ++i)
501  {
502  velocity[i] = Array<OneD, NekDouble>(nq);
503 
504  // Galerkin project solution back to original space
505  m_fields[0]->PhysGalerkinProjection1DScaled(
506  OneDptscale, physfield_interp[i+1], flux[0][i]);
507  }
508 
509  m_varConv->GetVelocityVector(physfield_interp, velocity);
510  m_varConv->GetPressure (physfield_interp, velocity, pressure);
511 
512  // Evaluation of flux vector for the velocity fields
513  for (i = 0; i < m_spacedim; ++i)
514  {
515  for (j = 0; j < m_spacedim; ++j)
516  {
517  Vmath::Vmul(nq, velocity[j], 1, physfield_interp[i+1], 1,
518  flux_interp[i+1][j], 1);
519  }
520 
521  // Add pressure to appropriate field
522  Vmath::Vadd(nq, flux_interp[i+1][i], 1, pressure,1,
523  flux_interp[i+1][i], 1);
524  }
525 
526  // Galerkin project solution back to origianl space
527  for (i = 0; i < m_spacedim; ++i)
528  {
529  for (j = 0; j < m_spacedim; ++j)
530  {
531  m_fields[0]->PhysGalerkinProjection1DScaled(
532  OneDptscale, flux_interp[i+1][j], flux[i+1][j]);
533  }
534  }
535 
536  // Evaluation of flux vector for energy
537  Vmath::Vadd(nq, physfield_interp[m_spacedim+1], 1, pressure, 1,
538  pressure, 1);
539 
540  for (j = 0; j < m_spacedim; ++j)
541  {
542  Vmath::Vmul(nq, velocity[j], 1, pressure, 1,
543  flux_interp[m_spacedim+1][j], 1);
544 
545  // Galerkin project solution back to origianl space
546  m_fields[0]->PhysGalerkinProjection1DScaled(
547  OneDptscale,
548  flux_interp[m_spacedim+1][j],
549  flux[m_spacedim+1][j]);
550  }
551  }
VariableConverterSharedPtr m_varConv
int m_spacedim
Spatial dimension (>= expansion dim).
double NekDouble
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:299
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:183
NekDouble Nektar::CompressibleFlowSystem::GetGamma ( )
inlineprotected

Definition at line 161 of file CompressibleFlowSystem.h.

References m_gamma.

Referenced by InitAdvection().

162  {
163  return m_gamma;
164  }
const Array<OneD, const Array<OneD, NekDouble> >& Nektar::CompressibleFlowSystem::GetNormals ( )
inlineprotected

Definition at line 171 of file CompressibleFlowSystem.h.

References Nektar::SolverUtils::EquationSystem::m_traceNormals.

Referenced by InitAdvection().

172  {
173  return m_traceNormals;
174  }
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
NekDouble Nektar::CompressibleFlowSystem::GetStabilityLimit ( int  n)

Function to calculate the stability limit for DG/CG.

Set the denominator to compute the time step when a cfl control is employed. This function is no longer used but is still here for being utilised in the future.

Parameters
nOrder of expansion element by element.

Definition at line 870 of file CompressibleFlowSystem.cpp.

References ASSERTL0, Nektar::MultiRegions::eDiscontinuous, and Nektar::SolverUtils::EquationSystem::m_projectionType.

Referenced by GetStabilityLimitVector().

871  {
872  ASSERTL0(n <= 20, "Illegal modes dimension for CFL calculation "
873  "(P has to be less then 20)");
874 
875  NekDouble CFLDG[21] = { 2.0000, 6.0000, 11.8424, 19.1569,
876  27.8419, 37.8247, 49.0518, 61.4815,
877  75.0797, 89.8181, 105.6700, 122.6200,
878  140.6400, 159.7300, 179.8500, 201.0100,
879  223.1800, 246.3600, 270.5300, 295.6900,
880  321.8300}; //CFLDG 1D [0-20]
881  NekDouble CFL = 0.0;
882 
884  {
885  CFL = CFLDG[n];
886  }
887  else
888  {
889  ASSERTL0(false, "Continuous Galerkin stability coefficients "
890  "not introduced yet.");
891  }
892 
893  return CFL;
894  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
double NekDouble
Array< OneD, NekDouble > Nektar::CompressibleFlowSystem::GetStabilityLimitVector ( const Array< OneD, int > &  ExpOrder)

Function to calculate the stability limit for DG/CG (a vector of them).

Compute the vector of denominators to compute the time step when a cfl control is employed. This function is no longer used but is still here for being utilised in the future.

Parameters
ExpOrderOrder of expansion element by element.

Definition at line 903 of file CompressibleFlowSystem.cpp.

References Nektar::SolverUtils::EquationSystem::GetExpSize(), GetStabilityLimit(), and Nektar::SolverUtils::EquationSystem::m_fields.

905  {
906  int i;
907  Array<OneD,NekDouble> returnval(m_fields[0]->GetExpSize(), 0.0);
908  for (i =0; i<m_fields[0]->GetExpSize(); i++)
909  {
910  returnval[i] = GetStabilityLimit(ExpOrder[i]);
911  }
912  return returnval;
913  }
NekDouble GetStabilityLimit(int n)
Function to calculate the stability limit for DG/CG.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetExpSize()
void Nektar::CompressibleFlowSystem::GetStdVelocity ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, NekDouble > &  stdV 
)
protected

Compute the advection velocity in the standard space for each element of the expansion.

Parameters
inarrayMomentum field.
stdVStandard velocity field.

Definition at line 777 of file CompressibleFlowSystem.cpp.

References Nektar::SpatialDomains::eDeformed, Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_spacedim, m_varConv, Vmath::Smul(), Vmath::Svtvp(), Vmath::Vmul(), Vmath::Vvtvp(), and Vmath::Zero().

Referenced by v_GetTimeStep().

780  {
781  int nTotQuadPoints = GetTotPoints();
782  int n_element = m_fields[0]->GetExpSize();
783  int nBCEdgePts = 0;
784 
785  // Getting the velocity vector on the 2D normal space
786  Array<OneD, Array<OneD, NekDouble> > velocity (m_spacedim);
787  Array<OneD, Array<OneD, NekDouble> > stdVelocity(m_spacedim);
788  Array<OneD, NekDouble> pressure (nTotQuadPoints);
789  Array<OneD, NekDouble> soundspeed (nTotQuadPoints);
791 
792  // Zero output array
793  Vmath::Zero(stdV.num_elements(), stdV, 1);
794 
795  for (int i = 0; i < m_spacedim; ++i)
796  {
797  velocity [i] = Array<OneD, NekDouble>(nTotQuadPoints);
798  stdVelocity[i] = Array<OneD, NekDouble>(nTotQuadPoints, 0.0);
799  }
800 
801  m_varConv->GetVelocityVector(inarray, velocity);
802  m_varConv->GetPressure (inarray, velocity, pressure);
803  m_varConv->GetSoundSpeed (inarray, pressure, soundspeed);
804 
805  for(int el = 0; el < n_element; ++el)
806  {
807  ptsKeys = m_fields[0]->GetExp(el)->GetPointsKeys();
808 
809  // Possible bug: not multiply by jacobian??
810  const SpatialDomains::GeomFactorsSharedPtr metricInfo =
811  m_fields[0]->GetExp(el)->GetGeom()->GetMetricInfo();
812  const Array<TwoD, const NekDouble> &gmat =
813  m_fields[0]->GetExp(el)->GetGeom()->GetMetricInfo()
814  ->GetDerivFactors(ptsKeys);
815 
816  int nq = m_fields[0]->GetExp(el)->GetTotPoints();
817 
818  if(metricInfo->GetGtype() == SpatialDomains::eDeformed)
819  {
820  // d xi/ dx = gmat = 1/J * d x/d xi
821  for (int i = 0; i < m_spacedim; ++i)
822  {
823  Vmath::Vmul(nq, gmat[i], 1, velocity[0], 1,
824  stdVelocity[i], 1);
825  for (int j = 1; j < m_spacedim; ++j)
826  {
827  Vmath::Vvtvp(nq, gmat[m_spacedim*j+i], 1, velocity[j],
828  1, stdVelocity[i], 1, stdVelocity[i], 1);
829  }
830  }
831  }
832  else
833  {
834  for (int i = 0; i < m_spacedim; ++i)
835  {
836  Vmath::Smul(nq, gmat[i][0], velocity[0], 1,
837  stdVelocity[i], 1);
838  for (int j = 1; j < m_spacedim; ++j)
839  {
840  Vmath::Svtvp(nq, gmat[m_spacedim*j+i][0], velocity[j],
841  1, stdVelocity[i], 1, stdVelocity[i], 1);
842  }
843  }
844  }
845 
846  for (int i = 0; i < nq; ++i)
847  {
848  NekDouble pntVelocity = 0.0;
849  for (int j = 0; j < m_spacedim; ++j)
850  {
851  pntVelocity += stdVelocity[j][i]*stdVelocity[j][i];
852  }
853  pntVelocity = sqrt(pntVelocity) + soundspeed[nBCEdgePts];
854  if (pntVelocity > stdV[el])
855  {
856  stdV[el] = pntVelocity;
857  }
858  nBCEdgePts++;
859  }
860  }
861  }
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:242
VariableConverterSharedPtr m_varConv
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
svtvp (scalar times vector plus vector): z = alpha*x + y
Definition: Vmath.cpp:485
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:442
SOLVER_UTILS_EXPORT int GetTotPoints()
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:213
int m_spacedim
Spatial dimension (>= expansion dim).
double NekDouble
boost::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
Definition: GeomFactors.h:62
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:373
Geometry is curved or has non-constant factors.
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:183
const Array<OneD, const Array<OneD, NekDouble> >& Nektar::CompressibleFlowSystem::GetVecLocs ( )
inlineprotected

Definition at line 166 of file CompressibleFlowSystem.h.

References m_vecLocs.

Referenced by InitAdvection().

167  {
168  return m_vecLocs;
169  }
Array< OneD, Array< OneD, NekDouble > > m_vecLocs
void Nektar::CompressibleFlowSystem::InitAdvection ( )
protected

Create advection and diffusion objects for CFS.

Definition at line 206 of file CompressibleFlowSystem.cpp.

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, >::CreateInstance(), Nektar::MultiRegions::eDiscontinuous, Nektar::SolverUtils::GetAdvectionFactory(), GetFluxVector(), GetFluxVectorDeAlias(), GetGamma(), GetNormals(), Nektar::SolverUtils::GetRiemannSolverFactory(), GetVecLocs(), m_advection, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_projectionType, Nektar::SolverUtils::EquationSystem::m_session, and Nektar::SolverUtils::EquationSystem::m_specHP_dealiasing.

Referenced by v_InitObject().

207  {
208  // Check if projection type is correct
210  "Unsupported projection type.");
211 
212  string advName, riemName;
213  m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
214 
216  .CreateInstance(advName, advName);
217 
219  {
220  m_advection->SetFluxVector(&CompressibleFlowSystem::
221  GetFluxVectorDeAlias, this);
222  }
223  else
224  {
225  m_advection->SetFluxVector (&CompressibleFlowSystem::
226  GetFluxVector, this);
227  }
228 
229  // Setting up Riemann solver for advection operator
230  m_session->LoadSolverInfo("UpwindType", riemName, "Average");
231 
233  riemannSolver = SolverUtils::GetRiemannSolverFactory()
234  .CreateInstance(riemName);
235 
236  // Setting up parameters for advection operator Riemann solver
237  riemannSolver->SetParam (
238  "gamma", &CompressibleFlowSystem::GetGamma, this);
239  riemannSolver->SetAuxVec(
240  "vecLocs", &CompressibleFlowSystem::GetVecLocs, this);
241  riemannSolver->SetVector(
243 
244  // Concluding initialisation of advection / diffusion operators
245  m_advection->SetRiemannSolver (riemannSolver);
246  m_advection->InitObject (m_session, m_fields);
247  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
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
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
bool m_specHP_dealiasing
Flag to determine if dealisising is usde for the Spectral/hp element discretisation.
boost::shared_ptr< RiemannSolver > RiemannSolverSharedPtr
A shared pointer to an EquationSystem object.
const Array< OneD, const Array< OneD, NekDouble > > & GetVecLocs()
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...
RiemannSolverFactory & GetRiemannSolverFactory()
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:46
SolverUtils::AdvectionSharedPtr m_advection
const Array< OneD, const Array< OneD, NekDouble > > & GetNormals()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
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 Nektar::CompressibleFlowSystem::InitialiseParameters ( )
protected

Load CFS parameters from the session file.

Definition at line 137 of file CompressibleFlowSystem.cpp.

References ASSERTL0, m_Cp, m_gamma, m_mu, m_pInf, m_Prandtl, m_rhoInf, Nektar::SolverUtils::EquationSystem::m_session, m_shockCaptureType, Nektar::SolverUtils::EquationSystem::m_spacedim, m_steadyStateTol, m_thermalConductivity, m_UInf, and m_ViscosityType.

Referenced by v_InitObject().

138  {
139  NekDouble velInf, gasConstant;
140 
141  // Get gamma parameter from session file.
142  m_session->LoadParameter("Gamma", m_gamma, 1.4);
143 
144  // Get gas constant from session file and compute Cp
145  m_session->LoadParameter ("GasConstant", gasConstant, 287.058);
146  m_Cp = m_gamma / (m_gamma - 1.0) * gasConstant;
147 
148  // Get pInf parameter from session file.
149  m_session->LoadParameter("pInf", m_pInf, 101325);
150 
151  // Get rhoInf parameter from session file.
152  m_session->LoadParameter("rhoInf", m_rhoInf, 1.225);
153 
154  // Get uInf parameter from session file.
155  m_session->LoadParameter("uInf", velInf, 0.1);
156 
157  m_UInf = velInf*velInf;
158 
159  // Get vInf parameter from session file.
160  if (m_spacedim == 2 || m_spacedim == 3)
161  {
162  m_session->LoadParameter("vInf", velInf, 0.0);
163  m_UInf += velInf*velInf;
164  }
165 
166  // Get wInf parameter from session file.
167  if (m_spacedim == 3)
168  {
169  m_session->LoadParameter("wInf", velInf, 0.0);
170  m_UInf += velInf*velInf;
171  }
172  m_UInf = sqrt(m_UInf);
173 
174  // Viscosity
175  m_session->LoadSolverInfo("ViscosityType", m_ViscosityType, "Constant");
176  m_session->LoadParameter ("mu", m_mu, 1.78e-05);
177 
178  // Thermal conductivity or Prandtl
179  if( m_session->DefinesParameter("thermalConductivity"))
180  {
181  ASSERTL0( !m_session->DefinesParameter("Pr"),
182  "Cannot define both Pr and thermalConductivity.");
183 
184  m_session->LoadParameter ("thermalConductivity",
187  }
188  else
189  {
190  m_session->LoadParameter ("Pr",
191  m_Prandtl, 0.72);
193  }
194 
195  // Steady state tolerance
196  m_session->LoadParameter("SteadyStateTol", m_steadyStateTol, 0.0);
197 
198  // Shock capture
199  m_session->LoadSolverInfo("ShockCaptureType",
200  m_shockCaptureType, "Off");
201  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
int m_spacedim
Spatial dimension (>= expansion dim).
double NekDouble
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
void Nektar::CompressibleFlowSystem::InitializeSteadyState ( )
protected

Definition at line 734 of file CompressibleFlowSystem.cpp.

References Nektar::SolverUtils::EquationSystem::m_comm, m_errFile, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_session, m_un, and Vmath::Vcopy().

Referenced by Nektar::IsentropicVortex::v_SetInitialConditions(), and v_SetInitialConditions().

735  {
736  if (m_session->DefinesParameter("SteadyStateTol"))
737  {
738  const int nPoints = m_fields[0]->GetTotPoints();
739  m_un = Array<OneD, Array<OneD, NekDouble> > (
740  m_fields.num_elements());
741 
742  for (int i = 0; i < m_fields.num_elements(); ++i)
743  {
744  m_un[i] = Array<OneD, NekDouble>(nPoints);
745  Vmath::Vcopy(nPoints, m_fields[i]->GetPhys(), 1, m_un[i], 1);
746  }
747 
748  if (m_comm->GetRank() == 0)
749  {
750  std::string fName = m_session->GetSessionName() +
751  std::string(".res");
752  m_errFile.open(fName.c_str());
753  m_errFile << "# "
754  << setw(15) << left << "Time"
755  << setw(22) << left << "rho";
756 
757  std::string velFields[3] = {"u", "v", "w"};
758 
759  for (int i = 0; i < m_fields.num_elements()-2; ++i)
760  {
761  m_errFile << setw(22) << "rho"+velFields[i];
762  }
763 
764  m_errFile << setw(22) << left << "E" << endl;
765  }
766  }
767  }
Array< OneD, Array< OneD, NekDouble > > m_un
LibUtilities::CommSharedPtr m_comm
Communicator.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
void Nektar::CompressibleFlowSystem::SetBoundaryConditions ( Array< OneD, Array< OneD, NekDouble > > &  physarray,
NekDouble  time 
)
protected

Definition at line 374 of file CompressibleFlowSystem.cpp.

References Nektar::SolverUtils::EquationSystem::GetTraceTotPoints(), Nektar::iterator, m_bndConds, and Nektar::SolverUtils::EquationSystem::m_fields.

Referenced by DoOdeProjection().

377  {
378  int nTracePts = GetTraceTotPoints();
379  int nvariables = physarray.num_elements();
380 
381  Array<OneD, Array<OneD, NekDouble> > Fwd(nvariables);
382  for (int i = 0; i < nvariables; ++i)
383  {
384  Fwd[i] = Array<OneD, NekDouble>(nTracePts);
385  m_fields[i]->ExtractTracePhys(physarray[i], Fwd[i]);
386  }
387 
388  if (m_bndConds.size())
389  {
390  // Loop over user-defined boundary conditions
392  for (x = m_bndConds.begin(); x != m_bndConds.end(); ++x)
393  {
394  (*x)->Apply(Fwd, physarray, time);
395  }
396  }
397  }
SOLVER_UTILS_EXPORT int GetTraceTotPoints()
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
std::vector< CFSBndCondSharedPtr > m_bndConds
virtual void Nektar::CompressibleFlowSystem::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 
)
inlineprotectedvirtual

Reimplemented in Nektar::NavierStokesCFE.

Definition at line 180 of file CompressibleFlowSystem.h.

Referenced by DoDiffusion().

185  {
186  // Do nothing by default
187  }
void Nektar::CompressibleFlowSystem::v_ExtraFldOutput ( std::vector< Array< OneD, NekDouble > > &  fieldcoeffs,
std::vector< std::string > &  variables 
)
protectedvirtual

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 915 of file CompressibleFlowSystem.cpp.

References m_artificialDiffusion, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_session, and m_varConv.

918  {
919  bool extraFields;
920  m_session->MatchSolverInfo("OutputExtraFields","True",
921  extraFields, true);
922  if (extraFields)
923  {
924  const int nPhys = m_fields[0]->GetNpoints();
925  const int nCoeffs = m_fields[0]->GetNcoeffs();
926  Array<OneD, Array<OneD, NekDouble> > tmp(m_fields.num_elements());
927 
928  for (int i = 0; i < m_fields.num_elements(); ++i)
929  {
930  tmp[i] = m_fields[i]->GetPhys();
931  }
932 
933  Array<OneD, NekDouble> pressure(nPhys), soundspeed(nPhys), mach(nPhys);
934  Array<OneD, NekDouble> sensor(nPhys), SensorKappa(nPhys);
935 
936  m_varConv->GetPressure (tmp, pressure);
937  m_varConv->GetSoundSpeed(tmp, pressure, soundspeed);
938  m_varConv->GetMach (tmp, soundspeed, mach);
939  m_varConv->GetSensor (m_fields[0], tmp, sensor, SensorKappa);
940 
941  Array<OneD, NekDouble> pFwd(nCoeffs), sFwd(nCoeffs), mFwd(nCoeffs);
942  Array<OneD, NekDouble> sensFwd(nCoeffs);
943 
944  m_fields[0]->FwdTrans_IterPerExp(pressure, pFwd);
945  m_fields[0]->FwdTrans_IterPerExp(soundspeed, sFwd);
946  m_fields[0]->FwdTrans_IterPerExp(mach, mFwd);
947  m_fields[0]->FwdTrans_IterPerExp(sensor, sensFwd);
948 
949  variables.push_back ("p");
950  variables.push_back ("a");
951  variables.push_back ("Mach");
952  variables.push_back ("Sensor");
953  fieldcoeffs.push_back(pFwd);
954  fieldcoeffs.push_back(sFwd);
955  fieldcoeffs.push_back(mFwd);
956  fieldcoeffs.push_back(sensFwd);
957 
959  {
960  // reuse pressure
961  m_artificialDiffusion->GetArtificialViscosity(tmp, pressure);
962  m_fields[0]->FwdTrans_IterPerExp(pressure, pFwd);
963 
964  variables.push_back ("ArtificialVisc");
965  fieldcoeffs.push_back(pFwd);
966  }
967  }
968  }
VariableConverterSharedPtr m_varConv
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
ArtificialDiffusionSharedPtr m_artificialDiffusion
NekDouble Nektar::CompressibleFlowSystem::v_GetTimeStep ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray)
protectedvirtual

Calculate the maximum timestep subject to CFL restrictions.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 650 of file CompressibleFlowSystem.cpp.

References Nektar::SolverUtils::EquationSystem::GetNumExpModesPerExp(), GetStdVelocity(), Nektar::SolverUtils::UnsteadySystem::m_cflSafetyFactor, Nektar::SolverUtils::EquationSystem::m_comm, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::UnsteadySystem::MaxTimeStepEstimator(), Nektar::LibUtilities::ReduceMin, and Vmath::Vmin().

652  {
653  int n;
654  int nElements = m_fields[0]->GetExpSize();
655  const Array<OneD, int> ExpOrder = GetNumExpModesPerExp();
656 
657  Array<OneD, NekDouble> tstep (nElements, 0.0);
658  Array<OneD, NekDouble> stdVelocity(nElements);
659 
660  // Get standard velocity to compute the time-step limit
661  GetStdVelocity(inarray, stdVelocity);
662 
663  // Factors to compute the time-step limit
664  NekDouble minLength = 0.0;
666  NekDouble cLambda = 0.2; // Spencer book-317
667 
668  // Loop over elements to compute the time-step limit for each element
669  for(n = 0; n < nElements; ++n)
670  {
671  int npoints = m_fields[0]->GetExp(n)->GetTotPoints();
672  Array<OneD, NekDouble> one2D(npoints, 1.0);
673  NekDouble Area = m_fields[0]->GetExp(n)->Integral(one2D);
674 
675  minLength = sqrt(Area);
676  if (m_fields[0]->GetExp(n)->as<LocalRegions::TriExp>())
677  {
678  minLength *= 2.0;
679  }
680 
681  tstep[n] = m_cflSafetyFactor * alpha * minLength
682  / (stdVelocity[n] * cLambda
683  * (ExpOrder[n] - 1) * (ExpOrder[n] - 1));
684  }
685 
686  // Get the minimum time-step limit and return the time-step
687  NekDouble TimeStep = Vmath::Vmin(nElements, tstep, 1);
688  m_comm->AllReduce(TimeStep, LibUtilities::ReduceMin);
689  return TimeStep;
690  }
void GetStdVelocity(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &stdV)
Compute the advection velocity in the standard space for each element of the expansion.
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:871
SOLVER_UTILS_EXPORT const Array< OneD, int > GetNumExpModesPerExp()
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator()
Get the maximum timestep estimator for cfl control.
LibUtilities::CommSharedPtr m_comm
Communicator.
double NekDouble
NekDouble m_cflSafetyFactor
CFL safety factor (comprise between 0 to 1).
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
void Nektar::CompressibleFlowSystem::v_InitObject ( )
protectedvirtual

Initialization object for CompressibleFlowSystem class.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Reimplemented in Nektar::NavierStokesCFE.

Definition at line 54 of file CompressibleFlowSystem.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, >::CreateInstance(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineOdeRhs(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineProjection(), DoOdeProjection(), DoOdeRhs(), Nektar::GetArtificialDiffusionFactory(), Nektar::GetCFSBndCondFactory(), InitAdvection(), InitialiseParameters(), Nektar::SolverUtils::Forcing::Load(), m_artificialDiffusion, m_bndConds, Nektar::SolverUtils::UnsteadySystem::m_explicitAdvection, Nektar::SolverUtils::EquationSystem::m_fields, m_forcing, Nektar::SolverUtils::UnsteadySystem::m_homoInitialFwd, Nektar::SolverUtils::UnsteadySystem::m_ode, Nektar::SolverUtils::EquationSystem::m_session, m_shockCaptureType, Nektar::SolverUtils::EquationSystem::m_spacedim, Nektar::SolverUtils::EquationSystem::m_traceNormals, m_varConv, m_vecLocs, and Nektar::SolverUtils::UnsteadySystem::v_InitObject().

Referenced by Nektar::NavierStokesCFE::v_InitObject().

55  {
57 
60 
61  ASSERTL0(m_session->DefinesSolverInfo("UPWINDTYPE"),
62  "No UPWINDTYPE defined in session.");
63 
64  // Do not forwards transform initial condition
65  m_homoInitialFwd = false;
66 
67  // Set up locations of velocity vector.
68  m_vecLocs = Array<OneD, Array<OneD, NekDouble> >(1);
69  m_vecLocs[0] = Array<OneD, NekDouble>(m_spacedim);
70  for (int i = 0; i < m_spacedim; ++i)
71  {
72  m_vecLocs[0][i] = 1 + i;
73  }
74 
75  // Loading parameters from session file
77 
78  // Setting up advection and diffusion operators
79  InitAdvection();
80 
81  // Create artificial diffusion
82  if (m_shockCaptureType != "Off")
83  {
86  m_session,
87  m_fields,
88  m_spacedim);
89  }
90 
91  // Forcing terms for the sponge region
93  m_fields.num_elements());
94 
95  // User-defined boundary conditions
96  int cnt = 0;
97  for (int n = 0; n < m_fields[0]->GetBndConditions().num_elements(); ++n)
98  {
99  std::string type =
100  m_fields[0]->GetBndConditions()[n]->GetUserDefined();
101  if(!type.empty())
102  {
103  m_bndConds.push_back(GetCFSBndCondFactory().CreateInstance(
104  type,
105  m_session,
106  m_fields,
108  m_spacedim,
109  n,
110  cnt));
111  }
112  cnt += m_fields[0]->GetBndCondExpansions()[n]->GetExpSize();
113  }
114 
116  {
119  }
120  else
121  {
122  ASSERTL0(false, "Implicit CFS not set up.");
123  }
124  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
bool m_homoInitialFwd
Flag to determine if simulation should start in homogeneous forward transformed state.
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
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
VariableConverterSharedPtr m_varConv
static SOLVER_UTILS_EXPORT std::vector< ForcingSharedPtr > Load(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields=0)
Definition: Forcing.cpp:86
Array< OneD, Array< OneD, NekDouble > > m_vecLocs
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 discontinu...
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
CFSBndCondFactory & GetCFSBndCondFactory()
Declaration of the boundary condition factory singleton.
Definition: CFSBndCond.cpp:42
bool m_explicitAdvection
Indicates if explicit or implicit treatment of advection is used.
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
int m_spacedim
Spatial dimension (>= expansion dim).
void InitialiseParameters()
Load CFS parameters from the session file.
ArtificialDiffusionFactory & GetArtificialDiffusionFactory()
Declaration of the artificial diffusion factory singleton.
virtual SOLVER_UTILS_EXPORT void v_InitObject()
Init object for UnsteadySystem class.
void InitAdvection()
Create advection and diffusion objects for CFS.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
ArtificialDiffusionSharedPtr m_artificialDiffusion
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
std::vector< CFSBndCondSharedPtr > m_bndConds
bool Nektar::CompressibleFlowSystem::v_PostIntegrate ( int  step)
protectedvirtual

Perform post-integration checks, presently just to check steady state behaviour.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 557 of file CompressibleFlowSystem.cpp.

References CalcSteadyState(), Nektar::SolverUtils::EquationSystem::m_comm, Nektar::SolverUtils::UnsteadySystem::m_infosteps, and m_steadyStateTol.

558  {
559  if (m_steadyStateTol > 0.0)
560  {
561  bool doOutput = step % m_infosteps == 0;
562  if (CalcSteadyState(doOutput))
563  {
564  if (m_comm->GetRank() == 0)
565  {
566  cout << "Reached Steady State to tolerance "
567  << m_steadyStateTol << endl;
568  }
569  return true;
570  }
571  }
572  return false;
573  }
LibUtilities::CommSharedPtr m_comm
Communicator.
bool CalcSteadyState(bool output)
Calculate whether the system has reached a steady state by observing residuals to a user-defined tole...
int m_infosteps
Number of time steps between outputting status information.
void Nektar::CompressibleFlowSystem::v_SetInitialConditions ( NekDouble  initialtime = 0.0,
bool  dumpInitialConditions = true,
const int  domain = 0 
)
protectedvirtual

Set up logic for residual calculation.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Reimplemented in Nektar::IsentropicVortex.

Definition at line 695 of file CompressibleFlowSystem.cpp.

References Nektar::SolverUtils::EquationSystem::Checkpoint_Output(), Vmath::FillWhiteNoise(), InitializeSteadyState(), Nektar::SolverUtils::EquationSystem::m_comm, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::v_SetInitialConditions(), and Vmath::Vadd().

699  {
700  EquationSystem::v_SetInitialConditions(initialtime, false);
701 
702  // insert white noise in initial condition
703  NekDouble Noise;
704  int phystot = m_fields[0]->GetTotPoints();
705  Array<OneD, NekDouble> noise(phystot);
706 
707  m_session->LoadParameter("Noise", Noise,0.0);
708  int m_nConvectiveFields = m_fields.num_elements();
709 
710  if (Noise > 0.0)
711  {
712  int seed = - m_comm->GetRank()*m_nConvectiveFields;
713  for (int i = 0; i < m_nConvectiveFields; i++)
714  {
715  Vmath::FillWhiteNoise(phystot, Noise, noise, 1,
716  seed);
717  --seed;
718  Vmath::Vadd(phystot, m_fields[i]->GetPhys(), 1,
719  noise, 1, m_fields[i]->UpdatePhys(), 1);
720  m_fields[i]->FwdTrans_IterPerExp(m_fields[i]->GetPhys(),
721  m_fields[i]->UpdateCoeffs());
722  }
723  }
724 
726 
727  if (dumpInitialConditions)
728  {
729  // Dump initial conditions to file
731  }
732  }
SOLVER_UTILS_EXPORT void Checkpoint_Output(const int n)
Write checkpoint file of m_fields.
LibUtilities::CommSharedPtr m_comm
Communicator.
double NekDouble
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
void FillWhiteNoise(int n, const T eps, T *x, const int incx, int outseed)
Fills a vector with white noise.
Definition: Vmath.cpp:138
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:299

Friends And Related Function Documentation

friend class MemoryManager< CompressibleFlowSystem >
friend

Definition at line 57 of file CompressibleFlowSystem.h.

Member Data Documentation

SolverUtils::AdvectionSharedPtr Nektar::CompressibleFlowSystem::m_advection
protected

Definition at line 70 of file CompressibleFlowSystem.h.

Referenced by DoAdvection(), and InitAdvection().

ArtificialDiffusionSharedPtr Nektar::CompressibleFlowSystem::m_artificialDiffusion
protected

Definition at line 72 of file CompressibleFlowSystem.h.

Referenced by DoDiffusion(), v_ExtraFldOutput(), and v_InitObject().

std::vector<CFSBndCondSharedPtr> Nektar::CompressibleFlowSystem::m_bndConds
protected

Definition at line 89 of file CompressibleFlowSystem.h.

Referenced by SetBoundaryConditions(), and v_InitObject().

NekDouble Nektar::CompressibleFlowSystem::m_Cp
protected
SolverUtils::DiffusionSharedPtr Nektar::CompressibleFlowSystem::m_diffusion
protected
std::ofstream Nektar::CompressibleFlowSystem::m_errFile
protected

Definition at line 92 of file CompressibleFlowSystem.h.

Referenced by CalcSteadyState(), and InitializeSteadyState().

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

Definition at line 98 of file CompressibleFlowSystem.h.

Referenced by DoOdeRhs(), and v_InitObject().

NekDouble Nektar::CompressibleFlowSystem::m_gamma
protected
NekDouble Nektar::CompressibleFlowSystem::m_mu
protected
NekDouble Nektar::CompressibleFlowSystem::m_pInf
protected

Definition at line 75 of file CompressibleFlowSystem.h.

Referenced by CalcSteadyState(), and InitialiseParameters().

NekDouble Nektar::CompressibleFlowSystem::m_Prandtl
protected
NekDouble Nektar::CompressibleFlowSystem::m_rhoInf
protected

Definition at line 76 of file CompressibleFlowSystem.h.

Referenced by CalcSteadyState(), and InitialiseParameters().

std::string Nektar::CompressibleFlowSystem::m_shockCaptureType
protected

Definition at line 79 of file CompressibleFlowSystem.h.

Referenced by DoDiffusion(), InitialiseParameters(), and v_InitObject().

NekDouble Nektar::CompressibleFlowSystem::m_steadyStateTol
protected

Definition at line 95 of file CompressibleFlowSystem.h.

Referenced by CalcSteadyState(), InitialiseParameters(), and v_PostIntegrate().

NekDouble Nektar::CompressibleFlowSystem::m_thermalConductivity
protected
NekDouble Nektar::CompressibleFlowSystem::m_UInf
protected

Definition at line 77 of file CompressibleFlowSystem.h.

Referenced by CalcSteadyState(), and InitialiseParameters().

Array<OneD, Array<OneD, NekDouble> > Nektar::CompressibleFlowSystem::m_un
protected

Definition at line 101 of file CompressibleFlowSystem.h.

Referenced by CalcSteadyState(), and InitializeSteadyState().

VariableConverterSharedPtr Nektar::CompressibleFlowSystem::m_varConv
protected
Array<OneD, Array<OneD, NekDouble> > Nektar::CompressibleFlowSystem::m_vecLocs
protected

Definition at line 73 of file CompressibleFlowSystem.h.

Referenced by GetVecLocs(), and v_InitObject().

std::string Nektar::CompressibleFlowSystem::m_ViscosityType
protected