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

#include <NavierStokesCFE.h>

Inheritance diagram for Nektar::NavierStokesCFE:
[legend]

Public Member Functions

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...
 
virtual void v_GetPressure (const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &pressure) override
 
virtual void v_GetDensity (const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &density) override
 
virtual bool v_HasConstantDensity () override
 
virtual void v_GetVelocity (const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, Array< OneD, NekDouble >> &velocity) override
 
- 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 (const bool FlagAcousticCFL=true)
 
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...
 
SOLVER_UTILS_EXPORT void SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
- 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 (bool DeclareField=true)
 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, const Array< OneD, const NekDouble > &input)
 
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 int GetInfoSteps ()
 
SOLVER_UTILS_EXPORT void SetInfoSteps (int num)
 
SOLVER_UTILS_EXPORT int GetPararealIterationNumber ()
 
SOLVER_UTILS_EXPORT void SetPararealIterationNumber (int num)
 
SOLVER_UTILS_EXPORT bool GetUseInitialCondition ()
 
SOLVER_UTILS_EXPORT void SetUseInitialCondition (bool num)
 
SOLVER_UTILS_EXPORT Array< OneD, const Array< OneD, NekDouble > > GetTraceNormals ()
 
SOLVER_UTILS_EXPORT void SetTime (const NekDouble time)
 
SOLVER_UTILS_EXPORT void SetTimeStep (const NekDouble timestep)
 
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...
 
SOLVER_UTILS_EXPORT bool ParallelInTime ()
 Check if solver use Parallel-in-Time. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::FluidInterface
virtual ~FluidInterface ()=default
 
SOLVER_UTILS_EXPORT void GetVelocity (const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, Array< OneD, NekDouble >> &velocity)
 Extract array with velocity from physfield. More...
 
SOLVER_UTILS_EXPORT bool HasConstantDensity ()
 
SOLVER_UTILS_EXPORT void GetDensity (const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &density)
 Extract array with density from physfield. More...
 
SOLVER_UTILS_EXPORT void GetPressure (const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &pressure)
 Extract array with pressure from physfield. More...
 
SOLVER_UTILS_EXPORT void SetMovingFrameVelocities (const Array< OneD, NekDouble > &vFrameVels)
 
SOLVER_UTILS_EXPORT void GetMovingFrameVelocities (Array< OneD, NekDouble > &vFrameVels)
 
SOLVER_UTILS_EXPORT void SetMovingFrameProjectionMat (const boost::numeric::ublas::matrix< NekDouble > &vProjMat)
 
SOLVER_UTILS_EXPORT void GetMovingFrameProjectionMat (boost::numeric::ublas::matrix< NekDouble > &vProjMat)
 
SOLVER_UTILS_EXPORT void SetMovingFrameAngles (const Array< OneD, NekDouble > &vFrameTheta)
 
SOLVER_UTILS_EXPORT void GetMovingFrameAngles (Array< OneD, NekDouble > &vFrameTheta)
 

Static Public Member Functions

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::SolverUtils::UnsteadySystem
static std::string cmdSetStartTime
 
static std::string cmdSetStartChkNum
 

Protected Member Functions

 NavierStokesCFE (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
void GetViscousFluxVectorConservVar (const size_t nDim, const Array< OneD, Array< OneD, NekDouble >> &inarray, const TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &outarray, Array< OneD, int > &nonZeroIndex=NullInt1DArray, const Array< OneD, Array< OneD, NekDouble >> &normal=NullNekDoubleArrayOfArray)
 
void GetViscousSymmtrFluxConservVar (const size_t nSpaceDim, const Array< OneD, Array< OneD, NekDouble >> &inaverg, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &outarray, Array< OneD, int > &nonZeroIndex, const Array< OneD, Array< OneD, NekDouble >> &normals)
 Calculate and return the Symmetric flux in IP method. More...
 
void SpecialBndTreat (Array< OneD, Array< OneD, NekDouble >> &consvar)
 For very special treatment. For general boundaries it does nothing But for WallViscous and WallAdiabatic, special treatment is needed because they get the same Bwd value, special treatment is needed for boundary treatment of diffusion flux Note: This special treatment could be removed by seperating WallViscous and WallAdiabatic into two different classes. More...
 
void GetArtificialViscosity (const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, NekDouble > &muav)
 
void CalcViscosity (const Array< OneD, const Array< OneD, NekDouble >> &inaverg, Array< OneD, NekDouble > &mu)
 
void InitObject_Explicit ()
 
virtual void v_InitObject (bool DeclareField=true) override
 Initialization object for CompressibleFlowSystem class. More...
 
virtual void v_DoDiffusion (const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd) override
 
virtual void v_GetViscousFluxVector (const Array< OneD, const Array< OneD, NekDouble >> &physfield, TensorOfArray3D< NekDouble > &derivatives, TensorOfArray3D< NekDouble > &viscousTensor)
 Return the flux vector for the LDG diffusion problem. More...
 
virtual void v_GetViscousFluxVectorDeAlias (const Array< OneD, const Array< OneD, NekDouble >> &physfield, TensorOfArray3D< NekDouble > &derivatives, TensorOfArray3D< NekDouble > &viscousTensor)
 Return the flux vector for the LDG diffusion problem. More...
 
void GetPhysicalAV (const Array< OneD, const Array< OneD, NekDouble >> &physfield)
 
void Ducros (Array< OneD, NekDouble > &field)
 
void C0Smooth (Array< OneD, NekDouble > &field)
 
virtual void v_GetFluxPenalty (const Array< OneD, const Array< OneD, NekDouble >> &uFwd, const Array< OneD, const 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...
 
void GetDivCurlSquared (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &cnsVar, Array< OneD, NekDouble > &div, Array< OneD, NekDouble > &curlSquare, const Array< OneD, Array< OneD, NekDouble >> &cnsVarFwd, const Array< OneD, Array< OneD, NekDouble >> &cnsVarBwd)
 Get divergence and curl squared. More...
 
void GetDivCurlFromDvelT (const TensorOfArray3D< NekDouble > &pVarDer, Array< OneD, NekDouble > &div, Array< OneD, NekDouble > &curlSquare)
 Get divergence and curl from velocity derivative tensor. More...
 
virtual void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble >> &fieldcoeffs, std::vector< std::string > &variables) override
 
template<class T , typename = typename std::enable_if< std::is_floating_point<T>::value || tinysimd::is_vector_floating_point<T>::value>::type>
void GetViscosityAndThermalCondFromTempKernel (const T &temperature, T &mu, T &thermalCond)
 
template<class T , typename = typename std::enable_if< std::is_floating_point<T>::value || tinysimd::is_vector_floating_point<T>::value>::type>
void GetViscosityFromTempKernel (const T &temperature, T &mu)
 
template<class T , typename = typename std::enable_if< std::is_floating_point<T>::value || tinysimd::is_vector_floating_point<T>::value>::type>
void GetViscousFluxBilinearFormKernel (const unsigned short nDim, const unsigned short FluxDirection, const unsigned short DerivDirection, const T *inaverg, const T *injumpp, const T &mu, T *outarray)
 Calculate diffusion flux using the Jacobian form. More...
 
template<bool IS_TRACE>
void GetViscousFluxVectorConservVar (const size_t nDim, const Array< OneD, Array< OneD, NekDouble >> &inarray, const TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &outarray, Array< OneD, int > &nonZeroIndex, const Array< OneD, Array< OneD, NekDouble >> &normal)
 Return the flux vector for the IP diffusion problem, based on conservative variables. More...
 
virtual bool v_SupportsShockCaptType (const std::string type) const override
 
- 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, 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, 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, const Array< OneD, NekDouble >> &physfield, TensorOfArray3D< NekDouble > &flux)
 Return the flux vector for the compressible Euler equations. More...
 
void GetFluxVectorDeAlias (const Array< OneD, const Array< OneD, NekDouble >> &physfield, TensorOfArray3D< 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 SetBoundaryConditionsBwdWeight ()
 Set up a weight on physical boundaries for boundary condition applications. More...
 
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) override
 Calculate the maximum timestep subject to CFL restrictions. More...
 
virtual void v_SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0) override
 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 MultiRegions::ExpListSharedPtr v_GetPressure () override
 
virtual Array< OneD, NekDoublev_GetMaxStdVelocity (const NekDouble SpeedSoundFactor) override
 Compute the advection velocity in the standard space for each element of the expansion. More...
 
virtual void v_SteadyStateResidual (int step, Array< OneD, NekDouble > &L2) override
 
- Protected Member Functions inherited from Nektar::SolverUtils::AdvectionSystem
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step) override
 
- 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 () override
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise () override
 Sets up initial conditions. More...
 
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary (SummaryList &s) override
 Print a summary of time stepping parameters. 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...
 
virtual SOLVER_UTILS_EXPORT bool v_UpdateTimeStepCheck ()
 
SOLVER_UTILS_EXPORT void DoDummyProjection (const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
 Perform dummy projection. 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)
 
- Protected Member Functions inherited from Nektar::SolverUtils::FluidInterface
virtual SOLVER_UTILS_EXPORT void v_SetMovingFrameVelocities (const Array< OneD, NekDouble > &vFrameVels)
 
virtual SOLVER_UTILS_EXPORT void v_GetMovingFrameVelocities (Array< OneD, NekDouble > &vFrameVels)
 
virtual SOLVER_UTILS_EXPORT void v_SetMovingFrameProjectionMat (const boost::numeric::ublas::matrix< NekDouble > &vProjMat)
 
virtual SOLVER_UTILS_EXPORT void v_GetMovingFrameProjectionMat (boost::numeric::ublas::matrix< NekDouble > &vProjMat)
 
virtual SOLVER_UTILS_EXPORT void v_SetMovingFrameAngles (const Array< OneD, NekDouble > &vFrameTheta)
 
virtual SOLVER_UTILS_EXPORT void v_GetMovingFrameAngles (Array< OneD, NekDouble > &vFrameTheta)
 

Protected Attributes

std::string m_ViscosityType
 
bool m_is_mu_variable {false}
 flag to switch between constant viscosity and Sutherland an enum could be added for more options More...
 
bool m_is_diffIP {false}
 flag to switch between IP and LDG an enum could be added for more options More...
 
bool m_is_shockCaptPhys {false}
 flag for shock capturing switch on/off an enum could be added for more options More...
 
NekDouble m_Cp
 
NekDouble m_Cv
 
NekDouble m_Prandtl
 
std::string m_physicalSensorType
 
std::string m_smoothing
 
MultiRegions::ContFieldSharedPtr m_C0ProjectExp
 
EquationOfStateSharedPtr m_eos
 Equation of system for computing temperature. More...
 
NekDouble m_Twall
 
NekDouble m_muRef
 
NekDouble m_thermalConductivityRef
 
- 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
 
Array< OneD, NekDoublem_muav
 
Array< OneD, NekDoublem_muavTrace
 
VariableConverterSharedPtr m_varConv
 
std::vector< CFSBndCondSharedPtrm_bndConds
 
NekDouble m_bndEvaluateTime
 
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_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::TimeIntegrationSchemeSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
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...
 
NekDouble m_steadyStateRes = 1.0
 
NekDouble m_steadyStateRes0 = 1.0
 
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...
 
NekDouble m_TimeIntegLambda = 0.0
 coefff of spacial derivatives(rhs or m_F in GLM) in calculating the residual of the whole equation(used in unsteady time integrations) More...
 
bool m_flagImplicitItsStatistics
 
bool m_flagImplicitSolver = false
 
Array< OneD, NekDoublem_magnitdEstimat
 estimate the magnitude of each conserved varibles More...
 
Array< OneD, NekDoublem_locTimeStep
 local time step(notice only for jfnk other see m_cflSafetyFactor) More...
 
NekDouble m_inArrayNorm = -1.0
 
int m_TotLinItePerStep = 0
 
int m_StagesPerStep = 1
 
bool m_flagUpdatePreconMat
 
int m_maxLinItePerNewton
 
int m_TotNewtonIts = 0
 
int m_TotLinIts = 0
 
int m_TotImpStages = 0
 
bool m_CalcPhysicalAV = true
 flag to update artificial viscosity More...
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
bool m_verbose
 
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_timestepMax = -1.0
 Time step size. More...
 
NekDouble m_lambda
 Lambda constant in real system if one required. More...
 
NekDouble m_checktime
 Time between checkpoints. More...
 
NekDouble m_lastCheckTime
 
NekDouble m_TimeIncrementFactor
 
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_infosteps
 Number of time steps between outputting status information. More...
 
int m_pararealIter
 Number of parareal time iteration. 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_useInitialCondition
 Flag to determine if IC are used. 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...
 
Array< OneD, NekDoublem_movingFrameVelsxyz
 Moving frame of reference velocities. More...
 
Array< OneD, NekDoublem_movingFrameTheta
 Moving frame of reference angles with respect to the. More...
 
boost::numeric::ublas::matrix< NekDoublem_movingFrameProjMat
 Projection matrix for transformation between inertial and moving. 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< NavierStokesCFE >
 

Additional Inherited Members

- Public Attributes inherited from Nektar::SolverUtils::UnsteadySystem
NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1). More...
 
NekDouble m_cflNonAcoustic
 
NekDouble m_CFLGrowth
 CFL growth rate. More...
 
NekDouble m_CFLEnd
 maximun cfl in cfl growth 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 49 of file NavierStokesCFE.h.

Constructor & Destructor Documentation

◆ ~NavierStokesCFE()

Nektar::NavierStokesCFE::~NavierStokesCFE ( )
virtual

Definition at line 53 of file NavierStokesCFE.cpp.

54 {
55 }

◆ NavierStokesCFE()

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

Definition at line 46 of file NavierStokesCFE.cpp.

49  : UnsteadySystem(pSession, pGraph), CompressibleFlowSystem(pSession, pGraph)
50 {
51 }
CompressibleFlowSystem(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

◆ C0Smooth()

void Nektar::NavierStokesCFE::C0Smooth ( Array< OneD, NekDouble > &  field)
protected

◆ CalcViscosity()

void Nektar::NavierStokesCFE::CalcViscosity ( const Array< OneD, const Array< OneD, NekDouble >> &  inaverg,
Array< OneD, NekDouble > &  mu 
)
protected

◆ create()

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

Definition at line 55 of file NavierStokesCFE.h.

58  {
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.

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

◆ Ducros()

void Nektar::NavierStokesCFE::Ducros ( Array< OneD, NekDouble > &  field)
protected

◆ GetArtificialViscosity()

void Nektar::NavierStokesCFE::GetArtificialViscosity ( const Array< OneD, Array< OneD, NekDouble >> &  inarray,
Array< OneD, NekDouble > &  muav 
)
protected

◆ GetDivCurlFromDvelT()

void Nektar::NavierStokesCFE::GetDivCurlFromDvelT ( const TensorOfArray3D< NekDouble > &  pVarDer,
Array< OneD, NekDouble > &  div,
Array< OneD, NekDouble > &  curlSquare 
)
protected

Get divergence and curl from velocity derivative tensor.

Definition at line 772 of file NavierStokesCFE.cpp.

775 {
776  size_t nDim = pVarDer.size();
777  size_t nPts = div.size();
778 
779  // div velocity
780  for (size_t p = 0; p < nPts; ++p)
781  {
782  NekDouble divTmp = 0;
783  for (unsigned short j = 0; j < nDim; ++j)
784  {
785  divTmp += pVarDer[j][j][p];
786  }
787  div[p] = divTmp;
788  }
789 
790  // |curl velocity| ** 2
791  if (nDim > 2)
792  {
793  for (size_t p = 0; p < nPts; ++p)
794  {
795  // curl[0] 3/2 - 2/3
796  NekDouble curl032 = pVarDer[2][1][p]; // load 1x
797  NekDouble curl023 = pVarDer[1][2][p]; // load 1x
798  NekDouble curl0 = curl032 - curl023;
799  // square curl[0]
800  NekDouble curl0sqr = curl0 * curl0;
801 
802  // curl[1] 3/1 - 1/3
803  NekDouble curl131 = pVarDer[2][0][p]; // load 1x
804  NekDouble curl113 = pVarDer[0][2][p]; // load 1x
805  NekDouble curl1 = curl131 - curl113;
806  // square curl[1]
807  NekDouble curl1sqr = curl1 * curl1;
808 
809  // curl[2] 1/2 - 2/1
810  NekDouble curl212 = pVarDer[0][1][p]; // load 1x
811  NekDouble curl221 = pVarDer[1][0][p]; // load 1x
812  NekDouble curl2 = curl212 - curl221;
813  // square curl[2]
814  NekDouble curl2sqr = curl2 * curl2;
815 
816  // reduce
817  curl0sqr += curl1sqr + curl2sqr;
818  // store
819  curlSquare[p] = curl0sqr; // store 1x
820  }
821  }
822  else if (nDim > 1)
823  {
824  for (size_t p = 0; p < nPts; ++p)
825  {
826  // curl[2] 1/2
827  NekDouble c212 = pVarDer[0][1][p]; // load 1x
828  // curl[2] 2/1
829  NekDouble c221 = pVarDer[1][0][p]; // load 1x
830  // curl[2] 1/2 - 2/1
831  NekDouble curl = c212 - c221;
832  // square curl[2]
833  curlSquare[p] = curl * curl; // store 1x
834  }
835  }
836  else
837  {
838  Vmath::Fill(nPts, 0.0, curlSquare, 1);
839  }
840 }
double NekDouble
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:45

References Vmath::Fill(), and CellMLToNektar.cellml_metadata::p.

Referenced by GetDivCurlSquared().

◆ GetDivCurlSquared()

void Nektar::NavierStokesCFE::GetDivCurlSquared ( const Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
const Array< OneD, Array< OneD, NekDouble >> &  cnsVar,
Array< OneD, NekDouble > &  div,
Array< OneD, NekDouble > &  curlSquare,
const Array< OneD, Array< OneD, NekDouble >> &  cnsVarFwd,
const Array< OneD, Array< OneD, NekDouble >> &  cnsVarBwd 
)
protected

Get divergence and curl squared.

Parameters
inputfields -> expansion list pointer cnsVar -> conservative variables cnsVarFwd -> forward trace of conservative variables cnsVarBwd -> backward trace of conservative variables @paran output divSquare -> divergence curlSquare -> curl squared

Definition at line 704 of file NavierStokesCFE.cpp.

710 {
711  size_t nDim = fields[0]->GetCoordim(0);
712  size_t nVar = cnsVar.size();
713  size_t nPts = cnsVar[0].size();
714  size_t nPtsTrc = cnsVarFwd[0].size();
715 
716  // These should be allocated once
717  Array<OneD, Array<OneD, NekDouble>> primVar(nVar - 1), primVarFwd(nVar - 1),
718  primVarBwd(nVar - 1);
719 
720  for (unsigned short d = 0; d < nVar - 2; ++d)
721  {
722  primVar[d] = Array<OneD, NekDouble>(nPts, 0.0);
723  primVarFwd[d] = Array<OneD, NekDouble>(nPtsTrc, 0.0);
724  primVarBwd[d] = Array<OneD, NekDouble>(nPtsTrc, 0.0);
725  }
726  size_t ergLoc = nVar - 2;
727  primVar[ergLoc] = Array<OneD, NekDouble>(nPts, 0.0);
728  primVarFwd[ergLoc] = Array<OneD, NekDouble>(nPtsTrc, 0.0);
729  primVarBwd[ergLoc] = Array<OneD, NekDouble>(nPtsTrc, 0.0);
730 
731  // Get primitive variables [u,v,w,T=0]
732  // Possibly should be changed to [rho,u,v,w,T] to make IP and LDGNS more
733  // consistent with each other
734  for (unsigned short d = 0; d < nVar - 2; ++d)
735  {
736  // Volume
737  for (size_t p = 0; p < nPts; ++p)
738  {
739  primVar[d][p] = cnsVar[d + 1][p] / cnsVar[0][p];
740  }
741  // Trace
742  for (size_t p = 0; p < nPtsTrc; ++p)
743  {
744  primVarFwd[d][p] = cnsVarFwd[d + 1][p] / cnsVarFwd[0][p];
745  primVarBwd[d][p] = cnsVarBwd[d + 1][p] / cnsVarBwd[0][p];
746  }
747  }
748 
749  // this should be allocated once
751  for (unsigned short j = 0; j < nDim; ++j)
752  {
753  primVarDer[j] = Array<OneD, Array<OneD, NekDouble>>(nVar - 1);
754  for (unsigned short i = 0; i < nVar - 1; ++i)
755  {
756  primVarDer[j][i] = Array<OneD, NekDouble>(nPts, 0.0);
757  }
758  }
759 
760  // Get derivative tensor
761  m_diffusion->DiffuseCalcDerivative(fields, primVar, primVarDer, primVarFwd,
762  primVarBwd);
763 
764  // Get div curl squared
765  GetDivCurlFromDvelT(primVarDer, div, curlSquare);
766 }
SolverUtils::DiffusionSharedPtr m_diffusion
void GetDivCurlFromDvelT(const TensorOfArray3D< NekDouble > &pVarDer, Array< OneD, NekDouble > &div, Array< OneD, NekDouble > &curlSquare)
Get divergence and curl from velocity derivative tensor.

References GetDivCurlFromDvelT(), Nektar::CompressibleFlowSystem::m_diffusion, and CellMLToNektar.cellml_metadata::p.

Referenced by v_DoDiffusion(), Nektar::NavierStokesImplicitCFE::v_DoDiffusionCoeff(), and v_ExtraFldOutput().

◆ GetPhysicalAV()

void Nektar::NavierStokesCFE::GetPhysicalAV ( const Array< OneD, const Array< OneD, NekDouble >> &  physfield)
protected

◆ GetViscosityAndThermalCondFromTemp()

void Nektar::NavierStokesCFE::GetViscosityAndThermalCondFromTemp ( const Array< OneD, NekDouble > &  temperature,
Array< OneD, NekDouble > &  mu,
Array< OneD, NekDouble > &  thermalCond 
)
protected

Update viscosity todo: add artificial viscosity here.

Definition at line 659 of file NavierStokesCFE.cpp.

662 {
663  size_t nPts = temperature.size();
664 
665  for (size_t p = 0; p < nPts; ++p)
666  {
667  GetViscosityAndThermalCondFromTempKernel(temperature[p], mu[p],
668  thermalCond[p]);
669  }
670 
671  // Add artificial viscosity if wanted
672  // move this above and add in kernel
673  if (m_is_shockCaptPhys)
674  {
675  size_t nTracePts = m_fields[0]->GetTrace()->GetTotPoints();
676  if (nPts != nTracePts)
677  {
678  Vmath::Vadd(nPts, mu, 1, m_varConv->GetAv(), 1, mu, 1);
679  }
680  else
681  {
682  Vmath::Vadd(nPts, mu, 1, m_varConv->GetAvTrace(), 1, mu, 1);
683  }
684 
685  // Update thermal conductivity
686  NekDouble tRa = m_Cp / m_Prandtl;
687  Vmath::Smul(nPts, tRa, mu, 1, thermalCond, 1);
688  }
689 }
VariableConverterSharedPtr m_varConv
bool m_is_shockCaptPhys
flag for shock capturing switch on/off an enum could be added for more options
void GetViscosityAndThermalCondFromTempKernel(const T &temperature, T &mu, T &thermalCond)
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:359
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:248

References GetViscosityAndThermalCondFromTempKernel(), m_Cp, Nektar::SolverUtils::EquationSystem::m_fields, m_is_shockCaptPhys, m_Prandtl, Nektar::CompressibleFlowSystem::m_varConv, CellMLToNektar.cellml_metadata::p, Vmath::Smul(), and Vmath::Vadd().

Referenced by GetViscousFluxVectorConservVar(), GetViscousSymmtrFluxConservVar(), Nektar::NavierStokesImplicitCFE::v_CalcMuDmuDT(), Nektar::NavierStokesImplicitCFE::v_GetFluxDerivJacDirctn(), v_GetFluxPenalty(), v_GetViscousFluxVector(), Nektar::NavierStokesCFEAxisym::v_GetViscousFluxVector(), and v_GetViscousFluxVectorDeAlias().

◆ GetViscosityAndThermalCondFromTempKernel()

template<class T , typename = typename std::enable_if< std::is_floating_point<T>::value || tinysimd::is_vector_floating_point<T>::value>::type>
void Nektar::NavierStokesCFE::GetViscosityAndThermalCondFromTempKernel ( const T &  temperature,
T &  mu,
T &  thermalCond 
)
inlineprotected

Definition at line 176 of file NavierStokesCFE.h.

178  {
179  GetViscosityFromTempKernel(temperature, mu);
180  NekDouble tRa = m_Cp / m_Prandtl;
181  thermalCond = tRa * mu;
182  }
void GetViscosityFromTempKernel(const T &temperature, T &mu)

References GetViscosityFromTempKernel(), m_Cp, and m_Prandtl.

Referenced by GetViscosityAndThermalCondFromTemp().

◆ GetViscosityFromTempKernel()

template<class T , typename = typename std::enable_if< std::is_floating_point<T>::value || tinysimd::is_vector_floating_point<T>::value>::type>
void Nektar::NavierStokesCFE::GetViscosityFromTempKernel ( const T &  temperature,
T &  mu 
)
inlineprotected

Definition at line 187 of file NavierStokesCFE.h.

188  {
189  // Variable viscosity through the Sutherland's law
190  if (m_is_mu_variable)
191  {
192  mu = m_varConv->GetDynamicViscosity(temperature);
193  }
194  else
195  {
196  mu = m_muRef;
197  }
198  }
bool m_is_mu_variable
flag to switch between constant viscosity and Sutherland an enum could be added for more options

References m_is_mu_variable, m_muRef, and Nektar::CompressibleFlowSystem::m_varConv.

Referenced by GetViscosityAndThermalCondFromTempKernel().

◆ GetViscousFluxBilinearFormKernel()

template<class T , typename = typename std::enable_if< std::is_floating_point<T>::value || tinysimd::is_vector_floating_point<T>::value>::type>
void Nektar::NavierStokesCFE::GetViscousFluxBilinearFormKernel ( const unsigned short  nDim,
const unsigned short  FluxDirection,
const unsigned short  DerivDirection,
const T *  inaverg,
const T *  injumpp,
const T &  mu,
T *  outarray 
)
inlineprotected

Calculate diffusion flux using the Jacobian form.

Parameters
in
outoutarray[nvars] flux

Definition at line 211 of file NavierStokesCFE.h.

216  {
217  // Constants
218  unsigned short nDim_plus_one = nDim + 1;
219  unsigned short FluxDirection_plus_one = FluxDirection + 1;
220  unsigned short DerivDirection_plus_one = DerivDirection + 1;
221 
222  NekDouble gammaoPr = m_gamma / m_Prandtl;
223  NekDouble one_minus_gammaoPr = 1.0 - gammaoPr;
224 
225  constexpr NekDouble OneThird = 1. / 3.;
226  constexpr NekDouble TwoThird = 2. * OneThird;
227  constexpr NekDouble FourThird = 4. * OneThird;
228 
229  if (DerivDirection == FluxDirection)
230  {
231  // rho flux always zero
232  outarray[0] = 0.0; // store 1x
233 
234  // load 1/rho
235  T oneOrho = 1.0 / inaverg[0]; // load 1x
236  // get vel, vel^2, sum of vel^2
237  std::array<T, 3> u = {{0.0, 0.0, 0.0}};
238  std::array<T, 3> u2 = {{0.0, 0.0, 0.0}};
239  T u2sum{};
240  for (unsigned short d = 0; d < nDim; ++d)
241  {
242  u[d] = inaverg[d + 1] * oneOrho; // load 1x
243  u2[d] = u[d] * u[d];
244  u2sum += u2[d];
245  }
246 
247  // get E - sum v^2
248  T E_minus_u2sum = inaverg[nDim_plus_one]; // load 1x
249  E_minus_u2sum *= oneOrho;
250  E_minus_u2sum -= u2sum;
251 
252  // get nu = mu/rho
253  T nu = mu * oneOrho; // load 1x
254 
255  // ^^^^ above is almost the same for both loops
256 
257  T tmp1 = OneThird * u2[FluxDirection] + u2sum;
258  tmp1 += gammaoPr * E_minus_u2sum;
259  tmp1 *= injumpp[0]; // load 1x
260 
261  T tmp2 = gammaoPr * injumpp[nDim_plus_one] - tmp1; // load 1x
262 
263  // local var for energy output
264  T outTmpE = 0.0;
265  for (unsigned short d = 0; d < nDim; ++d)
266  {
267  unsigned short d_plus_one = d + 1;
268  // flux[rhou, rhov, rhow]
269  T outTmpD = injumpp[d_plus_one] - u[d] * injumpp[0];
270  outTmpD *= nu;
271  // flux rhoE
272  T tmp3 = one_minus_gammaoPr * u[d];
273  outTmpE += tmp3 * injumpp[d_plus_one];
274 
275  if (d == FluxDirection)
276  {
277  outTmpD *= FourThird;
278  T tmp4 = OneThird * u[FluxDirection];
279  outTmpE += tmp4 * injumpp[FluxDirection_plus_one];
280  }
281 
282  outarray[d_plus_one] = outTmpD; // store 1x
283  }
284 
285  outTmpE += tmp2;
286  outTmpE *= nu;
287  outarray[nDim_plus_one] = outTmpE; // store 1x
288  }
289  else
290  {
291  // rho flux always zero
292  outarray[0] = 0.0; // store 1x
293 
294  // load 1/rho
295  T oneOrho = 1.0 / inaverg[0]; // load 1x
296  // get vel, vel^2, sum of vel^2
297  std::array<T, 3> u, u2;
298  T u2sum{};
299  for (unsigned short d = 0; d < nDim; ++d)
300  {
301  unsigned short d_plus_one = d + 1;
302  u[d] = inaverg[d_plus_one] * oneOrho; // load 1x
303  u2[d] = u[d] * u[d];
304  u2sum += u2[d];
305  // Not all directions are set
306  // one could work out the one that is not set
307  outarray[d_plus_one] = 0.0; // store 1x
308  }
309 
310  // get E - sum v^2
311  T E_minus_u2sum = inaverg[nDim_plus_one]; // load 1x
312  E_minus_u2sum *= oneOrho;
313  E_minus_u2sum -= u2sum;
314 
315  // get nu = mu/rho
316  T nu = mu * oneOrho; // load 1x
317 
318  // ^^^^ above is almost the same for both loops
319 
320  T tmp1 = u[DerivDirection] * injumpp[0] -
321  injumpp[DerivDirection_plus_one]; // load 2x
322  tmp1 *= TwoThird;
323  outarray[FluxDirection_plus_one] = nu * tmp1; // store 1x
324 
325  tmp1 =
326  injumpp[FluxDirection_plus_one] - u[FluxDirection] * injumpp[0];
327  outarray[DerivDirection_plus_one] = nu * tmp1; // store 1x
328 
329  tmp1 = OneThird * u[FluxDirection] * u[DerivDirection];
330  tmp1 *= injumpp[0];
331 
332  T tmp2 =
333  TwoThird * u[FluxDirection] * injumpp[DerivDirection_plus_one];
334 
335  tmp1 += tmp2;
336 
337  tmp1 = u[DerivDirection] * injumpp[FluxDirection_plus_one] - tmp1;
338  outarray[nDim_plus_one] = nu * tmp1; // store 1x
339  }
340  }

References Nektar::CompressibleFlowSystem::m_gamma, and m_Prandtl.

Referenced by GetViscousFluxVectorConservVar(), and GetViscousSymmtrFluxConservVar().

◆ GetViscousFluxVectorConservVar() [1/2]

template<bool IS_TRACE>
void Nektar::NavierStokesCFE::GetViscousFluxVectorConservVar ( const size_t  nDim,
const Array< OneD, Array< OneD, NekDouble >> &  inarray,
const TensorOfArray3D< NekDouble > &  qfields,
TensorOfArray3D< NekDouble > &  outarray,
Array< OneD, int > &  nonZeroIndex,
const Array< OneD, Array< OneD, NekDouble >> &  normal 
)
inlineprotected

Return the flux vector for the IP diffusion problem, based on conservative variables.

Definition at line 347 of file NavierStokesCFE.h.

352  {
353  size_t nConvectiveFields = inarray.size();
354  size_t nPts = inarray[0].size();
355  size_t n_nonZero = nConvectiveFields - 1;
356 
357  // max outfield dimensions
358  constexpr unsigned short nOutMax = 3 - 2 * IS_TRACE;
359  constexpr unsigned short nVarMax = 5;
360  constexpr unsigned short nDimMax = 3;
361 
362  // Update viscosity and thermal conductivity
363  // unfortunately the artificial viscosity is difficult to vectorize
364  // with the current implementation
365  Array<OneD, NekDouble> temperature(nPts, 0.0);
366  Array<OneD, NekDouble> mu(nPts, 0.0);
367  Array<OneD, NekDouble> thermalConductivity(nPts, 0.0);
368  m_varConv->GetTemperature(inarray, temperature);
369  GetViscosityAndThermalCondFromTemp(temperature, mu,
370  thermalConductivity);
371 
372  // vector loop
373  using namespace tinysimd;
374  using vec_t = simd<NekDouble>;
375  size_t sizeVec = (nPts / vec_t::width) * vec_t::width;
376  size_t p = 0;
377 
378  for (; p < sizeVec; p += vec_t::width)
379  {
380  // there is a significant penalty to use std::vector
381  alignas(vec_t::alignment) std::array<vec_t, nVarMax> inTmp,
382  qfieldsTmp, outTmp;
383  alignas(vec_t::alignment) std::array<vec_t, nDimMax> normalTmp;
384  alignas(vec_t::alignment) std::array<vec_t, nVarMax * nOutMax>
385  outArrTmp{{}};
386 
387  // rearrenge and load data
388  for (size_t f = 0; f < nConvectiveFields; ++f)
389  {
390  inTmp[f].load(&(inarray[f][p]), is_not_aligned);
391  // zero output vector
392  if (IS_TRACE)
393  {
394  outArrTmp[f] = 0.0;
395  }
396  else
397  {
398  for (size_t d = 0; d < nDim; ++d)
399  {
400  outArrTmp[f + nConvectiveFields * d] = 0.0;
401  }
402  }
403  }
404  if (IS_TRACE)
405  {
406  for (size_t d = 0; d < nDim; ++d)
407  {
408  normalTmp[d].load(&(normal[d][p]), is_not_aligned);
409  }
410  }
411 
412  // get viscosity
413  vec_t muV{};
414  muV.load(&(mu[p]), is_not_aligned);
415 
416  for (size_t nderiv = 0; nderiv < nDim; ++nderiv)
417  {
418  // rearrenge and load data
419  for (size_t f = 0; f < nConvectiveFields; ++f)
420  {
421  qfieldsTmp[f].load(&(qfields[nderiv][f][p]),
423  }
424 
425  for (size_t d = 0; d < nDim; ++d)
426  {
428  nDim, d, nderiv, inTmp.data(), qfieldsTmp.data(), muV,
429  outTmp.data());
430 
431  if (IS_TRACE)
432  {
433  for (size_t f = 0; f < nConvectiveFields; ++f)
434  {
435  outArrTmp[f] += normalTmp[d] * outTmp[f];
436  }
437  }
438  else
439  {
440  for (size_t f = 0; f < nConvectiveFields; ++f)
441  {
442  outArrTmp[f + nConvectiveFields * d] += outTmp[f];
443  }
444  }
445  }
446  }
447 
448  // store data
449  if (IS_TRACE)
450  {
451  for (size_t f = 0; f < nConvectiveFields; ++f)
452  {
453  outArrTmp[f].store(&(outarray[0][f][p]), is_not_aligned);
454  }
455  }
456  else
457  {
458  for (size_t d = 0; d < nDim; ++d)
459  {
460  for (size_t f = 0; f < nConvectiveFields; ++f)
461  {
462  outArrTmp[f + nConvectiveFields * d].store(
463  &(outarray[d][f][p]), is_not_aligned);
464  }
465  }
466  }
467  }
468 
469  // scalar loop
470  for (; p < nPts; ++p)
471  {
472  std::array<NekDouble, nVarMax> inTmp, qfieldsTmp, outTmp;
473  std::array<NekDouble, nDimMax> normalTmp;
474  std::array<NekDouble, nVarMax * nOutMax> outArrTmp{{}};
475  // rearrenge and load data
476  for (size_t f = 0; f < nConvectiveFields; ++f)
477  {
478  inTmp[f] = inarray[f][p];
479  // zero output vector
480  if (IS_TRACE)
481  {
482  outArrTmp[f] = 0.0;
483  }
484  else
485  {
486  for (size_t d = 0; d < nDim; ++d)
487  {
488  outArrTmp[f + nConvectiveFields * d] = 0.0;
489  }
490  }
491  }
492 
493  if (IS_TRACE)
494  {
495  for (size_t d = 0; d < nDim; ++d)
496  {
497  normalTmp[d] = normal[d][p];
498  }
499  }
500 
501  // get viscosity
502  NekDouble muS = mu[p];
503 
504  for (size_t nderiv = 0; nderiv < nDim; ++nderiv)
505  {
506  // rearrenge and load data
507  for (size_t f = 0; f < nConvectiveFields; ++f)
508  {
509  qfieldsTmp[f] = qfields[nderiv][f][p];
510  }
511 
512  for (size_t d = 0; d < nDim; ++d)
513  {
515  nDim, d, nderiv, inTmp.data(), qfieldsTmp.data(), muS,
516  outTmp.data());
517 
518  if (IS_TRACE)
519  {
520  for (size_t f = 0; f < nConvectiveFields; ++f)
521  {
522  outArrTmp[f] += normalTmp[d] * outTmp[f];
523  }
524  }
525  else
526  {
527  for (size_t f = 0; f < nConvectiveFields; ++f)
528  {
529  outArrTmp[f + nConvectiveFields * d] += outTmp[f];
530  }
531  }
532  }
533  }
534 
535  // store data
536  if (IS_TRACE)
537  {
538  for (size_t f = 0; f < nConvectiveFields; ++f)
539  {
540  outarray[0][f][p] = outArrTmp[f];
541  }
542  }
543  else
544  {
545  for (size_t d = 0; d < nDim; ++d)
546  {
547  for (size_t f = 0; f < nConvectiveFields; ++f)
548  {
549  outarray[d][f][p] =
550  outArrTmp[f + nConvectiveFields * d];
551  }
552  }
553  }
554  }
555 
556  // this is always the same, it should be initialized with the IP class
557  nonZeroIndex = Array<OneD, int>{n_nonZero, 0};
558  for (int i = 1; i < n_nonZero + 1; ++i)
559  {
560  nonZeroIndex[n_nonZero - i] = nConvectiveFields - i;
561  }
562  }
void GetViscosityAndThermalCondFromTemp(const Array< OneD, NekDouble > &temperature, Array< OneD, NekDouble > &mu, Array< OneD, NekDouble > &thermalCond)
Update viscosity todo: add artificial viscosity here.
void GetViscousFluxBilinearFormKernel(const unsigned short nDim, const unsigned short FluxDirection, const unsigned short DerivDirection, const T *inaverg, const T *injumpp, const T &mu, T *outarray)
Calculate diffusion flux using the Jacobian form.
tinysimd::simd< NekDouble > vec_t
static constexpr struct tinysimd::is_not_aligned_t is_not_aligned
typename abi< ScalarType, width >::type simd
Definition: tinysimd.hpp:80

References GetViscosityAndThermalCondFromTemp(), GetViscousFluxBilinearFormKernel(), tinysimd::is_not_aligned, Nektar::CompressibleFlowSystem::m_varConv, and CellMLToNektar.cellml_metadata::p.

◆ GetViscousFluxVectorConservVar() [2/2]

void Nektar::NavierStokesCFE::GetViscousFluxVectorConservVar ( const size_t  nDim,
const Array< OneD, Array< OneD, NekDouble >> &  inarray,
const TensorOfArray3D< NekDouble > &  qfields,
TensorOfArray3D< NekDouble > &  outarray,
Array< OneD, int > &  nonZeroIndex = NullInt1DArray,
const Array< OneD, Array< OneD, NekDouble >> &  normal = NullNekDoubleArrayOfArray 
)
protected

◆ GetViscousSymmtrFluxConservVar()

void Nektar::NavierStokesCFE::GetViscousSymmtrFluxConservVar ( const size_t  nSpaceDim,
const Array< OneD, Array< OneD, NekDouble >> &  inaverg,
const Array< OneD, Array< OneD, NekDouble >> &  inarray,
TensorOfArray3D< NekDouble > &  outarray,
Array< OneD, int > &  nonZeroIndex,
const Array< OneD, Array< OneD, NekDouble >> &  normals 
)
protected

Calculate and return the Symmetric flux in IP method.

Definition at line 565 of file NavierStokesCFE.cpp.

570 {
571  size_t nConvectiveFields = inarray.size();
572  size_t nPts = inaverg[nConvectiveFields - 1].size();
573  nonZeroIndex = Array<OneD, int>{nConvectiveFields - 1, 0};
574  for (size_t i = 0; i < nConvectiveFields - 1; ++i)
575  {
576  nonZeroIndex[i] = i + 1;
577  }
578 
579  Array<OneD, NekDouble> mu(nPts, 0.0);
580  Array<OneD, NekDouble> thermalConductivity(nPts, 0.0);
581  Array<OneD, NekDouble> temperature(nPts, 0.0);
582  m_varConv->GetTemperature(inaverg, temperature);
583  GetViscosityAndThermalCondFromTemp(temperature, mu, thermalConductivity);
584 
585  std::vector<NekDouble> inAvgTmp(nConvectiveFields);
586  std::vector<NekDouble> inTmp(nConvectiveFields);
587  std::vector<NekDouble> outTmp(nConvectiveFields);
588  for (size_t d = 0; d < nDim; ++d)
589  {
590  for (size_t nderiv = 0; nderiv < nDim; ++nderiv)
591  {
592  for (size_t p = 0; p < nPts; ++p)
593  {
594  // rearrenge data
595  for (size_t f = 0; f < nConvectiveFields; ++f)
596  {
597  inAvgTmp[f] = inaverg[f][p];
598  inTmp[f] = inarray[f][p];
599  }
600 
601  GetViscousFluxBilinearFormKernel(nDim, d, nderiv,
602  inAvgTmp.data(), inTmp.data(),
603  mu[p], outTmp.data());
604 
605  for (size_t f = 0; f < nConvectiveFields; ++f)
606  {
607  outarray[d][f][p] += normal[nderiv][p] * outTmp[f];
608  }
609  }
610  }
611  }
612 }

References GetViscosityAndThermalCondFromTemp(), GetViscousFluxBilinearFormKernel(), Nektar::CompressibleFlowSystem::m_varConv, and CellMLToNektar.cellml_metadata::p.

Referenced by InitObject_Explicit().

◆ InitObject_Explicit()

void Nektar::NavierStokesCFE::InitObject_Explicit ( )
protected

Definition at line 69 of file NavierStokesCFE.cpp.

70 {
71  // Get gas constant from session file and compute Cp
72  NekDouble gasConstant;
73  m_session->LoadParameter("GasConstant", gasConstant, 287.058);
74  m_Cp = m_gamma / (m_gamma - 1.0) * gasConstant;
75  m_Cv = m_Cp / m_gamma;
76 
77  m_session->LoadParameter("Twall", m_Twall, 300.15);
78 
79  // Viscosity
80  m_session->LoadSolverInfo("ViscosityType", m_ViscosityType, "Constant");
81  m_session->LoadParameter("mu", m_muRef, 1.78e-05);
82  if (m_ViscosityType == "Variable")
83  {
84  m_is_mu_variable = true;
85  }
86 
87  // Thermal conductivity or Prandtl
88  if (m_session->DefinesParameter("thermalConductivity"))
89  {
90  ASSERTL0(!m_session->DefinesParameter("Pr"),
91  "Cannot define both Pr and thermalConductivity.");
92 
93  m_session->LoadParameter("thermalConductivity",
96  }
97  else
98  {
99  m_session->LoadParameter("Pr", m_Prandtl, 0.72);
101  }
102 
103  if (m_shockCaptureType == "Physical")
104  {
105  m_is_shockCaptPhys = true;
106  }
107 
108  string diffName;
109  m_session->LoadSolverInfo("DiffusionType", diffName, "LDGNS");
110 
111  m_diffusion =
112  SolverUtils::GetDiffusionFactory().CreateInstance(diffName, diffName);
113  if ("InteriorPenalty" == diffName)
114  {
115  m_is_diffIP = true;
117  }
118 
119  if ("LDGNS" == diffName || "LDGNS3DHomogeneous1D" == diffName)
120  {
121  m_diffusion->SetFluxPenaltyNS(&NavierStokesCFE::v_GetFluxPenalty, this);
122  }
123 
125  {
126  m_diffusion->SetFluxVectorNS(
128  }
129  else
130  {
132  this);
133  }
134 
135  m_diffusion->SetDiffusionFluxCons(
136  &NavierStokesCFE::GetViscousFluxVectorConservVar<false>, this);
137 
138  m_diffusion->SetDiffusionFluxConsTrace(
139  &NavierStokesCFE::GetViscousFluxVectorConservVar<true>, this);
140 
141  m_diffusion->SetSpecialBndTreat(&NavierStokesCFE::SpecialBndTreat, this);
142 
143  m_diffusion->SetDiffusionSymmFluxCons(
145 
146  // Concluding initialisation of diffusion operator
147  m_diffusion->InitObject(m_session, m_fields);
148 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
void SetBoundaryConditionsBwdWeight()
Set up a weight on physical boundaries for boundary condition applications.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
virtual void v_GetViscousFluxVector(const Array< OneD, const Array< OneD, NekDouble >> &physfield, TensorOfArray3D< NekDouble > &derivatives, TensorOfArray3D< NekDouble > &viscousTensor)
Return the flux vector for the LDG diffusion problem.
void SpecialBndTreat(Array< OneD, Array< OneD, NekDouble >> &consvar)
For very special treatment. For general boundaries it does nothing But for WallViscous and WallAdiaba...
void GetViscousSymmtrFluxConservVar(const size_t nSpaceDim, const Array< OneD, Array< OneD, NekDouble >> &inaverg, const Array< OneD, Array< OneD, NekDouble >> &inarray, TensorOfArray3D< NekDouble > &outarray, Array< OneD, int > &nonZeroIndex, const Array< OneD, Array< OneD, NekDouble >> &normals)
Calculate and return the Symmetric flux in IP method.
NekDouble m_thermalConductivityRef
virtual void v_GetViscousFluxVectorDeAlias(const Array< OneD, const Array< OneD, NekDouble >> &physfield, TensorOfArray3D< NekDouble > &derivatives, TensorOfArray3D< NekDouble > &viscousTensor)
Return the flux vector for the LDG diffusion problem.
bool m_is_diffIP
flag to switch between IP and LDG an enum could be added for more options
virtual void v_GetFluxPenalty(const Array< OneD, const Array< OneD, NekDouble >> &uFwd, const Array< OneD, const Array< OneD, NekDouble >> &uBwd, Array< OneD, Array< OneD, NekDouble >> &penaltyCoeff)
Return the penalty vector for the LDGNS diffusion problem.
bool m_specHP_dealiasing
Flag to determine if dealisising is usde for the Spectral/hp element discretisation.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:41

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::SolverUtils::GetDiffusionFactory(), GetViscousSymmtrFluxConservVar(), m_Cp, m_Cv, Nektar::CompressibleFlowSystem::m_diffusion, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::CompressibleFlowSystem::m_gamma, m_is_diffIP, m_is_mu_variable, m_is_shockCaptPhys, m_muRef, m_Prandtl, Nektar::SolverUtils::EquationSystem::m_session, Nektar::CompressibleFlowSystem::m_shockCaptureType, Nektar::SolverUtils::EquationSystem::m_specHP_dealiasing, m_thermalConductivityRef, m_Twall, m_ViscosityType, Nektar::CompressibleFlowSystem::SetBoundaryConditionsBwdWeight(), SpecialBndTreat(), v_GetFluxPenalty(), v_GetViscousFluxVector(), and v_GetViscousFluxVectorDeAlias().

Referenced by v_InitObject(), and Nektar::NavierStokesImplicitCFE::v_InitObject().

◆ SpecialBndTreat()

void Nektar::NavierStokesCFE::SpecialBndTreat ( Array< OneD, Array< OneD, NekDouble >> &  consvar)
protected

For very special treatment. For general boundaries it does nothing But for WallViscous and WallAdiabatic, special treatment is needed because they get the same Bwd value, special treatment is needed for boundary treatment of diffusion flux Note: This special treatment could be removed by seperating WallViscous and WallAdiabatic into two different classes.

Definition at line 476 of file NavierStokesCFE.cpp.

478 {
479  size_t nConvectiveFields = consvar.size();
480  size_t ndens = 0;
481  size_t nengy = nConvectiveFields - 1;
482 
483  Array<OneD, Array<OneD, NekDouble>> bndCons{nConvectiveFields};
484  Array<OneD, NekDouble> bndTotEngy;
485  Array<OneD, NekDouble> bndPressure;
486  Array<OneD, NekDouble> bndRho;
487  Array<OneD, NekDouble> bndIntEndy;
488  size_t nLengthArray = 0;
489 
490  size_t cnt = 0;
491  size_t nBndRegions = m_fields[nengy]->GetBndCondExpansions().size();
492  for (size_t j = 0; j < nBndRegions; ++j)
493  {
494  if (m_fields[nengy]
495  ->GetBndConditions()[j]
496  ->GetBoundaryConditionType() == SpatialDomains::ePeriodic)
497  {
498  continue;
499  }
500 
501  size_t nBndEdges =
502  m_fields[nengy]->GetBndCondExpansions()[j]->GetExpSize();
503  for (size_t e = 0; e < nBndEdges; ++e)
504  {
505  size_t nBndEdgePts = m_fields[nengy]
506  ->GetBndCondExpansions()[j]
507  ->GetExp(e)
508  ->GetTotPoints();
509 
510  int id2 = m_fields[0]->GetTrace()->GetPhys_Offset(
511  m_fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt++));
512 
513  // Imposing Temperature Twall at the wall
514  if (boost::iequals(
515  m_fields[nengy]->GetBndConditions()[j]->GetUserDefined(),
516  "WallViscous"))
517  {
518  if (nBndEdgePts != nLengthArray)
519  {
520  for (size_t i = 0; i < nConvectiveFields; ++i)
521  {
522  bndCons[i] = Array<OneD, NekDouble>{nBndEdgePts, 0.0};
523  }
524  bndTotEngy = Array<OneD, NekDouble>{nBndEdgePts, 0.0};
525  bndPressure = Array<OneD, NekDouble>{nBndEdgePts, 0.0};
526  bndRho = Array<OneD, NekDouble>{nBndEdgePts, 0.0};
527  bndIntEndy = Array<OneD, NekDouble>{nBndEdgePts, 0.0};
528  nLengthArray = nBndEdgePts;
529  }
530  else
531  {
532  Vmath::Zero(nLengthArray, bndPressure, 1);
533  Vmath::Zero(nLengthArray, bndRho, 1);
534  Vmath::Zero(nLengthArray, bndIntEndy, 1);
535  }
536 
538 
539  for (size_t k = 0; k < nConvectiveFields; ++k)
540  {
541  Vmath::Vcopy(nBndEdgePts, tmp = consvar[k] + id2, 1,
542  bndCons[k], 1);
543  }
544 
545  m_varConv->GetPressure(bndCons, bndPressure);
546  Vmath::Fill(nLengthArray, m_Twall, bndTotEngy, 1);
547  m_varConv->GetRhoFromPT(bndPressure, bndTotEngy, bndRho);
548  m_varConv->GetEFromRhoP(bndRho, bndPressure, bndIntEndy);
549  m_varConv->GetDynamicEnergy(bndCons, bndTotEngy);
550 
551  Vmath::Vvtvp(nBndEdgePts, bndIntEndy, 1, bndCons[ndens], 1,
552  bndTotEngy, 1, bndTotEngy, 1);
553 
554  Vmath::Vcopy(nBndEdgePts, bndTotEngy, 1,
555  tmp = consvar[nengy] + id2, 1);
556  }
557  }
558  }
559 }
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:574
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255

References Nektar::SpatialDomains::ePeriodic, Vmath::Fill(), Nektar::SolverUtils::EquationSystem::m_fields, m_Twall, Nektar::CompressibleFlowSystem::m_varConv, Vmath::Vcopy(), Vmath::Vvtvp(), and Vmath::Zero().

Referenced by InitObject_Explicit().

◆ v_DoDiffusion()

void Nektar::NavierStokesCFE::v_DoDiffusion ( const Array< OneD, Array< OneD, NekDouble >> &  inarray,
Array< OneD, Array< OneD, NekDouble >> &  outarray,
const Array< OneD, Array< OneD, NekDouble >> &  pFwd,
const Array< OneD, Array< OneD, NekDouble >> &  pBwd 
)
overrideprotectedvirtual

Implements Nektar::CompressibleFlowSystem.

Reimplemented in Nektar::NavierStokesImplicitCFE, and Nektar::NavierStokesCFEAxisym.

Definition at line 150 of file NavierStokesCFE.cpp.

155 {
156  size_t nvariables = inarray.size();
157  size_t npoints = GetNpoints();
158  size_t nTracePts = GetTraceTotPoints();
159 
160  // this should be preallocated
161  Array<OneD, Array<OneD, NekDouble>> outarrayDiff(nvariables);
162  for (size_t i = 0; i < nvariables; ++i)
163  {
164  outarrayDiff[i] = Array<OneD, NekDouble>(npoints, 0.0);
165  }
166 
167  // Set artificial viscosity based on NS viscous tensor
168  if (m_is_shockCaptPhys)
169  {
170  if (m_varConv->GetFlagCalcDivCurl())
171  {
172  Array<OneD, NekDouble> div(npoints), curlSquare(npoints);
173  GetDivCurlSquared(m_fields, inarray, div, curlSquare, pFwd, pBwd);
174 
175  // Set volume and trace artificial viscosity
176  m_varConv->SetAv(m_fields, inarray, div, curlSquare);
177  }
178  else
179  {
180  m_varConv->SetAv(m_fields, inarray);
181  }
182  }
183 
184  if (m_is_diffIP)
185  {
186  if (m_bndEvaluateTime < 0.0)
187  {
188  NEKERROR(ErrorUtil::efatal, "m_bndEvaluateTime not setup");
189  }
190  m_diffusion->Diffuse(nvariables, m_fields, inarray, outarrayDiff,
191  m_bndEvaluateTime, pFwd, pBwd);
192  for (size_t i = 0; i < nvariables; ++i)
193  {
194  Vmath::Vadd(npoints, outarrayDiff[i], 1, outarray[i], 1,
195  outarray[i], 1);
196  }
197  }
198  else
199  {
200  // Get primitive variables [u,v,w,T]
201  Array<OneD, Array<OneD, NekDouble>> inarrayDiff(nvariables - 1);
202  Array<OneD, Array<OneD, NekDouble>> inFwd(nvariables - 1);
203  Array<OneD, Array<OneD, NekDouble>> inBwd(nvariables - 1);
204 
205  for (size_t i = 0; i < nvariables - 1; ++i)
206  {
207  inarrayDiff[i] = Array<OneD, NekDouble>{npoints};
208  inFwd[i] = Array<OneD, NekDouble>{nTracePts};
209  inBwd[i] = Array<OneD, NekDouble>{nTracePts};
210  }
211 
212  // Extract pressure
213  // (use inarrayDiff[0] as a temporary storage for the pressure)
214  m_varConv->GetPressure(inarray, inarrayDiff[0]);
215 
216  // Extract temperature
217  m_varConv->GetTemperature(inarray, inarrayDiff[nvariables - 2]);
218 
219  // Extract velocities
220  m_varConv->GetVelocityVector(inarray, inarrayDiff);
221 
222  // Repeat calculation for trace space
223  if (pFwd == NullNekDoubleArrayOfArray ||
225  {
228  }
229  else
230  {
231  m_varConv->GetPressure(pFwd, inFwd[0]);
232  m_varConv->GetPressure(pBwd, inBwd[0]);
233 
234  m_varConv->GetTemperature(pFwd, inFwd[nvariables - 2]);
235  m_varConv->GetTemperature(pBwd, inBwd[nvariables - 2]);
236 
237  m_varConv->GetVelocityVector(pFwd, inFwd);
238  m_varConv->GetVelocityVector(pBwd, inBwd);
239  }
240 
241  // Diffusion term in physical rhs form
242  m_diffusion->Diffuse(nvariables, m_fields, inarrayDiff, outarrayDiff,
243  inFwd, inBwd);
244 
245  for (size_t i = 0; i < nvariables; ++i)
246  {
247  Vmath::Vadd(npoints, outarrayDiff[i], 1, outarray[i], 1,
248  outarray[i], 1);
249  }
250  }
251 
252  // Add artificial diffusion through Laplacian operator
254  {
255  m_artificialDiffusion->DoArtificialDiffusion(inarray, outarray);
256  }
257 }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
ArtificialDiffusionSharedPtr m_artificialDiffusion
void GetDivCurlSquared(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble >> &cnsVar, Array< OneD, NekDouble > &div, Array< OneD, NekDouble > &curlSquare, const Array< OneD, Array< OneD, NekDouble >> &cnsVarFwd, const Array< OneD, Array< OneD, NekDouble >> &cnsVarBwd)
Get divergence and curl squared.
SOLVER_UTILS_EXPORT int GetTraceTotPoints()
SOLVER_UTILS_EXPORT int GetNpoints()
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray

References Nektar::ErrorUtil::efatal, GetDivCurlSquared(), Nektar::SolverUtils::EquationSystem::GetNpoints(), Nektar::SolverUtils::EquationSystem::GetTraceTotPoints(), Nektar::CompressibleFlowSystem::m_artificialDiffusion, Nektar::CompressibleFlowSystem::m_bndEvaluateTime, Nektar::CompressibleFlowSystem::m_diffusion, Nektar::SolverUtils::EquationSystem::m_fields, m_is_diffIP, m_is_shockCaptPhys, Nektar::CompressibleFlowSystem::m_varConv, NEKERROR, Nektar::NullNekDoubleArrayOfArray, and Vmath::Vadd().

Referenced by Nektar::NavierStokesCFEAxisym::v_DoDiffusion().

◆ v_ExtraFldOutput()

void Nektar::NavierStokesCFE::v_ExtraFldOutput ( std::vector< Array< OneD, NekDouble >> &  fieldcoeffs,
std::vector< std::string > &  variables 
)
overrideprotectedvirtual

Reimplemented from Nektar::CompressibleFlowSystem.

Definition at line 842 of file NavierStokesCFE.cpp.

845 {
846  bool extraFields;
847  m_session->MatchSolverInfo("OutputExtraFields", "True", extraFields, true);
848  if (extraFields)
849  {
850  const int nPhys = m_fields[0]->GetNpoints();
851  const int nCoeffs = m_fields[0]->GetNcoeffs();
853 
854  for (size_t i = 0; i < m_fields.size(); ++i)
855  {
856  cnsVar[i] = m_fields[i]->GetPhys();
857  }
858 
861  for (int i = 0; i < m_spacedim; ++i)
862  {
863  velocity[i] = Array<OneD, NekDouble>(nPhys);
864  velFwd[i] = Array<OneD, NekDouble>(nCoeffs);
865  }
866 
867  Array<OneD, NekDouble> pressure(nPhys), temperature(nPhys);
868  Array<OneD, NekDouble> entropy(nPhys);
869  Array<OneD, NekDouble> soundspeed(nPhys), mach(nPhys);
870  Array<OneD, NekDouble> sensor(nPhys), SensorKappa(nPhys);
871 
872  m_varConv->GetVelocityVector(cnsVar, velocity);
873  m_varConv->GetPressure(cnsVar, pressure);
874  m_varConv->GetTemperature(cnsVar, temperature);
875  m_varConv->GetEntropy(cnsVar, entropy);
876  m_varConv->GetSoundSpeed(cnsVar, soundspeed);
877  m_varConv->GetMach(cnsVar, soundspeed, mach);
878 
879  int sensorOffset;
880  m_session->LoadParameter("SensorOffset", sensorOffset, 1);
881  m_varConv->GetSensor(m_fields[0], cnsVar, sensor, SensorKappa,
882  sensorOffset);
883 
884  Array<OneD, NekDouble> pFwd(nCoeffs), TFwd(nCoeffs);
885  Array<OneD, NekDouble> sFwd(nCoeffs);
886  Array<OneD, NekDouble> aFwd(nCoeffs), mFwd(nCoeffs);
887  Array<OneD, NekDouble> sensFwd(nCoeffs);
888 
889  string velNames[3] = {"u", "v", "w"};
890  for (int i = 0; i < m_spacedim; ++i)
891  {
892  m_fields[0]->FwdTransLocalElmt(velocity[i], velFwd[i]);
893  variables.push_back(velNames[i]);
894  fieldcoeffs.push_back(velFwd[i]);
895  }
896 
897  m_fields[0]->FwdTransLocalElmt(pressure, pFwd);
898  m_fields[0]->FwdTransLocalElmt(temperature, TFwd);
899  m_fields[0]->FwdTransLocalElmt(entropy, sFwd);
900  m_fields[0]->FwdTransLocalElmt(soundspeed, aFwd);
901  m_fields[0]->FwdTransLocalElmt(mach, mFwd);
902  m_fields[0]->FwdTransLocalElmt(sensor, sensFwd);
903 
904  variables.push_back("p");
905  variables.push_back("T");
906  variables.push_back("s");
907  variables.push_back("a");
908  variables.push_back("Mach");
909  variables.push_back("Sensor");
910  fieldcoeffs.push_back(pFwd);
911  fieldcoeffs.push_back(TFwd);
912  fieldcoeffs.push_back(sFwd);
913  fieldcoeffs.push_back(aFwd);
914  fieldcoeffs.push_back(mFwd);
915  fieldcoeffs.push_back(sensFwd);
916 
918  {
919  // reuse pressure
920  Array<OneD, NekDouble> sensorFwd(nCoeffs);
921  m_artificialDiffusion->GetArtificialViscosity(cnsVar, pressure);
922  m_fields[0]->FwdTransLocalElmt(pressure, sensorFwd);
923 
924  variables.push_back("ArtificialVisc");
925  fieldcoeffs.push_back(sensorFwd);
926  }
927 
928  if (m_is_shockCaptPhys)
929  {
930 
931  Array<OneD, Array<OneD, NekDouble>> cnsVarFwd(m_fields.size()),
932  cnsVarBwd(m_fields.size());
933 
934  for (size_t i = 0; i < m_fields.size(); ++i)
935  {
936  cnsVarFwd[i] = Array<OneD, NekDouble>(GetTraceTotPoints());
937  cnsVarBwd[i] = Array<OneD, NekDouble>(GetTraceTotPoints());
938  m_fields[i]->GetFwdBwdTracePhys(cnsVar[i], cnsVarFwd[i],
939  cnsVarBwd[i]);
940  }
941 
942  Array<OneD, NekDouble> div(nPhys), curlSquare(nPhys);
943  GetDivCurlSquared(m_fields, cnsVar, div, curlSquare, cnsVarFwd,
944  cnsVarBwd);
945 
946  Array<OneD, NekDouble> divFwd(nCoeffs, 0.0);
947  m_fields[0]->FwdTransLocalElmt(div, divFwd);
948  variables.push_back("div");
949  fieldcoeffs.push_back(divFwd);
950 
951  Array<OneD, NekDouble> curlFwd(nCoeffs, 0.0);
952  m_fields[0]->FwdTransLocalElmt(curlSquare, curlFwd);
953  variables.push_back("curl^2");
954  fieldcoeffs.push_back(curlFwd);
955 
956  m_varConv->SetAv(m_fields, cnsVar, div, curlSquare);
957 
958  Array<OneD, NekDouble> muavFwd(nCoeffs);
959  m_fields[0]->FwdTransLocalElmt(m_varConv->GetAv(), muavFwd);
960  variables.push_back("ArtificialVisc");
961  fieldcoeffs.push_back(muavFwd);
962  }
963  }
964 }
int m_spacedim
Spatial dimension (>= expansion dim).

References GetDivCurlSquared(), Nektar::SolverUtils::EquationSystem::GetTraceTotPoints(), Nektar::CompressibleFlowSystem::m_artificialDiffusion, Nektar::SolverUtils::EquationSystem::m_fields, m_is_shockCaptPhys, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_spacedim, Nektar::CompressibleFlowSystem::m_varConv, and CG_Iterations::pressure.

◆ v_GetFluxPenalty()

void Nektar::NavierStokesCFE::v_GetFluxPenalty ( const Array< OneD, const Array< OneD, NekDouble >> &  uFwd,
const Array< OneD, const Array< OneD, NekDouble >> &  uBwd,
Array< OneD, Array< OneD, NekDouble >> &  penaltyCoeff 
)
protectedvirtual

Return the penalty vector for the LDGNS diffusion problem.

Definition at line 617 of file NavierStokesCFE.cpp.

621 {
622  size_t nTracePts = uFwd[0].size();
623 
624  // Compute average temperature
625  size_t nVariables = uFwd.size();
626  Array<OneD, NekDouble> tAve{nTracePts, 0.0};
627  Vmath::Svtsvtp(nTracePts, 0.5, uFwd[nVariables - 1], 1, 0.5,
628  uBwd[nVariables - 1], 1, tAve, 1);
629 
630  // Get average viscosity and thermal conductivity
631  Array<OneD, NekDouble> muAve{nTracePts, 0.0};
632  Array<OneD, NekDouble> tcAve{nTracePts, 0.0};
633 
634  GetViscosityAndThermalCondFromTemp(tAve, muAve, tcAve);
635 
636  // Compute penalty term
637  for (size_t i = 0; i < nVariables; ++i)
638  {
639  // Get jump of u variables
640  Vmath::Vsub(nTracePts, uFwd[i], 1, uBwd[i], 1, penaltyCoeff[i], 1);
641  // Multiply by variable coefficient = {coeff} ( u^+ - u^- )
642  if (i < nVariables - 1)
643  {
644  Vmath::Vmul(nTracePts, muAve, 1, penaltyCoeff[i], 1,
645  penaltyCoeff[i], 1);
646  }
647  else
648  {
649  Vmath::Vmul(nTracePts, tcAve, 1, penaltyCoeff[i], 1,
650  penaltyCoeff[i], 1);
651  }
652  }
653 }
void Svtsvtp(int n, const T alpha, const T *x, int incx, const T beta, const T *y, int incy, T *z, int incz)
svtvvtp (scalar times vector plus scalar times vector):
Definition: Vmath.cpp:751
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:209
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:419

References GetViscosityAndThermalCondFromTemp(), Vmath::Svtsvtp(), Vmath::Vmul(), and Vmath::Vsub().

Referenced by InitObject_Explicit().

◆ v_GetViscousFluxVector()

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

Return the flux vector for the LDG diffusion problem.

Todo:
Complete the viscous flux vector

Reimplemented in Nektar::NavierStokesCFEAxisym.

Definition at line 263 of file NavierStokesCFE.cpp.

267 {
268  // Auxiliary variables
269  size_t nScalar = physfield.size();
270  size_t nPts = physfield[0].size();
271  Array<OneD, NekDouble> divVel(nPts, 0.0);
272 
273  // Stokes hypothesis
274  const NekDouble lambda = -2.0 / 3.0;
275 
276  // Update viscosity and thermal conductivity
277  Array<OneD, NekDouble> mu(nPts, 0.0);
278  Array<OneD, NekDouble> thermalConductivity(nPts, 0.0);
279  GetViscosityAndThermalCondFromTemp(physfield[nScalar - 1], mu,
280  thermalConductivity);
281 
282  // Velocity divergence
283  for (int j = 0; j < m_spacedim; ++j)
284  {
285  Vmath::Vadd(nPts, divVel, 1, derivativesO1[j][j], 1, divVel, 1);
286  }
287 
288  // Velocity divergence scaled by lambda * mu
289  Vmath::Smul(nPts, lambda, divVel, 1, divVel, 1);
290  Vmath::Vmul(nPts, mu, 1, divVel, 1, divVel, 1);
291 
292  // Viscous flux vector for the rho equation = 0
293  for (int i = 0; i < m_spacedim; ++i)
294  {
295  Vmath::Zero(nPts, viscousTensor[i][0], 1);
296  }
297 
298  // Viscous stress tensor (for the momentum equations)
299  for (int i = 0; i < m_spacedim; ++i)
300  {
301  for (int j = i; j < m_spacedim; ++j)
302  {
303  Vmath::Vadd(nPts, derivativesO1[i][j], 1, derivativesO1[j][i], 1,
304  viscousTensor[i][j + 1], 1);
305 
306  Vmath::Vmul(nPts, mu, 1, viscousTensor[i][j + 1], 1,
307  viscousTensor[i][j + 1], 1);
308 
309  if (i == j)
310  {
311  // Add divergence term to diagonal
312  Vmath::Vadd(nPts, viscousTensor[i][j + 1], 1, divVel, 1,
313  viscousTensor[i][j + 1], 1);
314  }
315  else
316  {
317  // Copy to make symmetric
318  Vmath::Vcopy(nPts, viscousTensor[i][j + 1], 1,
319  viscousTensor[j][i + 1], 1);
320  }
321  }
322  }
323 
324  // Terms for the energy equation
325  for (int i = 0; i < m_spacedim; ++i)
326  {
327  Vmath::Zero(nPts, viscousTensor[i][m_spacedim + 1], 1);
328  // u_j * tau_ij
329  for (int j = 0; j < m_spacedim; ++j)
330  {
331  Vmath::Vvtvp(nPts, physfield[j], 1, viscousTensor[i][j + 1], 1,
332  viscousTensor[i][m_spacedim + 1], 1,
333  viscousTensor[i][m_spacedim + 1], 1);
334  }
335  // Add k*T_i
336  Vmath::Vvtvp(nPts, thermalConductivity, 1, derivativesO1[i][m_spacedim],
337  1, viscousTensor[i][m_spacedim + 1], 1,
338  viscousTensor[i][m_spacedim + 1], 1);
339  }
340 }

References GetViscosityAndThermalCondFromTemp(), Nektar::SolverUtils::EquationSystem::m_spacedim, Vmath::Smul(), Vmath::Vadd(), Vmath::Vcopy(), Vmath::Vmul(), Vmath::Vvtvp(), and Vmath::Zero().

Referenced by InitObject_Explicit().

◆ v_GetViscousFluxVectorDeAlias()

void Nektar::NavierStokesCFE::v_GetViscousFluxVectorDeAlias ( const Array< OneD, const Array< OneD, NekDouble >> &  physfield,
TensorOfArray3D< NekDouble > &  derivativesO1,
TensorOfArray3D< NekDouble > &  viscousTensor 
)
protectedvirtual

Return the flux vector for the LDG diffusion problem.

Todo:
Complete the viscous flux vector

Reimplemented in Nektar::NavierStokesCFEAxisym.

Definition at line 346 of file NavierStokesCFE.cpp.

350 {
351  // Factor to rescale 1d points in dealiasing.
352  NekDouble OneDptscale = 2;
353  // Get number of points to dealias a cubic non-linearity
354  size_t nScalar = physfield.size();
355  size_t nPts = m_fields[0]->Get1DScaledTotPoints(OneDptscale);
356  size_t nPts_orig = physfield[0].size();
357 
358  // Auxiliary variables
359  Array<OneD, NekDouble> divVel(nPts, 0.0);
360 
361  // Stokes hypothesis
362  const NekDouble lambda = -2.0 / 3.0;
363 
364  // Update viscosity and thermal conductivity
365  Array<OneD, NekDouble> mu(nPts, 0.0);
366  Array<OneD, NekDouble> thermalConductivity(nPts, 0.0);
367  GetViscosityAndThermalCondFromTemp(physfield[nScalar - 1], mu,
368  thermalConductivity);
369 
370  // Interpolate inputs and initialise interpolated output
374  for (int i = 0; i < m_spacedim; ++i)
375  {
376  // Interpolate velocity
377  vel_interp[i] = Array<OneD, NekDouble>(nPts);
378  m_fields[0]->PhysInterp1DScaled(OneDptscale, physfield[i],
379  vel_interp[i]);
380 
381  // Interpolate derivatives
382  deriv_interp[i] = Array<OneD, Array<OneD, NekDouble>>(m_spacedim + 1);
383  for (int j = 0; j < m_spacedim + 1; ++j)
384  {
385  deriv_interp[i][j] = Array<OneD, NekDouble>(nPts);
386  m_fields[0]->PhysInterp1DScaled(OneDptscale, derivativesO1[i][j],
387  deriv_interp[i][j]);
388  }
389 
390  // Output (start from j=1 since flux is zero for rho)
391  out_interp[i] = Array<OneD, Array<OneD, NekDouble>>(m_spacedim + 2);
392  for (int j = 1; j < m_spacedim + 2; ++j)
393  {
394  out_interp[i][j] = Array<OneD, NekDouble>(nPts);
395  }
396  }
397 
398  // Velocity divergence
399  for (int j = 0; j < m_spacedim; ++j)
400  {
401  Vmath::Vadd(nPts, divVel, 1, deriv_interp[j][j], 1, divVel, 1);
402  }
403 
404  // Velocity divergence scaled by lambda * mu
405  Vmath::Smul(nPts, lambda, divVel, 1, divVel, 1);
406  Vmath::Vmul(nPts, mu, 1, divVel, 1, divVel, 1);
407 
408  // Viscous flux vector for the rho equation = 0 (no need to dealias)
409  for (int i = 0; i < m_spacedim; ++i)
410  {
411  Vmath::Zero(nPts_orig, viscousTensor[i][0], 1);
412  }
413 
414  // Viscous stress tensor (for the momentum equations)
415  for (int i = 0; i < m_spacedim; ++i)
416  {
417  for (int j = i; j < m_spacedim; ++j)
418  {
419  Vmath::Vadd(nPts, deriv_interp[i][j], 1, deriv_interp[j][i], 1,
420  out_interp[i][j + 1], 1);
421 
422  Vmath::Vmul(nPts, mu, 1, out_interp[i][j + 1], 1,
423  out_interp[i][j + 1], 1);
424 
425  if (i == j)
426  {
427  // Add divergence term to diagonal
428  Vmath::Vadd(nPts, out_interp[i][j + 1], 1, divVel, 1,
429  out_interp[i][j + 1], 1);
430  }
431  else
432  {
433  // Make symmetric
434  out_interp[j][i + 1] = out_interp[i][j + 1];
435  }
436  }
437  }
438 
439  // Terms for the energy equation
440  for (int i = 0; i < m_spacedim; ++i)
441  {
442  Vmath::Zero(nPts, out_interp[i][m_spacedim + 1], 1);
443  // u_j * tau_ij
444  for (int j = 0; j < m_spacedim; ++j)
445  {
446  Vmath::Vvtvp(nPts, vel_interp[j], 1, out_interp[i][j + 1], 1,
447  out_interp[i][m_spacedim + 1], 1,
448  out_interp[i][m_spacedim + 1], 1);
449  }
450  // Add k*T_i
451  Vmath::Vvtvp(nPts, thermalConductivity, 1, deriv_interp[i][m_spacedim],
452  1, out_interp[i][m_spacedim + 1], 1,
453  out_interp[i][m_spacedim + 1], 1);
454  }
455 
456  // Project to original space
457  for (int i = 0; i < m_spacedim; ++i)
458  {
459  for (int j = 1; j < m_spacedim + 2; ++j)
460  {
461  m_fields[0]->PhysGalerkinProjection1DScaled(
462  OneDptscale, out_interp[i][j], viscousTensor[i][j]);
463  }
464  }
465 }

References GetViscosityAndThermalCondFromTemp(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_spacedim, Vmath::Smul(), Vmath::Vadd(), Vmath::Vmul(), Vmath::Vvtvp(), and Vmath::Zero().

Referenced by InitObject_Explicit().

◆ v_InitObject()

void Nektar::NavierStokesCFE::v_InitObject ( bool  DeclareField = true)
overrideprotectedvirtual

Initialization object for CompressibleFlowSystem class.

Reimplemented from Nektar::CompressibleFlowSystem.

Reimplemented in Nektar::NavierStokesImplicitCFE, and Nektar::NavierStokesCFEAxisym.

Definition at line 60 of file NavierStokesCFE.cpp.

61 {
63 
64  // rest of initialisation is in this routine so it can also be called
65  // in NavierStokesImplicitCFE initialisation
67 }
virtual void v_InitObject(bool DeclareFields=true) override
Initialization object for CompressibleFlowSystem class.

References InitObject_Explicit(), and Nektar::CompressibleFlowSystem::v_InitObject().

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

◆ v_SupportsShockCaptType()

bool Nektar::NavierStokesCFE::v_SupportsShockCaptType ( const std::string  type) const
overrideprotectedvirtual

Implements Nektar::CompressibleFlowSystem.

Reimplemented in Nektar::NavierStokesImplicitCFE.

Definition at line 966 of file NavierStokesCFE.cpp.

967 {
968  if (type == "NonSmooth" || type == "Physical" || type == "Off")
969  {
970  return true;
971  }
972  else
973  {
974  return false;
975  }
976 }

Friends And Related Function Documentation

◆ MemoryManager< NavierStokesCFE >

friend class MemoryManager< NavierStokesCFE >
friend

Definition at line 1 of file NavierStokesCFE.h.

Member Data Documentation

◆ className

string Nektar::NavierStokesCFE::className
static
Initial value:
=
"NavierStokesCFE", NavierStokesCFE::create,
"NavierStokes equations in conservative variables.")
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
EquationSystemFactory & GetEquationSystemFactory()

Definition at line 65 of file NavierStokesCFE.h.

◆ m_C0ProjectExp

MultiRegions::ContFieldSharedPtr Nektar::NavierStokesCFE::m_C0ProjectExp
protected

Definition at line 87 of file NavierStokesCFE.h.

◆ m_Cp

NekDouble Nektar::NavierStokesCFE::m_Cp
protected

◆ m_Cv

NekDouble Nektar::NavierStokesCFE::m_Cv
protected

◆ m_eos

EquationOfStateSharedPtr Nektar::NavierStokesCFE::m_eos
protected

Equation of system for computing temperature.

Definition at line 90 of file NavierStokesCFE.h.

◆ m_is_diffIP

bool Nektar::NavierStokesCFE::m_is_diffIP {false}
protected

flag to switch between IP and LDG an enum could be added for more options

Definition at line 76 of file NavierStokesCFE.h.

Referenced by InitObject_Explicit(), v_DoDiffusion(), and Nektar::NavierStokesImplicitCFE::v_DoDiffusionCoeff().

◆ m_is_mu_variable

bool Nektar::NavierStokesCFE::m_is_mu_variable {false}
protected

flag to switch between constant viscosity and Sutherland an enum could be added for more options

Definition at line 73 of file NavierStokesCFE.h.

Referenced by GetViscosityFromTempKernel(), and InitObject_Explicit().

◆ m_is_shockCaptPhys

bool Nektar::NavierStokesCFE::m_is_shockCaptPhys {false}
protected

flag for shock capturing switch on/off an enum could be added for more options

Definition at line 79 of file NavierStokesCFE.h.

Referenced by GetViscosityAndThermalCondFromTemp(), InitObject_Explicit(), v_DoDiffusion(), Nektar::NavierStokesImplicitCFE::v_DoDiffusionCoeff(), and v_ExtraFldOutput().

◆ m_muRef

NekDouble Nektar::NavierStokesCFE::m_muRef
protected

Definition at line 93 of file NavierStokesCFE.h.

Referenced by GetViscosityFromTempKernel(), and InitObject_Explicit().

◆ m_physicalSensorType

std::string Nektar::NavierStokesCFE::m_physicalSensorType
protected

Definition at line 85 of file NavierStokesCFE.h.

◆ m_Prandtl

NekDouble Nektar::NavierStokesCFE::m_Prandtl
protected

◆ m_smoothing

std::string Nektar::NavierStokesCFE::m_smoothing
protected

Definition at line 86 of file NavierStokesCFE.h.

◆ m_thermalConductivityRef

NekDouble Nektar::NavierStokesCFE::m_thermalConductivityRef
protected

Definition at line 94 of file NavierStokesCFE.h.

Referenced by InitObject_Explicit().

◆ m_Twall

NekDouble Nektar::NavierStokesCFE::m_Twall
protected

Definition at line 92 of file NavierStokesCFE.h.

Referenced by InitObject_Explicit(), and SpecialBndTreat().

◆ m_ViscosityType

std::string Nektar::NavierStokesCFE::m_ViscosityType
protected