Nektar++
Public Member Functions | Protected Member Functions | Protected Attributes | Static Protected Attributes | List of all members
Nektar::IncNavierStokes Class Referenceabstract

This class is the base class for Navier Stokes problems. More...

#include <IncNavierStokes.h>

Inheritance diagram for Nektar::IncNavierStokes:
[legend]

Public Member Functions

virtual ~IncNavierStokes ()
 
virtual void v_InitObject (bool DeclareField=true) override
 Init object for UnsteadySystem class. More...
 
int GetNConvectiveFields (void)
 
void AddForcing (const SolverUtils::ForcingSharedPtr &pForce)
 
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
 
virtual void v_SetMovingFrameVelocities (const Array< OneD, NekDouble > &vFrameVels) override
 
virtual void v_GetMovingFrameVelocities (Array< OneD, NekDouble > &vFrameVels) override
 
virtual void v_SetMovingFrameAngles (const Array< OneD, NekDouble > &vFrameTheta) override
 
virtual void v_GetMovingFrameAngles (Array< OneD, NekDouble > &vFrameTheta) override
 
virtual void v_SetMovingFrameProjectionMat (const bnu::matrix< NekDouble > &vProjMat) override
 
virtual void v_GetMovingFrameProjectionMat (bnu::matrix< NekDouble > &vProjMat) override
 
bool DefinedForcing (const std::string &sForce)
 
void GetPivotPoint (Array< OneD, NekDouble > &vPivotPoint)
 
- 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)
 

Protected Member Functions

 IncNavierStokes (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Constructor. More...
 
EquationType GetEquationType (void)
 
void EvaluateAdvectionTerms (const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
 
void WriteModalEnergy (void)
 
void SetBoundaryConditions (NekDouble time)
 time dependent boundary conditions updating More...
 
void SetRadiationBoundaryForcing (int fieldid)
 Set Radiation forcing term. More...
 
void SetZeroNormalVelocity ()
 Set Normal Velocity Component to Zero. More...
 
void SetWomersleyBoundary (const int fldid, const int bndid)
 Set Womersley Profile if specified. More...
 
void SetUpWomersley (const int fldid, const int bndid, std::string womstr)
 Set Up Womersley details. More...
 
void SetMovingReferenceFrameBCs (const NekDouble &time)
 Set the moving reference frame boundary conditions. More...
 
void SetMRFWallBCs (const NekDouble &time)
 
void SetMRFDomainVelBCs (const NekDouble &time)
 
virtual MultiRegions::ExpListSharedPtr v_GetPressure () override
 
virtual void v_TransCoeffToPhys (void) override
 Virtual function for transformation to physical space. More...
 
virtual void v_TransPhysToCoeff (void) override
 Virtual function for transformation to coefficient space. More...
 
virtual int v_GetForceDimension ()=0
 
virtual Array< OneD, NekDoublev_GetMaxStdVelocity (const NekDouble SpeedSoundFactor) override
 
virtual bool v_PreIntegrate (int step) 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 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_RequireFwdTrans ()
 
virtual SOLVER_UTILS_EXPORT void v_SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
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_SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 
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 void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble >> &fieldcoeffs, std::vector< std::string > &variables)
 
- Protected Member Functions inherited from Nektar::SolverUtils::FluidInterface
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)
 

Protected Attributes

ExtrapolateSharedPtr m_extrapolation
 
std::ofstream m_mdlFile
 modal energy file More...
 
bool m_SmoothAdvection
 bool to identify if advection term smoothing is requested More...
 
std::vector< SolverUtils::ForcingSharedPtrm_forcing
 Forcing terms. More...
 
int m_nConvectiveFields
 Number of fields to be convected;. More...
 
Array< OneD, int > m_velocity
 int which identifies which components of m_fields contains the velocity (u,v,w); More...
 
MultiRegions::ExpListSharedPtr m_pressure
 Pointer to field holding pressure field. More...
 
NekDouble m_kinvis
 Kinematic viscosity. More...
 
int m_energysteps
 dump energy to file at steps time More...
 
EquationType m_equationType
 equation type; More...
 
Array< OneD, Array< OneD, int > > m_fieldsBCToElmtID
 Mapping from BCs to Elmt IDs. More...
 
Array< OneD, Array< OneD, int > > m_fieldsBCToTraceID
 Mapping from BCs to Elmt Edge IDs. More...
 
Array< OneD, Array< OneD, NekDouble > > m_fieldsRadiationFactor
 RHS Factor for Radiation Condition. More...
 
int m_intSteps
 Number of time integration steps AND Order of extrapolation for pressure boundary conditions. More...
 
Array< OneD, NekDoublem_pivotPoint
 
std::map< int, std::map< int, WomersleyParamsSharedPtr > > m_womersleyParams
 Womersley parameters if required. More...
 
- 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...
 

Static Protected Attributes

static std::string eqTypeLookupIds []
 
- Static Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
static std::string equationSystemTypeLookupIds []
 

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...
 
- Static Public Attributes inherited from Nektar::SolverUtils::UnsteadySystem
static std::string cmdSetStartTime
 
static std::string cmdSetStartChkNum
 
- Protected Types inherited from Nektar::SolverUtils::EquationSystem
enum  HomogeneousType { eHomogeneous1D , eHomogeneous2D , eHomogeneous3D , eNotHomogeneous }
 Parameter for homogeneous expansions. More...
 

Detailed Description

This class is the base class for Navier Stokes problems.

Definition at line 133 of file IncNavierStokes.h.

Constructor & Destructor Documentation

◆ ~IncNavierStokes()

Nektar::IncNavierStokes::~IncNavierStokes ( void  )
virtual

Destructor

Definition at line 392 of file IncNavierStokes.cpp.

393 {
394 }

◆ IncNavierStokes()

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

Constructor.

Constructor. Creates ...

Parameters

Definition at line 83 of file IncNavierStokes.cpp.

86  : UnsteadySystem(pSession, pGraph), AdvectionSystem(pSession, pGraph),
87  m_SmoothAdvection(false)
88 {
89 }
bool m_SmoothAdvection
bool to identify if advection term smoothing is requested
SOLVER_UTILS_EXPORT AdvectionSystem(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

◆ AddForcing()

void Nektar::IncNavierStokes::AddForcing ( const SolverUtils::ForcingSharedPtr pForce)

Add an additional forcing term programmatically.

Definition at line 1131 of file IncNavierStokes.cpp.

1132 {
1133  m_forcing.push_back(pForce);
1134 }
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Forcing terms.

References m_forcing.

Referenced by Nektar::VortexWaveInteraction::ExecuteRoll().

◆ DefinedForcing()

bool Nektar::IncNavierStokes::DefinedForcing ( const std::string &  sForce)

Function to check the type of forcing

Definition at line 1308 of file IncNavierStokes.cpp.

1309 {
1310  vector<std::string> vForceList;
1311  bool hasForce{false};
1312 
1313  if (!m_session->DefinesElement("Nektar/Forcing"))
1314  {
1315  return hasForce;
1316  }
1317 
1318  TiXmlElement *vForcing = m_session->GetElement("Nektar/Forcing");
1319  if (vForcing)
1320  {
1321  TiXmlElement *vForce = vForcing->FirstChildElement("FORCE");
1322  while (vForce)
1323  {
1324  string vType = vForce->Attribute("TYPE");
1325 
1326  vForceList.push_back(vType);
1327  vForce = vForce->NextSiblingElement("FORCE");
1328  }
1329  }
1330 
1331  for (auto &f : vForceList)
1332  {
1333  if (boost::iequals(f, sForce))
1334  {
1335  hasForce = true;
1336  }
1337  }
1338 
1339  return hasForce;
1340 }
LibUtilities::SessionReaderSharedPtr m_session
The session reader.

References Nektar::SolverUtils::EquationSystem::m_session.

Referenced by GetPivotPoint(), and v_InitObject().

◆ EvaluateAdvectionTerms()

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

Evaluation -N(V) for all fields except pressure using m_velocity

Definition at line 399 of file IncNavierStokes.cpp.

402 {
403  size_t VelDim = m_velocity.size();
404  Array<OneD, Array<OneD, NekDouble>> velocity(VelDim);
405 
406  size_t npoints = m_fields[0]->GetNpoints();
407  for (size_t i = 0; i < VelDim; ++i)
408  {
409  velocity[i] = Array<OneD, NekDouble>(npoints);
410  Vmath::Vcopy(npoints, inarray[m_velocity[i]], 1, velocity[i], 1);
411  }
412  for (auto &x : m_forcing)
413  {
414  x->PreApply(m_fields, velocity, velocity, time);
415  }
416 
417  m_advObject->Advect(m_nConvectiveFields, m_fields, velocity, inarray,
418  outarray, time);
419 }
Array< OneD, int > m_velocity
int which identifies which components of m_fields contains the velocity (u,v,w);
int m_nConvectiveFields
Number of fields to be convected;.
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255

References Nektar::SolverUtils::AdvectionSystem::m_advObject, Nektar::SolverUtils::EquationSystem::m_fields, m_forcing, m_nConvectiveFields, m_velocity, and Vmath::Vcopy().

Referenced by Nektar::CoupledLinearNS::EvaluateAdvection(), Nektar::CoupledLinearNS::EvaluateNewtonRHS(), Nektar::VelocityCorrectionScheme::v_EvaluateAdvection_SetPressureBCs(), and Nektar::VCSMapping::v_EvaluateAdvection_SetPressureBCs().

◆ GetEquationType()

EquationType Nektar::IncNavierStokes::GetEquationType ( void  )
inlineprotected

Definition at line 232 of file IncNavierStokes.h.

233  {
234  return m_equationType;
235  }
EquationType m_equationType
equation type;

References m_equationType.

◆ GetNConvectiveFields()

int Nektar::IncNavierStokes::GetNConvectiveFields ( void  )
inline

Definition at line 142 of file IncNavierStokes.h.

143  {
144  return m_nConvectiveFields;
145  }

References m_nConvectiveFields.

◆ GetPivotPoint()

void Nektar::IncNavierStokes::GetPivotPoint ( Array< OneD, NekDouble > &  vPivotPoint)

Get the pivot point for moving reference

Definition at line 1345 of file IncNavierStokes.cpp.

1346 {
1347  std::string sMRFForcingType = "MovingReferenceFrame";
1348 
1349  // only if we use moving reference frame formulation
1350  if (!DefinedForcing(sMRFForcingType))
1351  {
1352  return;
1353  }
1354 
1355  TiXmlElement *vForcing = m_session->GetElement("Nektar/Forcing");
1356  if (vForcing)
1357  {
1358  TiXmlElement *vForce = vForcing->FirstChildElement("FORCE");
1359  while (vForce)
1360  {
1361  string vType = vForce->Attribute("TYPE");
1362 
1363  // if it is moving reference frame
1364  if (boost::iequals(vType, sMRFForcingType))
1365  {
1366  TiXmlElement *pivotElmt =
1367  vForce->FirstChildElement("PivotPoint");
1368 
1369  // if not defined, zero would be the default
1370  if (pivotElmt)
1371  {
1372  std::stringstream pivotPointStream;
1373  std::string pivotPointStr = pivotElmt->GetText();
1374  pivotPointStream.str(pivotPointStr);
1375 
1376  for (int j = 0; j < m_spacedim; ++j)
1377  {
1378  pivotPointStream >> pivotPointStr;
1379  if (!pivotPointStr.empty())
1380  {
1381  LibUtilities::Equation equ(
1382  m_session->GetInterpreter(), pivotPointStr);
1383  vPivotPoint[j] = equ.Evaluate();
1384  }
1385  }
1386  }
1387 
1388  break;
1389  }
1390  }
1391  }
1392 }
bool DefinedForcing(const std::string &sForce)
int m_spacedim
Spatial dimension (>= expansion dim).

References DefinedForcing(), Nektar::LibUtilities::Equation::Evaluate(), Nektar::SolverUtils::EquationSystem::m_session, and Nektar::SolverUtils::EquationSystem::m_spacedim.

Referenced by v_InitObject().

◆ SetBoundaryConditions()

void Nektar::IncNavierStokes::SetBoundaryConditions ( NekDouble  time)
protected

time dependent boundary conditions updating

Time dependent boundary conditions updating

Definition at line 424 of file IncNavierStokes.cpp.

425 {
426  size_t i, n;
427  std::string varName;
428  size_t nvariables = m_fields.size();
429 
430  for (i = 0; i < nvariables; ++i)
431  {
432  for (n = 0; n < m_fields[i]->GetBndConditions().size(); ++n)
433  {
434  if (m_fields[i]->GetBndConditions()[n]->IsTimeDependent())
435  {
436  varName = m_session->GetVariable(i);
437  m_fields[i]->EvaluateBoundaryConditions(time, varName);
438  }
439  else if (boost::istarts_with(
440  m_fields[i]->GetBndConditions()[n]->GetUserDefined(),
441  "Womersley"))
442  {
443  SetWomersleyBoundary(i, n);
444  }
445  }
446 
447  // Set Radiation conditions if required
449  }
450 
451  // Enforcing the boundary conditions (Inlet and wall) for the
452  // Moving reference frame
455 }
void SetWomersleyBoundary(const int fldid, const int bndid)
Set Womersley Profile if specified.
void SetZeroNormalVelocity()
Set Normal Velocity Component to Zero.
void SetRadiationBoundaryForcing(int fieldid)
Set Radiation forcing term.
void SetMovingReferenceFrameBCs(const NekDouble &time)
Set the moving reference frame boundary conditions.

References Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_session, SetMovingReferenceFrameBCs(), SetRadiationBoundaryForcing(), SetWomersleyBoundary(), and SetZeroNormalVelocity().

Referenced by Nektar::VelocityCorrectionScheme::v_DoInitialise(), and v_PreIntegrate().

◆ SetMovingReferenceFrameBCs()

void Nektar::IncNavierStokes::SetMovingReferenceFrameBCs ( const NekDouble time)
protected

Set the moving reference frame boundary conditions.

Set boundary conditions for moving frame of reference

Definition at line 906 of file IncNavierStokes.cpp.

907 {
908  SetMRFWallBCs(time);
909  SetMRFDomainVelBCs(time);
910 }
void SetMRFDomainVelBCs(const NekDouble &time)
void SetMRFWallBCs(const NekDouble &time)

References SetMRFDomainVelBCs(), and SetMRFWallBCs().

Referenced by SetBoundaryConditions().

◆ SetMRFDomainVelBCs()

void Nektar::IncNavierStokes::SetMRFDomainVelBCs ( const NekDouble time)
protected

Set inlet boundary conditions for moving frame of reference

Definition at line 1032 of file IncNavierStokes.cpp.

1033 {
1034  // The inlet conditions for the velocities given in the session xml file
1035  // however, those are in the inertial reference coordinate XYZ and
1036  // needed to be converted to the moving reference coordinate xyz
1037 
1038  // define arrays to calculate the prescribed velocities and modifed ones
1039  Array<OneD, Array<OneD, NekDouble>> definedVels(m_velocity.size());
1040  Array<OneD, Array<OneD, NekDouble>> velocities(m_velocity.size());
1041  Array<OneD, Array<OneD, NekDouble>> unitArray(m_velocity.size());
1042  Array<OneD, Array<OneD, NekDouble>> coords(3);
1043 
1044  Array<OneD, Array<OneD, const SpatialDomains::BoundaryConditionShPtr>>
1045  BndConds(m_spacedim);
1046  Array<OneD, Array<OneD, MultiRegions::ExpListSharedPtr>> BndExp(m_spacedim);
1047 
1048  for (int i = 0; i < m_spacedim; ++i)
1049  {
1050  BndConds[i] = m_fields[m_velocity[i]]->GetBndConditions();
1051  BndExp[i] = m_fields[m_velocity[i]]->GetBndCondExpansions();
1052  }
1053 
1054  int npoints;
1055  Array<OneD, NekDouble> Bphys, Bcoeffs;
1056 
1057  // loop over the boundary regions
1058  for (size_t n = 0; n < BndExp[0].size(); ++n)
1059  {
1060  if (BndConds[0][n]->GetBoundaryConditionType() ==
1062  (boost::iequals(BndConds[0][n]->GetUserDefined(),
1063  "MovingFrameDomainVel")))
1064  {
1065  npoints = BndExp[0][n]->GetNpoints();
1066  for (size_t k = 0; k < m_velocity.size(); ++k)
1067  {
1068  definedVels[k] = Array<OneD, NekDouble>(npoints, 0.0);
1069  velocities[k] = Array<OneD, NekDouble>(npoints, 0.0);
1070  unitArray[k] = Array<OneD, NekDouble>(npoints, 1.0);
1071  }
1072  Array<OneD, NekDouble> tmp(npoints, 0.0);
1073  for (size_t k = 0; k < 3; ++k)
1074  {
1075  coords[k] = Array<OneD, NekDouble>(npoints, 0.0);
1076  }
1077  BndExp[0][n]->GetCoords(coords[0], coords[1], coords[2]);
1078 
1079  // loop over the velocity fields and compute the boundary
1080  // condition
1081  for (size_t k = 0; k < m_velocity.size(); ++k)
1082  {
1083  LibUtilities::Equation condition =
1084  std::static_pointer_cast<
1085  SpatialDomains::DirichletBoundaryCondition>(
1086  BndConds[k][n])
1087  ->m_dirichletCondition;
1088  // Evaluate
1089  condition.Evaluate(coords[0], coords[1], coords[2], time,
1090  definedVels[k]);
1091  }
1092 
1093  // We have all velocity components
1094  // transform them to the moving refernce frame
1095  for (int l = 0; l < m_spacedim; ++l)
1096  {
1097  Array<OneD, NekDouble> tmp0(npoints, 0.0);
1098  Array<OneD, NekDouble> tmp1(npoints, 0.0);
1099  Array<OneD, NekDouble> tmp2(npoints, 0.0);
1100  for (int m = 0; m < m_spacedim; ++m)
1101  {
1102  Vmath::Svtvp(npoints, m_movingFrameProjMat(l, m),
1103  tmp0 = definedVels[m], 1, tmp1 = velocities[l],
1104  1, tmp2 = velocities[l], 1);
1105  }
1106  }
1107 
1108  // update the boundary values
1109  for (int k = 0; k < m_spacedim; ++k)
1110  {
1111  Vmath::Vmul(npoints, unitArray[k], 1, velocities[k], 1,
1112  BndExp[k][n]->UpdatePhys(), 1);
1113  }
1114  // update the coefficients
1115  for (int k = 0; k < m_spacedim; ++k)
1116  {
1117  if (m_fields[k]->GetExpType() == MultiRegions::e3DH1D)
1118  {
1119  BndExp[k][n]->SetWaveSpace(false);
1120  }
1121  BndExp[k][n]->FwdTransBndConstrained(
1122  BndExp[k][n]->GetPhys(), BndExp[k][n]->UpdateCoeffs());
1123  }
1124  }
1125  }
1126 }
boost::numeric::ublas::matrix< NekDouble > m_movingFrameProjMat
Projection matrix for transformation between inertial and moving.
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 Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
svtvp (scalar times vector plus vector): z = alpha*x + y
Definition: Vmath.cpp:622

References Nektar::MultiRegions::e3DH1D, Nektar::SpatialDomains::eDirichlet, Nektar::LibUtilities::Equation::Evaluate(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_movingFrameProjMat, Nektar::SolverUtils::EquationSystem::m_spacedim, m_velocity, Vmath::Svtvp(), and Vmath::Vmul().

Referenced by SetMovingReferenceFrameBCs().

◆ SetMRFWallBCs()

void Nektar::IncNavierStokes::SetMRFWallBCs ( const NekDouble time)
protected

Set Wall boundary conditions for moving frame of reference

Definition at line 915 of file IncNavierStokes.cpp.

916 {
917  boost::ignore_unused(time);
918 
919  // for the wall we need to calculate:
920  // [V_wall]_xyz = [V_frame]_xyz + [Omega X r]_xyz
921  // Note all vectors must be in moving frame coordinates xyz
922  // not in inertial frame XYZ
923 
924  Array<OneD, Array<OneD, const SpatialDomains::BoundaryConditionShPtr>>
925  BndConds(m_spacedim);
926  Array<OneD, Array<OneD, MultiRegions::ExpListSharedPtr>> BndExp(m_spacedim);
927 
928  for (int i = 0; i < m_spacedim; ++i)
929  {
930  BndConds[i] = m_fields[m_velocity[i]]->GetBndConditions();
931  BndExp[i] = m_fields[m_velocity[i]]->GetBndCondExpansions();
932  }
933 
934  int npoints;
935  Array<OneD, NekDouble> Bphys, Bcoeffs;
936 
937  // loop over the boundary regions
938  for (size_t n = 0; n < BndExp[0].size(); ++n)
939  {
940  if (BndConds[0][n]->GetBoundaryConditionType() ==
942  (boost::iequals(BndConds[0][n]->GetUserDefined(),
943  "MovingFrameWall")))
944  {
945  npoints = BndExp[0][n]->GetNpoints();
946  //
947  // define arrays to calculate the prescribed velocities and modifed
948  // ones
949  Array<OneD, Array<OneD, NekDouble>> velocities(m_velocity.size());
950  Array<OneD, Array<OneD, NekDouble>> unitArray(m_velocity.size());
951  Array<OneD, Array<OneD, NekDouble>> coords(3);
952  for (size_t k = 0; k < m_velocity.size(); ++k)
953  {
954  velocities[k] = Array<OneD, NekDouble>(npoints, 0.0);
955  unitArray[k] = Array<OneD, NekDouble>(npoints, 1.0);
956  }
957  Array<OneD, NekDouble> tmp(npoints, 0.0);
958  for (size_t k = 0; k < 3; ++k)
959  {
960  coords[k] = Array<OneD, NekDouble>(npoints, 0.0);
961  }
962  BndExp[0][n]->GetCoords(coords[0], coords[1], coords[2]);
963 
964  // move the centre to the location of pivot
965  for (int i = 0; i < m_spacedim; ++i)
966  {
967  Vmath::Sadd(npoints, -m_pivotPoint[i], coords[i], 1, coords[i],
968  1);
969  }
970  /////////////////////////////////////
971  // Note that both Omega and the absolute velocities have been
972  // expressed in the moving reference frame unit vectors
973  // therefore the result is in moving ref frame and no furhter
974  // transformaton is required
975  //
976  // compute Omega X r = vx ex + vy ey + vz ez
977  // Note OmegaX : movingFrameVelsxyz[3]
978  // Note OmegaY : movingFrameVelsxyz[4]
979  // Note OmegaZ : movingFrameVelsxyz[5]
980  //
981  // vx = OmegaY*z-OmegaZ*y
982  Vmath::Smul(npoints, -1 * m_movingFrameVelsxyz[5], coords[1], 1,
983  velocities[0], 1);
984  // vy = OmegaZ*x-OmegaX*z
985  Vmath::Smul(npoints, m_movingFrameVelsxyz[5], coords[0], 1,
986  velocities[1], 1);
987  if (m_spacedim == 3)
988  {
989  // add the OmegaY*z to vx
990  Vmath::Svtvp(npoints, m_movingFrameVelsxyz[4], coords[2], 1,
991  velocities[0], 1, velocities[0], 1);
992  // add the -OmegaX*z to vy
993  Vmath::Svtvp(npoints, -1 * m_movingFrameVelsxyz[3], coords[2],
994  1, velocities[1], 1, velocities[1], 1);
995 
996  // vz = OmegaX*y-OmegaY*x
997  Vmath::Svtsvtp(npoints, m_movingFrameVelsxyz[3], coords[1], 1,
998  -1.0 * m_movingFrameVelsxyz[4], coords[0], 1,
999  velocities[2], 1);
1000  }
1001 
1002  // add the translation velocity
1003  for (int k = 0; k < m_spacedim; ++k)
1004  {
1005  Vmath::Sadd(npoints, m_movingFrameVelsxyz[k], velocities[k], 1,
1006  velocities[k], 1);
1007  }
1008 
1009  // update the boundary values
1010  for (int k = 0; k < m_spacedim; ++k)
1011  {
1012  Vmath::Vmul(npoints, unitArray[k], 1, velocities[k], 1,
1013  BndExp[k][n]->UpdatePhys(), 1);
1014  }
1015  // update the coefficients
1016  for (int k = 0; k < m_spacedim; ++k)
1017  {
1018  if (m_fields[k]->GetExpType() == MultiRegions::e3DH1D)
1019  {
1020  BndExp[k][n]->SetWaveSpace(false);
1021  }
1022  BndExp[k][n]->FwdTransBndConstrained(
1023  BndExp[k][n]->GetPhys(), BndExp[k][n]->UpdateCoeffs());
1024  }
1025  }
1026  }
1027 }
Array< OneD, NekDouble > m_pivotPoint
Array< OneD, NekDouble > m_movingFrameVelsxyz
Moving frame of reference velocities.
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 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
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add scalar y = alpha + x.
Definition: Vmath.cpp:384

References Nektar::MultiRegions::e3DH1D, Nektar::SpatialDomains::eDirichlet, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_movingFrameVelsxyz, m_pivotPoint, Nektar::SolverUtils::EquationSystem::m_spacedim, m_velocity, Vmath::Sadd(), Vmath::Smul(), Vmath::Svtsvtp(), Vmath::Svtvp(), and Vmath::Vmul().

Referenced by SetMovingReferenceFrameBCs().

◆ SetRadiationBoundaryForcing()

void Nektar::IncNavierStokes::SetRadiationBoundaryForcing ( int  fieldid)
protected

Set Radiation forcing term.

Probably should be pushed back into ContField?

Definition at line 460 of file IncNavierStokes.cpp.

461 {
462  size_t i, n;
463 
464  Array<OneD, const SpatialDomains::BoundaryConditionShPtr> BndConds;
465  Array<OneD, MultiRegions::ExpListSharedPtr> BndExp;
466 
467  BndConds = m_fields[fieldid]->GetBndConditions();
468  BndExp = m_fields[fieldid]->GetBndCondExpansions();
469 
472 
473  size_t cnt;
474  size_t elmtid, nq, offset, boundary;
475  Array<OneD, NekDouble> Bvals, U;
476  size_t cnt1 = 0;
477 
478  for (cnt = n = 0; n < BndConds.size(); ++n)
479  {
480  std::string type = BndConds[n]->GetUserDefined();
481 
482  if ((BndConds[n]->GetBoundaryConditionType() ==
484  (boost::iequals(type, "Radiation")))
485  {
486  size_t nExp = BndExp[n]->GetExpSize();
487  for (i = 0; i < nExp; ++i, cnt++)
488  {
489  elmtid = m_fieldsBCToElmtID[m_velocity[fieldid]][cnt];
490  elmt = m_fields[fieldid]->GetExp(elmtid);
491  offset = m_fields[fieldid]->GetPhys_Offset(elmtid);
492 
493  U = m_fields[fieldid]->UpdatePhys() + offset;
494  Bc = BndExp[n]->GetExp(i);
495 
496  boundary = m_fieldsBCToTraceID[fieldid][cnt];
497 
498  // Get edge values and put into ubc
499  nq = Bc->GetTotPoints();
500  Array<OneD, NekDouble> ubc(nq);
501  elmt->GetTracePhysVals(boundary, Bc, U, ubc);
502 
503  Vmath::Vmul(nq,
505  [fieldid][cnt1 + BndExp[n]->GetPhys_Offset(i)],
506  1, &ubc[0], 1, &ubc[0], 1);
507 
508  Bvals =
509  BndExp[n]->UpdateCoeffs() + BndExp[n]->GetCoeff_Offset(i);
510 
511  Bc->IProductWRTBase(ubc, Bvals);
512  }
513  cnt1 += BndExp[n]->GetTotPoints();
514  }
515  else
516  {
517  cnt += BndExp[n]->GetExpSize();
518  }
519  }
520 }
Array< OneD, Array< OneD, int > > m_fieldsBCToTraceID
Mapping from BCs to Elmt Edge IDs.
Array< OneD, Array< OneD, NekDouble > > m_fieldsRadiationFactor
RHS Factor for Radiation Condition.
Array< OneD, Array< OneD, int > > m_fieldsBCToElmtID
Mapping from BCs to Elmt IDs.
SOLVER_UTILS_EXPORT int GetPhys_Offset(int n)
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
std::shared_ptr< StdExpansion > StdExpansionSharedPtr

References Nektar::SpatialDomains::eRobin, Nektar::SolverUtils::EquationSystem::GetPhys_Offset(), Nektar::SolverUtils::EquationSystem::m_fields, m_fieldsBCToElmtID, m_fieldsBCToTraceID, m_fieldsRadiationFactor, m_velocity, and Vmath::Vmul().

Referenced by SetBoundaryConditions().

◆ SetUpWomersley()

void Nektar::IncNavierStokes::SetUpWomersley ( const int  fldid,
const int  bndid,
std::string  womstr 
)
protected

Set Up Womersley details.

Error value returned by TinyXML.

Definition at line 679 of file IncNavierStokes.cpp.

681 {
682  std::string::size_type indxBeg = womStr.find_first_of(':') + 1;
683  string filename = womStr.substr(indxBeg, string::npos);
684 
685  TiXmlDocument doc(filename);
686 
687  bool loadOkay = doc.LoadFile();
688  ASSERTL0(loadOkay,
689  (std::string("Failed to load file: ") + filename).c_str());
690 
691  TiXmlHandle docHandle(&doc);
692 
693  int err; /// Error value returned by TinyXML.
694 
695  TiXmlElement *nektar = doc.FirstChildElement("NEKTAR");
696  ASSERTL0(nektar, "Unable to find NEKTAR tag in file.");
697 
698  TiXmlElement *wombc = nektar->FirstChildElement("WOMERSLEYBC");
699  ASSERTL0(wombc, "Unable to find WOMERSLEYBC tag in file.");
700 
701  // read womersley parameters
702  TiXmlElement *womparam = wombc->FirstChildElement("WOMPARAMS");
703  ASSERTL0(womparam, "Unable to find WOMPARAMS tag in file.");
704 
705  // Input coefficients
706  TiXmlElement *params = womparam->FirstChildElement("W");
707  map<std::string, std::string> Wparams;
708 
709  // read parameter list
710  while (params)
711  {
712 
713  std::string propstr;
714  propstr = params->Attribute("PROPERTY");
715 
716  ASSERTL0(!propstr.empty(),
717  "Failed to read PROPERTY value Womersley BC Parameter");
718 
719  std::string valstr;
720  valstr = params->Attribute("VALUE");
721 
722  ASSERTL0(!valstr.empty(),
723  "Failed to read VALUE value Womersley BC Parameter");
724 
725  std::transform(propstr.begin(), propstr.end(), propstr.begin(),
726  ::toupper);
727  Wparams[propstr] = valstr;
728 
729  params = params->NextSiblingElement("W");
730  }
731  bool parseGood;
732 
733  // Read parameters
734 
735  ASSERTL0(
736  Wparams.count("RADIUS") == 1,
737  "Failed to find Radius parameter in Womersley boundary conditions");
738  std::vector<NekDouble> rad;
739  ParseUtils::GenerateVector(Wparams["RADIUS"], rad);
740  m_womersleyParams[fldid][bndid]->m_radius = rad[0];
741 
742  ASSERTL0(
743  Wparams.count("PERIOD") == 1,
744  "Failed to find period parameter in Womersley boundary conditions");
745  std::vector<NekDouble> period;
746  parseGood = ParseUtils::GenerateVector(Wparams["PERIOD"], period);
747  m_womersleyParams[fldid][bndid]->m_period = period[0];
748 
749  ASSERTL0(
750  Wparams.count("AXISNORMAL") == 1,
751  "Failed to find axisnormal parameter in Womersley boundary conditions");
752  std::vector<NekDouble> anorm;
753  parseGood = ParseUtils::GenerateVector(Wparams["AXISNORMAL"], anorm);
754  m_womersleyParams[fldid][bndid]->m_axisnormal[0] = anorm[0];
755  m_womersleyParams[fldid][bndid]->m_axisnormal[1] = anorm[1];
756  m_womersleyParams[fldid][bndid]->m_axisnormal[2] = anorm[2];
757 
758  ASSERTL0(
759  Wparams.count("AXISPOINT") == 1,
760  "Failed to find axispoint parameter in Womersley boundary conditions");
761  std::vector<NekDouble> apt;
762  parseGood = ParseUtils::GenerateVector(Wparams["AXISPOINT"], apt);
763  m_womersleyParams[fldid][bndid]->m_axispoint[0] = apt[0];
764  m_womersleyParams[fldid][bndid]->m_axispoint[1] = apt[1];
765  m_womersleyParams[fldid][bndid]->m_axispoint[2] = apt[2];
766 
767  // Read Temporal Fourier Coefficients.
768 
769  // Find the FourierCoeff tag
770  TiXmlElement *coeff = wombc->FirstChildElement("FOURIERCOEFFS");
771 
772  // Input coefficients
773  TiXmlElement *fval = coeff->FirstChildElement("F");
774 
775  int indx;
776  int nextFourierCoeff = -1;
777 
778  while (fval)
779  {
780  nextFourierCoeff++;
781 
782  TiXmlAttribute *fvalAttr = fval->FirstAttribute();
783  std::string attrName(fvalAttr->Name());
784 
785  ASSERTL0(attrName == "ID",
786  (std::string("Unknown attribute name: ") + attrName).c_str());
787 
788  err = fvalAttr->QueryIntValue(&indx);
789  ASSERTL0(err == TIXML_SUCCESS, "Unable to read attribute ID.");
790 
791  std::string coeffStr = fval->FirstChild()->ToText()->ValueStr();
792  vector<NekDouble> coeffvals;
793 
794  parseGood = ParseUtils::GenerateVector(coeffStr, coeffvals);
795  ASSERTL0(
796  parseGood,
797  (std::string("Problem reading value of fourier coefficient, ID=") +
798  boost::lexical_cast<string>(indx))
799  .c_str());
800  ASSERTL1(
801  coeffvals.size() == 2,
802  (std::string(
803  "Have not read two entries of Fourier coefficicent from ID=" +
804  boost::lexical_cast<string>(indx))
805  .c_str()));
806 
807  m_womersleyParams[fldid][bndid]->m_wom_vel.push_back(
808  NekComplexDouble(coeffvals[0], coeffvals[1]));
809 
810  fval = fval->NextSiblingElement("F");
811  }
812 
813  // starting point of precalculation
814  size_t i, j, k;
815  // M fourier coefficients
816  size_t M_coeffs = m_womersleyParams[fldid][bndid]->m_wom_vel.size();
817  NekDouble R = m_womersleyParams[fldid][bndid]->m_radius;
818  NekDouble T = m_womersleyParams[fldid][bndid]->m_period;
819  Array<OneD, NekDouble> x0 = m_womersleyParams[fldid][bndid]->m_axispoint;
820 
821  NekComplexDouble rqR;
822  // Womersley Number
823  NekComplexDouble omega_c(2.0 * M_PI / T, 0.0);
824  NekComplexDouble alpha_c(R * sqrt(omega_c.real() / m_kinvis), 0.0);
825  NekComplexDouble z1(1.0, 0.0);
826  NekComplexDouble i_pow_3q2(-1.0 / sqrt(2.0), 1.0 / sqrt(2.0));
827 
829  BndCondExp = m_fields[fldid]->GetBndCondExpansions()[bndid];
830 
832  size_t cnt = 0;
833  size_t nfq;
834  Array<OneD, NekDouble> Bvals;
835 
836  size_t exp_npts = BndCondExp->GetExpSize();
837  Array<OneD, NekDouble> wbc(exp_npts, 0.0);
838 
839  // allocate time indepedent variables
840  m_womersleyParams[fldid][bndid]->m_poiseuille =
841  Array<OneD, Array<OneD, NekDouble>>(exp_npts);
842  m_womersleyParams[fldid][bndid]->m_zvel =
843  Array<OneD, Array<OneD, Array<OneD, NekComplexDouble>>>(exp_npts);
844  // could use M_coeffs - 1 but need to avoid complicating things
845  Array<OneD, NekComplexDouble> zJ0(M_coeffs);
846  Array<OneD, NekComplexDouble> lamda_n(M_coeffs);
847  Array<OneD, NekComplexDouble> k_c(M_coeffs);
848  NekComplexDouble zJ0r;
849 
850  for (k = 1; k < M_coeffs; k++)
851  {
852  k_c[k] = NekComplexDouble((NekDouble)k, 0.0);
853  lamda_n[k] = i_pow_3q2 * alpha_c * sqrt(k_c[k]);
854  zJ0[k] = Polylib::ImagBesselComp(0, lamda_n[k]);
855  }
856 
857  // Loop over each element in an expansion
858  for (i = 0; i < exp_npts; ++i, cnt++)
859  {
860  // Get Boundary and trace expansion
861  bc = BndCondExp->GetExp(i);
862  nfq = bc->GetTotPoints();
863 
864  Array<OneD, NekDouble> x(nfq, 0.0);
865  Array<OneD, NekDouble> y(nfq, 0.0);
866  Array<OneD, NekDouble> z(nfq, 0.0);
867  bc->GetCoords(x, y, z);
868 
869  m_womersleyParams[fldid][bndid]->m_poiseuille[i] =
870  Array<OneD, NekDouble>(nfq);
871  m_womersleyParams[fldid][bndid]->m_zvel[i] =
872  Array<OneD, Array<OneD, NekComplexDouble>>(nfq);
873 
874  // Compute coefficients
875  for (j = 0; j < nfq; j++)
876  {
877  rqR = NekComplexDouble(sqrt((x[j] - x0[0]) * (x[j] - x0[0]) +
878  (y[j] - x0[1]) * (y[j] - x0[1]) +
879  (z[j] - x0[2]) * (z[j] - x0[2])) /
880  R,
881  0.0);
882 
883  // Compute Poiseulle Flow
884  m_womersleyParams[fldid][bndid]->m_poiseuille[i][j] =
885  m_womersleyParams[fldid][bndid]->m_wom_vel[0].real() *
886  (1. - rqR.real() * rqR.real());
887 
888  m_womersleyParams[fldid][bndid]->m_zvel[i][j] =
889  Array<OneD, NekComplexDouble>(M_coeffs);
890 
891  // compute the velocity information
892  for (k = 1; k < M_coeffs; k++)
893  {
894  zJ0r = Polylib::ImagBesselComp(0, rqR * lamda_n[k]);
895  m_womersleyParams[fldid][bndid]->m_zvel[i][j][k] =
896  m_womersleyParams[fldid][bndid]->m_wom_vel[k] *
897  (z1 - (zJ0r / zJ0[k]));
898  }
899  }
900  }
901 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
std::map< int, std::map< int, WomersleyParamsSharedPtr > > m_womersleyParams
Womersley parameters if required.
NekDouble m_kinvis
Kinematic viscosity.
static bool GenerateVector(const std::string &str, std::vector< T > &out)
Takes a comma-separated string and converts it to entries in a vector.
Definition: ParseUtils.cpp:131
static NekDouble rad(NekDouble x, NekDouble y)
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::complex< double > NekComplexDouble
double NekDouble
std::complex< Nektar::NekDouble > ImagBesselComp(int n, std::complex< Nektar::NekDouble > y)
Calcualte the bessel function of the first kind with complex double input y. Taken from Numerical Rec...
Definition: Polylib.cpp:1837
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References ASSERTL0, ASSERTL1, Nektar::ParseUtils::GenerateVector(), Polylib::ImagBesselComp(), Nektar::SolverUtils::EquationSystem::m_fields, m_kinvis, m_womersleyParams, Nektar::LibUtilities::rad(), and tinysimd::sqrt().

Referenced by v_InitObject().

◆ SetWomersleyBoundary()

void Nektar::IncNavierStokes::SetWomersleyBoundary ( const int  fldid,
const int  bndid 
)
protected

Set Womersley Profile if specified.

Womersley boundary condition defintion

Definition at line 610 of file IncNavierStokes.cpp.

611 {
612  ASSERTL1(m_womersleyParams.count(bndid) == 1,
613  "Womersley parameters for this boundary have not been set up");
614 
615  WomersleyParamsSharedPtr WomParam = m_womersleyParams[fldid][bndid];
616  NekComplexDouble zvel;
617  size_t i, j, k;
618 
619  size_t M_coeffs = WomParam->m_wom_vel.size();
620 
621  NekDouble T = WomParam->m_period;
622  NekDouble axis_normal = WomParam->m_axisnormal[fldid];
623 
624  // Womersley Number
625  NekComplexDouble omega_c(2.0 * M_PI / T, 0.0);
626  NekComplexDouble k_c(0.0, 0.0);
627  NekComplexDouble m_time_c(m_time, 0.0);
628  NekComplexDouble zi(0.0, 1.0);
629  NekComplexDouble i_pow_3q2(-1.0 / sqrt(2.0), 1.0 / sqrt(2.0));
630 
632  BndCondExp = m_fields[fldid]->GetBndCondExpansions()[bndid];
633 
635  size_t cnt = 0;
636  size_t nfq;
637  Array<OneD, NekDouble> Bvals;
638  size_t exp_npts = BndCondExp->GetExpSize();
639  Array<OneD, NekDouble> wbc(exp_npts, 0.0);
640 
641  Array<OneD, NekComplexDouble> zt(M_coeffs);
642 
643  // preallocate the exponent
644  for (k = 1; k < M_coeffs; k++)
645  {
646  k_c = NekComplexDouble((NekDouble)k, 0.0);
647  zt[k] = std::exp(zi * omega_c * k_c * m_time_c);
648  }
649 
650  // Loop over each element in an expansion
651  for (i = 0; i < exp_npts; ++i, cnt++)
652  {
653  // Get Boundary and trace expansion
654  bc = BndCondExp->GetExp(i);
655  nfq = bc->GetTotPoints();
656  Array<OneD, NekDouble> wbc(nfq, 0.0);
657 
658  // Compute womersley solution
659  for (j = 0; j < nfq; j++)
660  {
661  wbc[j] = WomParam->m_poiseuille[i][j];
662  for (k = 1; k < M_coeffs; k++)
663  {
664  zvel = WomParam->m_zvel[i][j][k] * zt[k];
665  wbc[j] = wbc[j] + zvel.real();
666  }
667  }
668 
669  // Multiply w by normal to get u,v,w component of velocity
670  Vmath::Smul(nfq, axis_normal, wbc, 1, wbc, 1);
671  // get the offset
672  Bvals = BndCondExp->UpdateCoeffs() + BndCondExp->GetCoeff_Offset(i);
673 
674  // Push back to Coeff space
675  bc->FwdTrans(wbc, Bvals);
676  }
677 }
NekDouble m_time
Current time of simulation.
std::shared_ptr< WomersleyParams > WomersleyParamsSharedPtr

References ASSERTL1, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_time, m_womersleyParams, Vmath::Smul(), and tinysimd::sqrt().

Referenced by SetBoundaryConditions().

◆ SetZeroNormalVelocity()

void Nektar::IncNavierStokes::SetZeroNormalVelocity ( )
protected

Set Normal Velocity Component to Zero.

Definition at line 522 of file IncNavierStokes.cpp.

523 {
524  // use static trip since cannot use UserDefinedTag for zero
525  // velocity and have time dependent conditions
526  static bool Setup = false;
527 
528  if (Setup == true)
529  {
530  return;
531  }
532  Setup = true;
533 
534  size_t i, n;
535 
536  Array<OneD, Array<OneD, const SpatialDomains::BoundaryConditionShPtr>>
537  BndConds(m_spacedim);
538  Array<OneD, Array<OneD, MultiRegions::ExpListSharedPtr>> BndExp(m_spacedim);
539 
540  for (int i = 0; i < m_spacedim; ++i)
541  {
542  BndConds[i] = m_fields[m_velocity[i]]->GetBndConditions();
543  BndExp[i] = m_fields[m_velocity[i]]->GetBndCondExpansions();
544  }
545 
547 
548  size_t cnt;
549  size_t elmtid, nq, boundary;
550 
551  Array<OneD, Array<OneD, NekDouble>> normals;
552  Array<OneD, NekDouble> Bphys, Bcoeffs;
553 
554  size_t fldid = m_velocity[0];
555 
556  for (cnt = n = 0; n < BndConds[0].size(); ++n)
557  {
558  if ((BndConds[0][n]->GetBoundaryConditionType() ==
560  (boost::iequals(BndConds[0][n]->GetUserDefined(),
561  "ZeroNormalComponent")))
562  {
563  size_t nExp = BndExp[0][n]->GetExpSize();
564  for (i = 0; i < nExp; ++i, cnt++)
565  {
566  elmtid = m_fieldsBCToElmtID[fldid][cnt];
567  elmt = m_fields[0]->GetExp(elmtid);
568  boundary = m_fieldsBCToTraceID[fldid][cnt];
569 
570  normals = elmt->GetTraceNormal(boundary);
571 
572  nq = BndExp[0][n]->GetExp(i)->GetTotPoints();
573  Array<OneD, NekDouble> normvel(nq, 0.0);
574 
575  for (int k = 0; k < m_spacedim; ++k)
576  {
577  Bphys = BndExp[k][n]->UpdatePhys() +
578  BndExp[k][n]->GetPhys_Offset(i);
579  Bc = BndExp[k][n]->GetExp(i);
580  Vmath::Vvtvp(nq, normals[k], 1, Bphys, 1, normvel, 1,
581  normvel, 1);
582  }
583 
584  // negate normvel for next step
585  Vmath::Neg(nq, normvel, 1);
586 
587  for (int k = 0; k < m_spacedim; ++k)
588  {
589  Bphys = BndExp[k][n]->UpdatePhys() +
590  BndExp[k][n]->GetPhys_Offset(i);
591  Bcoeffs = BndExp[k][n]->UpdateCoeffs() +
592  BndExp[k][n]->GetCoeff_Offset(i);
593  Bc = BndExp[k][n]->GetExp(i);
594  Vmath::Vvtvp(nq, normvel, 1, normals[k], 1, Bphys, 1, Bphys,
595  1);
596  Bc->FwdTransBndConstrained(Bphys, Bcoeffs);
597  }
598  }
599  }
600  else
601  {
602  cnt += BndExp[0][n]->GetExpSize();
603  }
604  }
605 }
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:518
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

References Nektar::SpatialDomains::eDirichlet, Nektar::SolverUtils::EquationSystem::m_fields, m_fieldsBCToElmtID, m_fieldsBCToTraceID, Nektar::SolverUtils::EquationSystem::m_spacedim, m_velocity, Vmath::Neg(), and Vmath::Vvtvp().

Referenced by SetBoundaryConditions().

◆ v_GetDensity()

void Nektar::IncNavierStokes::v_GetDensity ( const Array< OneD, const Array< OneD, NekDouble >> &  physfield,
Array< OneD, NekDouble > &  density 
)
overridevirtual

Implements Nektar::SolverUtils::FluidInterface.

Definition at line 1186 of file IncNavierStokes.cpp.

1189 {
1190  int nPts = physfield[0].size();
1191  Vmath::Fill(nPts, 1.0, density, 1);
1192 }
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().

◆ v_GetForceDimension()

virtual int Nektar::IncNavierStokes::v_GetForceDimension ( )
protectedpure virtual

◆ v_GetMaxStdVelocity()

Array< OneD, NekDouble > Nektar::IncNavierStokes::v_GetMaxStdVelocity ( const NekDouble  SpeedSoundFactor)
overrideprotectedvirtual

Reimplemented from Nektar::SolverUtils::AdvectionSystem.

Definition at line 1139 of file IncNavierStokes.cpp.

1141 {
1142  boost::ignore_unused(SpeedSoundFactor);
1143  size_t nvel = m_velocity.size();
1144  size_t nelmt = m_fields[0]->GetExpSize();
1145 
1146  Array<OneD, NekDouble> stdVelocity(nelmt, 0.0);
1147  Array<OneD, Array<OneD, NekDouble>> velfields;
1148 
1149  if (m_HomogeneousType == eHomogeneous1D) // just do check on 2D info
1150  {
1151  velfields = Array<OneD, Array<OneD, NekDouble>>(2);
1152 
1153  for (size_t i = 0; i < 2; ++i)
1154  {
1155  velfields[i] = m_fields[m_velocity[i]]->UpdatePhys();
1156  }
1157  }
1158  else
1159  {
1160  velfields = Array<OneD, Array<OneD, NekDouble>>(nvel);
1161 
1162  for (size_t i = 0; i < nvel; ++i)
1163  {
1164  velfields[i] = m_fields[m_velocity[i]]->UpdatePhys();
1165  }
1166  }
1167 
1168  stdVelocity = m_extrapolation->GetMaxStdVelocity(velfields);
1169 
1170  return stdVelocity;
1171 }
ExtrapolateSharedPtr m_extrapolation
enum HomogeneousType m_HomogeneousType

References Nektar::SolverUtils::EquationSystem::eHomogeneous1D, m_extrapolation, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_HomogeneousType, and m_velocity.

◆ v_GetMovingFrameAngles()

void Nektar::IncNavierStokes::v_GetMovingFrameAngles ( Array< OneD, NekDouble > &  vFrameTheta)
overridevirtual

Function to get the angles between the moving frame of reference and stationary inertial reference frame

Reimplemented from Nektar::SolverUtils::FluidInterface.

Definition at line 1294 of file IncNavierStokes.cpp.

1296 {
1297  ASSERTL0(vFrameTheta.size() == m_movingFrameTheta.size(),
1298  "Arrays have different size, cannot get moving frame angles");
1299  for (size_t i = 0; i < m_movingFrameTheta.size(); ++i)
1300  {
1301  vFrameTheta[i] = m_movingFrameTheta[i];
1302  }
1303 }
Array< OneD, NekDouble > m_movingFrameTheta
Moving frame of reference angles with respect to the.

References ASSERTL0, and Nektar::SolverUtils::EquationSystem::m_movingFrameTheta.

◆ v_GetMovingFrameProjectionMat()

void Nektar::IncNavierStokes::v_GetMovingFrameProjectionMat ( bnu::matrix< NekDouble > &  vProjMat)
overridevirtual

Definition at line 1256 of file IncNavierStokes.cpp.

1258 {
1259  ASSERTL0(vProjMat.size1() == m_movingFrameProjMat.size1(),
1260  "Matrices have different numbers of rows, cannot Get the "
1261  "moving frame projection matrix");
1262  ASSERTL0(vProjMat.size2() == m_movingFrameProjMat.size2(),
1263  "Matrices have different numbers of columns, cannot Get the "
1264  "moving frame projection matrix");
1265 
1266  for (size_t i = 0; i < vProjMat.size1(); ++i)
1267  {
1268  for (size_t j = 0; j < vProjMat.size2(); ++j)
1269  {
1270  vProjMat(i, j) = m_movingFrameProjMat(i, j);
1271  }
1272  }
1273 }

References ASSERTL0, and Nektar::SolverUtils::EquationSystem::m_movingFrameProjMat.

◆ v_GetMovingFrameVelocities()

void Nektar::IncNavierStokes::v_GetMovingFrameVelocities ( Array< OneD, NekDouble > &  vFrameVels)
overridevirtual

Reimplemented from Nektar::SolverUtils::FluidInterface.

Definition at line 1222 of file IncNavierStokes.cpp.

1224 {
1225  ASSERTL0(vFrameVels.size() == m_movingFrameVelsxyz.size(),
1226  "Arrays have different dimensions, cannot get moving frame "
1227  "velocities");
1228  size_t size = m_movingFrameVelsxyz.size();
1229  Vmath::Vcopy(size, m_movingFrameVelsxyz, 1, vFrameVels, 1);
1230 }

References ASSERTL0, Nektar::SolverUtils::EquationSystem::m_movingFrameVelsxyz, and Vmath::Vcopy().

◆ v_GetPressure() [1/2]

virtual MultiRegions::ExpListSharedPtr Nektar::IncNavierStokes::v_GetPressure ( void  )
inlineoverrideprotectedvirtual

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 267 of file IncNavierStokes.h.

268  {
269  return m_pressure;
270  }
MultiRegions::ExpListSharedPtr m_pressure
Pointer to field holding pressure field.

References m_pressure.

◆ v_GetPressure() [2/2]

void Nektar::IncNavierStokes::v_GetPressure ( const Array< OneD, const Array< OneD, NekDouble >> &  physfield,
Array< OneD, NekDouble > &  pressure 
)
overridevirtual

Implements Nektar::SolverUtils::FluidInterface.

Definition at line 1176 of file IncNavierStokes.cpp.

1179 {
1180  pressure = physfield[m_nConvectiveFields];
1181 }

References m_nConvectiveFields, and CG_Iterations::pressure.

◆ v_GetVelocity()

void Nektar::IncNavierStokes::v_GetVelocity ( const Array< OneD, const Array< OneD, NekDouble >> &  physfield,
Array< OneD, Array< OneD, NekDouble >> &  velocity 
)
overridevirtual

Implements Nektar::SolverUtils::FluidInterface.

Definition at line 1197 of file IncNavierStokes.cpp.

1200 {
1201  for (int i = 0; i < m_spacedim; ++i)
1202  {
1203  velocity[i] = physfield[i];
1204  }
1205 }

References Nektar::SolverUtils::EquationSystem::m_spacedim.

◆ v_HasConstantDensity()

virtual bool Nektar::IncNavierStokes::v_HasConstantDensity ( )
inlineoverridevirtual

Implements Nektar::SolverUtils::FluidInterface.

Definition at line 157 of file IncNavierStokes.h.

158  {
159  return true;
160  }

◆ v_InitObject()

void Nektar::IncNavierStokes::v_InitObject ( bool  DeclareField = true)
overridevirtual

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::SolverUtils::AdvectionSystem.

Reimplemented in Nektar::VelocityCorrectionScheme, Nektar::VCSMapping, Nektar::SmoothedProfileMethod, and Nektar::CoupledLinearNS.

Definition at line 91 of file IncNavierStokes.cpp.

92 {
93  AdvectionSystem::v_InitObject(DeclareField);
94 
95  int i, j;
96  int numfields = m_fields.size();
97  std::string velids[] = {"u", "v", "w"};
98 
99  // Set up Velocity field to point to the first m_expdim of m_fields;
100  m_velocity = Array<OneD, int>(m_spacedim);
101 
102  for (i = 0; i < m_spacedim; ++i)
103  {
104  for (j = 0; j < numfields; ++j)
105  {
106  std::string var = m_boundaryConditions->GetVariable(j);
107  if (boost::iequals(velids[i], var))
108  {
109  m_velocity[i] = j;
110  break;
111  }
112 
113  ASSERTL0(j != numfields, "Failed to find field: " + var);
114  }
115  }
116 
117  // Set up equation type enum using kEquationTypeStr
118  for (i = 0; i < (int)eEquationTypeSize; ++i)
119  {
120  bool match;
121  m_session->MatchSolverInfo("EQTYPE", kEquationTypeStr[i], match, false);
122  if (match)
123  {
125  break;
126  }
127  }
128  ASSERTL0(i != eEquationTypeSize, "EQTYPE not found in SOLVERINFO section");
129 
130  m_session->LoadParameter("Kinvis", m_kinvis);
131 
132  // Default advection type per solver
133  std::string vConvectiveType;
134  switch (m_equationType)
135  {
136  case eUnsteadyStokes:
137  case eSteadyLinearisedNS:
138  vConvectiveType = "NoAdvection";
139  break;
141  case eSteadyNavierStokes:
142  vConvectiveType = "Convective";
143  break;
145  vConvectiveType = "Linearised";
146  break;
147  default:
148  break;
149  }
150 
151  // Check if advection type overridden
152  if (m_session->DefinesTag("AdvectiveType") &&
155  {
156  vConvectiveType = m_session->GetTag("AdvectiveType");
157  }
158 
159  // Initialise advection
161  vConvectiveType, vConvectiveType);
162  m_advObject->InitObject(m_session, m_fields);
163 
164  // Set up arrays for moving reference frame
165  // Note: this must be done before the forcing
166  m_movingFrameProjMat = bnu::identity_matrix<NekDouble>(3, 3);
167  if (DefinedForcing("MovingReferenceFrame"))
168  {
169  m_movingFrameVelsxyz = Array<OneD, NekDouble>(6, 0.0);
170  m_movingFrameTheta = Array<OneD, NekDouble>(3, 0.0);
171  m_pivotPoint = Array<OneD, NekDouble>(3, 0.0);
173  }
174 
175  // Forcing terms
176  m_forcing = SolverUtils::Forcing::Load(m_session, shared_from_this(),
178 
179  // check to see if any Robin boundary conditions and if so set
180  // up m_field to boundary condition maps;
181  m_fieldsBCToElmtID = Array<OneD, Array<OneD, int>>(numfields);
182  m_fieldsBCToTraceID = Array<OneD, Array<OneD, int>>(numfields);
183  m_fieldsRadiationFactor = Array<OneD, Array<OneD, NekDouble>>(numfields);
184 
185  for (size_t i = 0; i < m_fields.size(); ++i)
186  {
187  bool Set = false;
188 
189  Array<OneD, const SpatialDomains::BoundaryConditionShPtr> BndConds;
190  Array<OneD, MultiRegions::ExpListSharedPtr> BndExp;
191  int radpts = 0;
192 
193  BndConds = m_fields[i]->GetBndConditions();
194  BndExp = m_fields[i]->GetBndCondExpansions();
195  for (size_t n = 0; n < BndConds.size(); ++n)
196  {
197  if (boost::iequals(BndConds[n]->GetUserDefined(), "Radiation"))
198  {
199  ASSERTL0(
200  BndConds[n]->GetBoundaryConditionType() ==
202  "Radiation boundary condition must be of type Robin <R>");
203 
204  if (Set == false)
205  {
206  m_fields[i]->GetBoundaryToElmtMap(m_fieldsBCToElmtID[i],
208  Set = true;
209  }
210  radpts += BndExp[n]->GetTotPoints();
211  }
212  if (boost::iequals(BndConds[n]->GetUserDefined(),
213  "ZeroNormalComponent"))
214  {
215  ASSERTL0(BndConds[n]->GetBoundaryConditionType() ==
217  "Zero Normal Component boundary condition option must "
218  "be of type Dirichlet <D>");
219 
220  if (Set == false)
221  {
222  m_fields[i]->GetBoundaryToElmtMap(m_fieldsBCToElmtID[i],
224  Set = true;
225  }
226  }
227  }
228 
229  m_fieldsRadiationFactor[i] = Array<OneD, NekDouble>(radpts);
230 
231  radpts = 0; // reset to use as a counter
232 
233  for (size_t n = 0; n < BndConds.size(); ++n)
234  {
235  if (boost::iequals(BndConds[n]->GetUserDefined(), "Radiation"))
236  {
237 
238  int npoints = BndExp[n]->GetNpoints();
239  Array<OneD, NekDouble> x0(npoints, 0.0);
240  Array<OneD, NekDouble> x1(npoints, 0.0);
241  Array<OneD, NekDouble> x2(npoints, 0.0);
242  Array<OneD, NekDouble> tmpArray;
243 
244  BndExp[n]->GetCoords(x0, x1, x2);
245 
246  LibUtilities::Equation coeff =
247  std::static_pointer_cast<
248  SpatialDomains::RobinBoundaryCondition>(BndConds[n])
249  ->m_robinPrimitiveCoeff;
250 
251  coeff.Evaluate(x0, x1, x2, m_time,
252  tmpArray = m_fieldsRadiationFactor[i] + radpts);
253  // Vmath::Neg(npoints,tmpArray = m_fieldsRadiationFactor[i]+
254  // radpts,1);
255  radpts += npoints;
256  }
257  }
258  }
259 
260  // Set up maping for womersley BC - and load variables
261  for (size_t i = 0; i < m_fields.size(); ++i)
262  {
263  for (size_t n = 0; n < m_fields[i]->GetBndConditions().size(); ++n)
264  {
265  if (boost::istarts_with(
266  m_fields[i]->GetBndConditions()[n]->GetUserDefined(),
267  "Womersley"))
268  {
269  // assumes that boundary condition is applied in normal
270  // direction and is decomposed for each direction. There could
271  // be a unique file for each direction
272  m_womersleyParams[i][n] =
274  m_spacedim);
275  // Read in fourier coeffs and precompute coefficients
277  i, n, m_fields[i]->GetBndConditions()[n]->GetUserDefined());
278 
279  m_fields[i]->GetBoundaryToElmtMap(m_fieldsBCToElmtID[i],
281  }
282  }
283  }
284 
285  // Check if we have moving reference frame boundary conditions,if so, all
286  // velocity vectors at those boundaries should be Drichlet and have a same
287  // tag loop over the boundaries
288  for (size_t n = 0; n < m_fields[0]->GetBndConditions().size(); ++n)
289  {
290  // loop over the velocities
291  for (size_t i = 0; i < m_velocity.size(); ++i)
292  {
293  if (boost::iequals(
294  m_fields[i]->GetBndConditions()[n]->GetUserDefined(),
295  "MovingFrameWall"))
296  {
297  // all velocities must have the same userdefined tag and all
298  // must be Dirichlet
299  for (size_t j = 0; j < m_velocity.size(); ++j)
300  {
301  ASSERTL0(m_fields[j]
302  ->GetBndConditions()[n]
303  ->GetBoundaryConditionType() ==
305  "All velocities with userdefined tag: "
306  "\"MovingFrameWall\" must be Dirichlet boundary "
307  "condition");
308  ASSERTL0(boost::iequals(m_fields[j]
309  ->GetBndConditions()[n]
310  ->GetUserDefined(),
311  "MovingFrameWall"),
312  "If any of the velocity components at a "
313  "boundary is defined as \"MovingFrameWall\", "
314  "the rest of velocity components are also "
315  "must be defined as \"MovingFrameWall\" ");
316  }
317  }
318 
319  if (boost::iequals(
320  m_fields[i]->GetBndConditions()[n]->GetUserDefined(),
321  "MovingFrameDomainVel"))
322  {
323  // all velocities must have the same userdefined tag and all
324  // must be Dirichlet
325  for (size_t j = 0; j < m_velocity.size(); ++j)
326  {
327  ASSERTL0(
328  m_fields[j]
329  ->GetBndConditions()[n]
330  ->GetBoundaryConditionType() ==
332  "All velocities with userdefined tag: "
333  "\"MovingFrameDomainVel\" must be Dirichlet boundary "
334  "condition");
335  ASSERTL0(boost::iequals(m_fields[j]
336  ->GetBndConditions()[n]
337  ->GetUserDefined(),
338  "MovingFrameDomainVel"),
339  "If any of the velocity components at a "
340  "boundary is defined as \"MovingFrameDomainVel\", "
341  "the rest of velocity components are also "
342  "must be defined as \"MovingFrameDomainVel\" ");
343  }
344  }
345  }
346  }
347 
348  // Set up Field Meta Data for output files
349  m_fieldMetaDataMap["Kinvis"] = boost::lexical_cast<std::string>(m_kinvis);
350  m_fieldMetaDataMap["TimeStep"] =
351  boost::lexical_cast<std::string>(m_timestep);
352 
353  // Check to see if the metadata for moving reference frame is defined
354  if (DefinedForcing("MovingReferenceFrame"))
355  {
357  {
358  std::vector<std::string> vSuffix = {"_x", "_y", "_z"};
359  for (size_t i = 0; i < 3; ++i)
360  {
361  NekDouble dTheta{0};
362  std::string sTheta = "Theta" + vSuffix[i];
363  auto it = m_fieldMetaDataMap.find(sTheta);
364  if (it != m_fieldMetaDataMap.end())
365  {
366  dTheta = boost::lexical_cast<NekDouble>(it->second);
367  m_movingFrameTheta[i] = dTheta;
368  }
369  else
370  {
371  m_fieldMetaDataMap[sTheta] =
372  boost::lexical_cast<std::string>(m_movingFrameTheta[i]);
373  }
374  }
375  }
376  else
377  {
378  std::vector<std::string> vSuffix = {"_x", "_y", "_z"};
379  for (size_t i = 0; i < 3; ++i)
380  {
381  std::string sTheta = "Theta" + vSuffix[i];
382  m_fieldMetaDataMap[sTheta] =
383  boost::lexical_cast<std::string>(m_movingFrameTheta[i]);
384  }
385  }
386  }
387 }
virtual int v_GetForceDimension()=0
void GetPivotPoint(Array< OneD, NekDouble > &vPivotPoint)
void SetUpWomersley(const int fldid, const int bndid, std::string womstr)
Set Up Womersley details.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
virtual SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareField=true) override
Init object for UnsteadySystem class.
NekDouble m_timestep
Time step size.
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
Map to identify relevant solver info to dump in output fields.
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
Pointer to boundary conditions object.
static SOLVER_UTILS_EXPORT std::vector< ForcingSharedPtr > Load(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields=0)
Definition: Forcing.cpp:120
static FieldMetaDataMap NullFieldMetaDataMap
Definition: FieldIO.h:53
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47
@ eSteadyNavierStokes
@ eUnsteadyStokes
@ eUnsteadyNavierStokes
@ eSteadyLinearisedNS
@ eUnsteadyLinearisedNS
@ eEquationTypeSize
const std::string kEquationTypeStr[]

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), DefinedForcing(), Nektar::SpatialDomains::eDirichlet, Nektar::eEquationTypeSize, Nektar::SpatialDomains::eRobin, Nektar::eSteadyLinearisedNS, Nektar::eSteadyNavierStokes, Nektar::eUnsteadyLinearisedNS, Nektar::eUnsteadyNavierStokes, Nektar::eUnsteadyStokes, Nektar::LibUtilities::Equation::Evaluate(), Nektar::SolverUtils::GetAdvectionFactory(), GetPivotPoint(), Nektar::kEquationTypeStr, Nektar::SolverUtils::Forcing::Load(), Nektar::SolverUtils::AdvectionSystem::m_advObject, Nektar::SolverUtils::EquationSystem::m_boundaryConditions, m_equationType, Nektar::SolverUtils::EquationSystem::m_fieldMetaDataMap, Nektar::SolverUtils::EquationSystem::m_fields, m_fieldsBCToElmtID, m_fieldsBCToTraceID, m_fieldsRadiationFactor, m_forcing, m_kinvis, Nektar::SolverUtils::EquationSystem::m_movingFrameProjMat, Nektar::SolverUtils::EquationSystem::m_movingFrameTheta, Nektar::SolverUtils::EquationSystem::m_movingFrameVelsxyz, m_pivotPoint, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_spacedim, Nektar::SolverUtils::EquationSystem::m_time, Nektar::SolverUtils::EquationSystem::m_timestep, m_velocity, m_womersleyParams, Nektar::LibUtilities::NullFieldMetaDataMap, SetUpWomersley(), v_GetForceDimension(), and Nektar::SolverUtils::AdvectionSystem::v_InitObject().

Referenced by Nektar::CoupledLinearNS::v_InitObject(), and Nektar::VelocityCorrectionScheme::v_InitObject().

◆ v_PreIntegrate()

bool Nektar::IncNavierStokes::v_PreIntegrate ( int  step)
overrideprotectedvirtual

Perform the extrapolation.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 1397 of file IncNavierStokes.cpp.

1398 {
1399  m_extrapolation->SubStepSaveFields(step);
1400  m_extrapolation->SubStepAdvance(step, m_time);
1402  return false;
1403 }
void SetBoundaryConditions(NekDouble time)
time dependent boundary conditions updating

References m_extrapolation, Nektar::SolverUtils::EquationSystem::m_time, Nektar::SolverUtils::EquationSystem::m_timestep, and SetBoundaryConditions().

◆ v_SetMovingFrameAngles()

void Nektar::IncNavierStokes::v_SetMovingFrameAngles ( const Array< OneD, NekDouble > &  vFrameTheta)
overridevirtual

Function to set the angles between the moving frame of reference and stationary inertial reference frame

Reimplemented from Nektar::SolverUtils::FluidInterface.

Definition at line 1279 of file IncNavierStokes.cpp.

1281 {
1282  ASSERTL0(vFrameTheta.size() == m_movingFrameTheta.size(),
1283  "Arrays have different size, cannot set moving frame angles");
1284  for (size_t i = 0; i < vFrameTheta.size(); ++i)
1285  {
1286  m_movingFrameTheta[i] = vFrameTheta[i];
1287  }
1288 }

References ASSERTL0, and Nektar::SolverUtils::EquationSystem::m_movingFrameTheta.

◆ v_SetMovingFrameProjectionMat()

void Nektar::IncNavierStokes::v_SetMovingFrameProjectionMat ( const bnu::matrix< NekDouble > &  vProjMat)
overridevirtual

Function to set the rotation matrix computed in the forcing this gives access to the moving reference forcing to set the Projection matrix to be used later in IncNavierStokes calss for enforcing the boundary conditions

Definition at line 1238 of file IncNavierStokes.cpp.

1240 {
1241  ASSERTL0(vProjMat.size1() == m_movingFrameProjMat.size1(),
1242  "Matrices have different numbers of rows, cannot Set the "
1243  "moving frame projection matrix");
1244  ASSERTL0(vProjMat.size2() == m_movingFrameProjMat.size2(),
1245  "Matrices have different numbers of columns, cannot Set the "
1246  "moving frame projection matrix");
1247  for (size_t i = 0; i < vProjMat.size1(); ++i)
1248  {
1249  for (size_t j = 0; j < vProjMat.size2(); ++j)
1250  {
1251  m_movingFrameProjMat(i, j) = vProjMat(i, j);
1252  }
1253  }
1254 }

References ASSERTL0, and Nektar::SolverUtils::EquationSystem::m_movingFrameProjMat.

◆ v_SetMovingFrameVelocities()

void Nektar::IncNavierStokes::v_SetMovingFrameVelocities ( const Array< OneD, NekDouble > &  vFrameVels)
overridevirtual

Function to set the moving frame velocities calucated in the forcing this gives access to the moving reference forcing to set the velocities to be later used in enforcing the boundary condition in IncNavierStokes class

Reimplemented from Nektar::SolverUtils::FluidInterface.

Definition at line 1213 of file IncNavierStokes.cpp.

1215 {
1216  ASSERTL0(vFrameVels.size() == m_movingFrameVelsxyz.size(),
1217  "Arrays have different dimensions, cannot set moving frame "
1218  "velocities");
1219  Vmath::Vcopy(vFrameVels.size(), vFrameVels, 1, m_movingFrameVelsxyz, 1);
1220 }

References ASSERTL0, Nektar::SolverUtils::EquationSystem::m_movingFrameVelsxyz, and Vmath::Vcopy().

◆ v_TransCoeffToPhys()

virtual void Nektar::IncNavierStokes::v_TransCoeffToPhys ( void  )
inlineoverrideprotectedvirtual

Virtual function for transformation to physical space.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Reimplemented in Nektar::VelocityCorrectionScheme, and Nektar::CoupledLinearNS.

Definition at line 272 of file IncNavierStokes.h.

273  {
274  ASSERTL0(false, "This method is not defined in this class");
275  }

References ASSERTL0.

◆ v_TransPhysToCoeff()

virtual void Nektar::IncNavierStokes::v_TransPhysToCoeff ( void  )
inlineoverrideprotectedvirtual

Virtual function for transformation to coefficient space.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Reimplemented in Nektar::VelocityCorrectionScheme, and Nektar::CoupledLinearNS.

Definition at line 277 of file IncNavierStokes.h.

278  {
279  ASSERTL0(false, "This method is not defined in this class");
280  }

References ASSERTL0.

◆ WriteModalEnergy()

void Nektar::IncNavierStokes::WriteModalEnergy ( void  )
protected

Member Data Documentation

◆ eqTypeLookupIds

std::string Nektar::IncNavierStokes::eqTypeLookupIds
staticprotected
Initial value:
= {
"EqType", "SteadyNavierStokes", eSteadyNavierStokes),
"EqType", "SteadyLinearisedNS", eSteadyLinearisedNS),
"EqType", "UnsteadyNavierStokes", eUnsteadyNavierStokes),
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.

Definition at line 236 of file IncNavierStokes.h.

◆ m_energysteps

int Nektar::IncNavierStokes::m_energysteps
protected

dump energy to file at steps time

Definition at line 208 of file IncNavierStokes.h.

◆ m_equationType

EquationType Nektar::IncNavierStokes::m_equationType
protected

◆ m_extrapolation

ExtrapolateSharedPtr Nektar::IncNavierStokes::m_extrapolation
protected

◆ m_fieldsBCToElmtID

Array<OneD, Array<OneD, int> > Nektar::IncNavierStokes::m_fieldsBCToElmtID
protected

Mapping from BCs to Elmt IDs.

Definition at line 214 of file IncNavierStokes.h.

Referenced by SetRadiationBoundaryForcing(), SetZeroNormalVelocity(), and v_InitObject().

◆ m_fieldsBCToTraceID

Array<OneD, Array<OneD, int> > Nektar::IncNavierStokes::m_fieldsBCToTraceID
protected

Mapping from BCs to Elmt Edge IDs.

Definition at line 216 of file IncNavierStokes.h.

Referenced by SetRadiationBoundaryForcing(), SetZeroNormalVelocity(), and v_InitObject().

◆ m_fieldsRadiationFactor

Array<OneD, Array<OneD, NekDouble> > Nektar::IncNavierStokes::m_fieldsRadiationFactor
protected

RHS Factor for Radiation Condition.

Definition at line 218 of file IncNavierStokes.h.

Referenced by SetRadiationBoundaryForcing(), and v_InitObject().

◆ m_forcing

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

◆ m_intSteps

int Nektar::IncNavierStokes::m_intSteps
protected

Number of time integration steps AND Order of extrapolation for pressure boundary conditions.

Definition at line 222 of file IncNavierStokes.h.

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

◆ m_kinvis

NekDouble Nektar::IncNavierStokes::m_kinvis
protected

◆ m_mdlFile

std::ofstream Nektar::IncNavierStokes::m_mdlFile
protected

modal energy file

Definition at line 188 of file IncNavierStokes.h.

◆ m_nConvectiveFields

int Nektar::IncNavierStokes::m_nConvectiveFields
protected

◆ m_pivotPoint

Array<OneD, NekDouble> Nektar::IncNavierStokes::m_pivotPoint
protected

Definition at line 226 of file IncNavierStokes.h.

Referenced by SetMRFWallBCs(), and v_InitObject().

◆ m_pressure

MultiRegions::ExpListSharedPtr Nektar::IncNavierStokes::m_pressure
protected

◆ m_SmoothAdvection

bool Nektar::IncNavierStokes::m_SmoothAdvection
protected

◆ m_velocity

Array<OneD, int> Nektar::IncNavierStokes::m_velocity
protected

◆ m_womersleyParams

std::map<int, std::map<int, WomersleyParamsSharedPtr> > Nektar::IncNavierStokes::m_womersleyParams
protected

Womersley parameters if required.

Definition at line 265 of file IncNavierStokes.h.

Referenced by SetUpWomersley(), SetWomersleyBoundary(), and v_InitObject().