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

#include <NavierStokesCFEAxisym.h>

Inheritance diagram for Nektar::NavierStokesCFEAxisym:
[legend]

Public Member Functions

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

Static Public Member Functions

static SolverUtils::EquationSystemSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
- Static Public Member Functions inherited from Nektar::NavierStokesCFE
static SolverUtils::EquationSystemSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 

Static Public Attributes

static std::string className
 
- Static Public Attributes inherited from Nektar::NavierStokesCFE
static std::string className
 

Protected Member Functions

 NavierStokesCFEAxisym (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
virtual void v_InitObject ()
 Initialization object for CompressibleFlowSystem class. More...
 
virtual void v_DoDiffusion (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
 
virtual void v_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. More...
 
virtual void v_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. More...
 
- Protected Member Functions inherited from Nektar::NavierStokesCFE
 NavierStokesCFE (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
virtual void v_GetFluxPenalty (const Array< OneD, Array< OneD, NekDouble > > &uFwd, const Array< OneD, Array< OneD, NekDouble > > &uBwd, Array< OneD, Array< OneD, NekDouble > > &penaltyCoeff)
 Return the penalty vector for the LDGNS diffusion problem. More...
 
void GetViscosityAndThermalCondFromTemp (const Array< OneD, NekDouble > &temperature, Array< OneD, NekDouble > &mu, Array< OneD, NekDouble > &thermalCond)
 Update viscosity todo: add artificial viscosity here. More...
 
- Protected Member Functions inherited from Nektar::CompressibleFlowSystem
 CompressibleFlowSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
void InitialiseParameters ()
 Load CFS parameters from the session file. More...
 
void InitAdvection ()
 Create advection and diffusion objects for CFS. More...
 
void DoOdeRhs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 Compute the right-hand side. More...
 
void DoOdeProjection (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 Compute the projection and call the method for imposing the boundary conditions in case of discontinuous projection. More...
 
void DoAdvection (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
 Compute the advection terms for the right-hand side. More...
 
void DoDiffusion (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
 Add the diffusions terms to the right-hand side. More...
 
void GetFluxVector (const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
 Return the flux vector for the compressible Euler equations. More...
 
void GetFluxVectorDeAlias (const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
 Return the flux vector for the compressible Euler equations by using the de-aliasing technique. More...
 
void SetBoundaryConditions (Array< OneD, Array< OneD, NekDouble > > &physarray, NekDouble time)
 
void GetElmtTimeStep (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &tstep)
 Calculate the maximum timestep on each element subject to CFL restrictions. More...
 
virtual NekDouble v_GetTimeStep (const Array< OneD, const Array< OneD, NekDouble > > &inarray)
 Calculate the maximum timestep subject to CFL restrictions. More...
 
virtual void v_SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 Set up logic for residual calculation. More...
 
NekDouble GetGamma ()
 
const Array< OneD, const Array< OneD, NekDouble > > & GetVecLocs ()
 
const Array< OneD, const Array< OneD, NekDouble > > & GetNormals ()
 
virtual void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 
virtual Array< OneD, NekDoublev_GetMaxStdVelocity ()
 Compute the advection velocity in the standard space for each element of the expansion. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::AdvectionSystem
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step)
 
- Protected Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
SOLVER_UTILS_EXPORT UnsteadySystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises UnsteadySystem class members. More...
 
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve ()
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise ()
 Sets up initial conditions. More...
 
virtual SOLVER_UTILS_EXPORT void v_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 NekDouble v_GetTimeStep (const Array< OneD, const Array< OneD, NekDouble >> &inarray)
 Return the timestep to be used for the next step in the time-marching loop. More...
 
virtual SOLVER_UTILS_EXPORT bool v_PreIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT bool v_RequireFwdTrans ()
 
SOLVER_UTILS_EXPORT void CheckForRestartTime (NekDouble &time, int &nchk)
 
SOLVER_UTILS_EXPORT void SVVVarDiffCoeff (const Array< OneD, Array< OneD, NekDouble >> vel, StdRegions::VarCoeffMap &varCoeffMap)
 Evaluate the SVV diffusion coefficient according to Moura's paper where it should proportional to h time velocity. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises EquationSystem class members. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Virtual function for the L_inf error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false)
 Virtual function for the L_2 error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransCoeffToPhys ()
 Virtual function for transformation to physical space. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransPhysToCoeff ()
 Virtual function for transformation to coefficient space. More...
 
virtual SOLVER_UTILS_EXPORT void v_EvaluateExactSolution (unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 
virtual SOLVER_UTILS_EXPORT void v_Output (void)
 
virtual SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr v_GetPressure (void)
 

Protected Attributes

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

Friends

class MemoryManager< NavierStokesCFEAxisym >
 

Additional Inherited Members

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

Detailed Description

Definition at line 48 of file NavierStokesCFEAxisym.h.

Constructor & Destructor Documentation

◆ ~NavierStokesCFEAxisym()

Nektar::NavierStokesCFEAxisym::~NavierStokesCFEAxisym ( )
virtual

Definition at line 54 of file NavierStokesCFEAxisym.cpp.

55  {
56 
57  }

◆ NavierStokesCFEAxisym()

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

Definition at line 46 of file NavierStokesCFEAxisym.cpp.

49  : UnsteadySystem(pSession, pGraph),
50  NavierStokesCFE(pSession, pGraph)
51  {
52  }
NavierStokesCFE(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
SOLVER_UTILS_EXPORT UnsteadySystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises UnsteadySystem class members.

Member Function Documentation

◆ create()

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

Definition at line 54 of file NavierStokesCFEAxisym.h.

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

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

◆ v_DoDiffusion()

void Nektar::NavierStokesCFEAxisym::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 
)
protectedvirtual

Reimplemented from Nektar::NavierStokesCFE.

Definition at line 72 of file NavierStokesCFEAxisym.cpp.

References Nektar::SolverUtils::EquationSystem::GetNpoints(), m_viscousForcing, Nektar::NavierStokesCFE::v_DoDiffusion(), and Vmath::Vadd().

77  {
78  int npoints = GetNpoints();
79  int nvariables = inarray.num_elements();
80 
81  NavierStokesCFE::v_DoDiffusion(inarray, outarray, pFwd, pBwd);
82 
83  for (int i = 0; i < nvariables; ++i)
84  {
85  Vmath::Vadd(npoints,
86  m_viscousForcing[i], 1,
87  outarray[i], 1,
88  outarray[i], 1);
89  }
90  }
Array< OneD, Array< OneD, NekDouble > > m_viscousForcing
SOLVER_UTILS_EXPORT int GetNpoints()
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)
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:302

◆ v_GetViscousFluxVector()

void Nektar::NavierStokesCFEAxisym::v_GetViscousFluxVector ( const Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &  derivativesO1,
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &  viscousTensor 
)
protectedvirtual

Return the flux vector for the LDG diffusion problem.

Todo:
Complete the viscous flux vector

Reimplemented from Nektar::NavierStokesCFE.

Definition at line 96 of file NavierStokesCFEAxisym.cpp.

References Nektar::NavierStokesCFE::GetViscosityAndThermalCondFromTemp(), Nektar::NekConstants::kNekZeroTol, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::NavierStokesCFE::m_mu, Nektar::SolverUtils::EquationSystem::m_spacedim, Nektar::NavierStokesCFE::m_thermalConductivity, m_viscousForcing, Vmath::Smul(), Vmath::Vadd(), Vmath::Vcopy(), Vmath::Vmul(), Vmath::Vsub(), Vmath::Vvtvp(), and Vmath::Zero().

100  {
101  int i, j;
102  int nVariables = m_fields.num_elements();
103  int nPts = physfield[0].num_elements();
104 
105  // 1/r
106  Array<OneD, Array<OneD, NekDouble> > coords(3);
107  Array<OneD, NekDouble> invR (nPts,0.0);
108  for (int i = 0; i < 3; i++)
109  {
110  coords[i] = Array<OneD, NekDouble> (nPts);
111  }
112  m_fields[0]->GetCoords(coords[0], coords[1], coords[2]);
113  for (int i = 0; i < nPts; ++i)
114  {
115  if (coords[0][i] < NekConstants::kNekZeroTol)
116  {
117  invR[i] = 0;
118  }
119  else
120  {
121  invR[i] = 1.0/coords[0][i];
122  }
123  }
124 
125  // Stokes hypothesis
126  const NekDouble lambda = -2.0/3.0;
127 
128  // Auxiliary variables
129  Array<OneD, NekDouble > divVel (nPts, 0.0);
130  Array<OneD, NekDouble > tmp (nPts, 0.0);
131 
132  // Update viscosity and thermal conductivity
133  GetViscosityAndThermalCondFromTemp(physfield[nVariables-2], m_mu,
135 
136  // Velocity divergence = d(u_r)/dr + d(u_z)/dz + u_r/r
137  Vmath::Vadd(nPts, derivativesO1[0][0], 1, derivativesO1[1][1], 1,
138  divVel, 1);
139  Vmath::Vvtvp(nPts, physfield[0], 1 , invR, 1, divVel, 1, divVel, 1);
140 
141  // Velocity divergence scaled by lambda * mu
142  Vmath::Smul(nPts, lambda, divVel, 1, divVel, 1);
143  Vmath::Vmul(nPts, m_mu, 1, divVel, 1, divVel, 1);
144 
145  // Viscous flux vector for the rho equation = 0
146  for (i = 0; i < m_spacedim; ++i)
147  {
148  Vmath::Zero(nPts, viscousTensor[i][0], 1);
149  }
150 
151  // Viscous stress tensor (for the momentum equations)
152 
153  for (i = 0; i < 2; ++i)
154  {
155  for (j = i; j < 2; ++j)
156  {
157  Vmath::Vadd(nPts, derivativesO1[i][j], 1,
158  derivativesO1[j][i], 1,
159  viscousTensor[i][j+1], 1);
160 
161  Vmath::Vmul(nPts, m_mu, 1,
162  viscousTensor[i][j+1], 1,
163  viscousTensor[i][j+1], 1);
164 
165  if (i == j)
166  {
167  // Add divergence term to diagonal
168  Vmath::Vadd(nPts, viscousTensor[i][j+1], 1,
169  divVel, 1,
170  viscousTensor[i][j+1], 1);
171  }
172  else
173  {
174  // Copy to make symmetric
175  Vmath::Vcopy(nPts, viscousTensor[i][j+1], 1,
176  viscousTensor[j][i+1], 1);
177  }
178  }
179  }
180  // Swirl case
181  if(m_spacedim == 3)
182  {
183  // Tau_theta_theta = mu ( 2*u_r/r - 2/3*div(u))
184  Vmath::Vmul(nPts, physfield[0], 1 , invR, 1,
185  viscousTensor[2][3], 1);
186  Vmath::Smul(nPts, 2.0, viscousTensor[2][3], 1,
187  viscousTensor[2][3], 1);
188  Vmath::Vmul(nPts, m_mu, 1, viscousTensor[2][3], 1,
189  viscousTensor[2][3], 1);
190  Vmath::Vadd(nPts, viscousTensor[2][3], 1,
191  divVel, 1,
192  viscousTensor[2][3], 1);
193 
194  // Tau_r_theta = mu (-u_theta/r + d(u_theta)/dr )
195  Vmath::Vmul(nPts, physfield[2], 1 , invR, 1,
196  viscousTensor[2][1], 1);
197  Vmath::Smul(nPts, -1.0, viscousTensor[2][1], 1,
198  viscousTensor[2][1], 1);
199  Vmath::Vadd(nPts, derivativesO1[0][2], 1 , viscousTensor[2][1], 1,
200  viscousTensor[2][1], 1);
201  Vmath::Vmul(nPts, m_mu, 1, viscousTensor[2][1], 1,
202  viscousTensor[2][1], 1);
203  Vmath::Vcopy(nPts, viscousTensor[2][1], 1,
204  viscousTensor[0][3], 1);
205 
206  // Tau_z_theta = mu (d(u_theta)/dz )
207  Vmath::Vmul(nPts, m_mu, 1, derivativesO1[1][2], 1,
208  viscousTensor[2][2], 1);
209  Vmath::Vcopy(nPts, viscousTensor[2][2], 1,
210  viscousTensor[1][3], 1);
211  }
212 
213  // Terms for the energy equation
214  for (i = 0; i < m_spacedim; ++i)
215  {
216  Vmath::Zero(nPts, viscousTensor[i][m_spacedim+1], 1);
217  // u_j * tau_ij
218  for (j = 0; j < m_spacedim; ++j)
219  {
220  Vmath::Vvtvp(nPts, physfield[j], 1,
221  viscousTensor[i][j+1], 1,
222  viscousTensor[i][m_spacedim+1], 1,
223  viscousTensor[i][m_spacedim+1], 1);
224  }
225  // Add k*T_i
226  if (i != 2)
227  {
229  derivativesO1[i][m_spacedim], 1,
230  viscousTensor[i][m_spacedim+1], 1,
231  viscousTensor[i][m_spacedim+1], 1);
232  }
233  else
234  {
235  Vmath::Vmul(nPts, derivativesO1[i][m_spacedim], 1 ,
236  invR, 1, tmp, 1);
238  tmp, 1,
239  viscousTensor[i][m_spacedim+1], 1,
240  viscousTensor[i][m_spacedim+1], 1);
241  }
242  }
243 
244  // Update viscous forcing
245  // r-momentum: F = 1/r * (tau_rr - tau_theta_theta)
246  if(m_spacedim == 3)
247  {
248  Vmath::Vsub(nPts, viscousTensor[0][1], 1, viscousTensor[2][3], 1,
249  m_viscousForcing[1], 1);
250  Vmath::Vmul(nPts, m_viscousForcing[1], 1 ,
251  invR, 1, m_viscousForcing[1], 1);
252  }
253  else
254  {
255  Vmath::Vmul(nPts, viscousTensor[0][1], 1 ,
256  invR, 1, m_viscousForcing[1], 1);
257  }
258 
259  // z-momentum: F = 1/r * tau_r_z
260  Vmath::Vmul(nPts, viscousTensor[0][2], 1 ,
261  invR, 1, m_viscousForcing[2], 1);
262 
263  // Theta_momentum: F = 2* tau_r_theta
264  if(m_spacedim == 3)
265  {
266  Vmath::Vmul(nPts, viscousTensor[0][3], 1 ,
267  invR, 1, m_viscousForcing[3], 1);
268  Vmath::Smul(nPts, 2.0, m_viscousForcing[3], 1,
269  m_viscousForcing[3], 1);
270  }
271 
272  // Energy: F = 1/r* viscousTensor_T_r
273  Vmath::Vmul(nPts, viscousTensor[0][m_spacedim+1], 1 ,
274  invR, 1, m_viscousForcing[m_spacedim+1], 1);
275  }
Array< OneD, NekDouble > m_thermalConductivity
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:445
static const NekDouble kNekZeroTol
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:216
int m_spacedim
Spatial dimension (>= expansion dim).
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:346
Array< OneD, Array< OneD, NekDouble > > m_viscousForcing
void GetViscosityAndThermalCondFromTemp(const Array< OneD, NekDouble > &temperature, Array< OneD, NekDouble > &mu, Array< OneD, NekDouble > &thermalCond)
Update viscosity todo: add artificial viscosity here.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
Array< OneD, NekDouble > m_mu
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:376
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:302
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:186

◆ v_GetViscousFluxVectorDeAlias()

virtual void Nektar::NavierStokesCFEAxisym::v_GetViscousFluxVectorDeAlias ( const Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &  derivativesO1,
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &  viscousTensor 
)
inlineprotectedvirtual

Return the flux vector for the LDG diffusion problem.

Todo:
Complete the viscous flux vector

Reimplemented from Nektar::NavierStokesCFE.

Definition at line 88 of file NavierStokesCFEAxisym.h.

References Nektar::ErrorUtil::efatal, and NEKERROR.

92  {
93  boost::ignore_unused(physfield, derivatives, viscousTensor);
95  "Dealiased flux not implemented for axisymmetric case");
96  }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:209

◆ v_InitObject()

void Nektar::NavierStokesCFEAxisym::v_InitObject ( )
protectedvirtual

Initialization object for CompressibleFlowSystem class.

Reimplemented from Nektar::NavierStokesCFE.

Definition at line 59 of file NavierStokesCFEAxisym.cpp.

References Nektar::SolverUtils::EquationSystem::GetNpoints(), Nektar::SolverUtils::EquationSystem::m_fields, m_viscousForcing, and Nektar::NavierStokesCFE::v_InitObject().

60  {
62 
63  int nVariables = m_fields.num_elements();
64  int npoints = GetNpoints();
65  m_viscousForcing = Array<OneD, Array<OneD, NekDouble>> (nVariables);
66  for (int i = 0; i < nVariables; ++i)
67  {
68  m_viscousForcing[i] = Array<OneD, NekDouble>(npoints, 0.0);
69  }
70  }
virtual void v_InitObject()
Initialization object for CompressibleFlowSystem class.
Array< OneD, Array< OneD, NekDouble > > m_viscousForcing
SOLVER_UTILS_EXPORT int GetNpoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.

Friends And Related Function Documentation

◆ MemoryManager< NavierStokesCFEAxisym >

friend class MemoryManager< NavierStokesCFEAxisym >
friend

Definition at line 51 of file NavierStokesCFEAxisym.h.

Member Data Documentation

◆ className

string Nektar::NavierStokesCFEAxisym::className
static
Initial value:
=
"NavierStokesCFEAxisym", NavierStokesCFEAxisym::create,
"Axisymmetric NavierStokes equations in conservative variables.")

Definition at line 65 of file NavierStokesCFEAxisym.h.

◆ m_viscousForcing

Array<OneD, Array<OneD, NekDouble> > Nektar::NavierStokesCFEAxisym::m_viscousForcing
protected

Definition at line 70 of file NavierStokesCFEAxisym.h.

Referenced by v_DoDiffusion(), v_GetViscousFluxVector(), and v_InitObject().