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

 ~IncNavierStokes () override
 
void v_InitObject (bool DeclareField=true) override
 Initialisation object for EquationSystem. More...
 
int GetNConvectiveFields (void)
 
void AddForcing (const SolverUtils::ForcingSharedPtr &pForce)
 
void v_GetPressure (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &pressure) override
 
void v_GetDensity (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &density) override
 
bool v_HasConstantDensity () override
 
void v_GetVelocity (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &velocity) override
 
void v_SetMovingFrameVelocities (const Array< OneD, NekDouble > &vFrameVels, const int step) override
 
bool v_GetMovingFrameVelocities (Array< OneD, NekDouble > &vFrameVels, const int step) override
 
void v_SetMovingFrameDisp (const Array< OneD, NekDouble > &vFrameDisp, const int step) override
 
void v_SetMovingFramePivot (const Array< OneD, NekDouble > &vFramePivot) override
 
bool v_GetMovingFrameDisp (Array< OneD, NekDouble > &vFrameDisp, const int step) override
 
void v_SetAeroForce (Array< OneD, NekDouble > forces) override
 
void v_GetAeroForce (Array< OneD, NekDouble > forces) override
 
bool DefinedForcing (const std::string &sForce)
 
- Public Member Functions inherited from Nektar::SolverUtils::AdvectionSystem
SOLVER_UTILS_EXPORT AdvectionSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
SOLVER_UTILS_EXPORT ~AdvectionSystem () override
 
SOLVER_UTILS_EXPORT void v_InitObject (bool DeclareField=true) override
 Initialisation object for EquationSystem. More...
 
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
SOLVER_UTILS_EXPORT ~UnsteadySystem () override
 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 NekDouble GetTimeStep ()
 
SOLVER_UTILS_EXPORT void SetTimeStep (const NekDouble timestep)
 
SOLVER_UTILS_EXPORT void SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
SOLVER_UTILS_EXPORT LibUtilities::TimeIntegrationSchemeSharedPtrGetTimeIntegrationScheme ()
 Returns the time integration scheme. More...
 
SOLVER_UTILS_EXPORT LibUtilities::TimeIntegrationSchemeOperatorsGetTimeIntegrationSchemeOperators ()
 Returns the time integration scheme operators. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT void InitObject (bool DeclareField=true)
 Initialises the members of this object. More...
 
SOLVER_UTILS_EXPORT void DoInitialise (bool dumpInitialConditions=true)
 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 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 NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation. 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 GetTime ()
 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 Array< OneD, NekDouble > & UpdatePhysField (const int i)
 
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 void SetIterationNumberPIT (int num)
 
SOLVER_UTILS_EXPORT void SetWindowNumberPIT (int 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...
 
SOLVER_UTILS_EXPORT bool NegatedOp ()
 Identify if operator is negated in DoSolve. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::ALEHelper
virtual ~ALEHelper ()=default
 
virtual SOLVER_UTILS_EXPORT void v_ALEInitObject (int spaceDim, Array< OneD, MultiRegions::ExpListSharedPtr > &fields)
 
SOLVER_UTILS_EXPORT void InitObject (int spaceDim, Array< OneD, MultiRegions::ExpListSharedPtr > &fields)
 
virtual SOLVER_UTILS_EXPORT void v_UpdateGridVelocity (const NekDouble &time)
 
virtual SOLVER_UTILS_EXPORT void v_ALEPreMultiplyMass (Array< OneD, Array< OneD, NekDouble > > &fields)
 
SOLVER_UTILS_EXPORT void ALEDoElmtInvMass (Array< OneD, Array< OneD, NekDouble > > &traceNormals, Array< OneD, Array< OneD, NekDouble > > &fields, NekDouble time)
 Update m_fields with u^n by multiplying by inverse mass matrix. That's then used in e.g. checkpoint output and L^2 error calculation. More...
 
SOLVER_UTILS_EXPORT void ALEDoElmtInvMassBwdTrans (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
SOLVER_UTILS_EXPORT void MoveMesh (const NekDouble &time, Array< OneD, Array< OneD, NekDouble > > &traceNormals)
 
const Array< OneD, const Array< OneD, NekDouble > > & GetGridVelocity ()
 
SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & GetGridVelocityTrace ()
 
SOLVER_UTILS_EXPORT void ExtraFldOutputGridVelocity (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 
- 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, const int step)
 
SOLVER_UTILS_EXPORT bool GetMovingFrameVelocities (Array< OneD, NekDouble > &vFrameVels, const int step)
 
SOLVER_UTILS_EXPORT void SetMovingFrameDisp (const Array< OneD, NekDouble > &vFrameDisp, const int step)
 
SOLVER_UTILS_EXPORT void SetMovingFramePivot (const Array< OneD, NekDouble > &vFramePivot)
 
SOLVER_UTILS_EXPORT bool GetMovingFrameDisp (Array< OneD, NekDouble > &vFrameDisp, const int step)
 
SOLVER_UTILS_EXPORT void SetAeroForce (Array< OneD, NekDouble > forces)
 Set aerodynamic force and moment. More...
 
SOLVER_UTILS_EXPORT void GetAeroForce (Array< OneD, NekDouble > forces)
 Get aerodynamic force and moment. More...
 

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...
 
MultiRegions::ExpListSharedPtr v_GetPressure () override
 
void v_TransCoeffToPhys (void) override
 Virtual function for transformation to physical space. More...
 
void v_TransPhysToCoeff (void) override
 Virtual function for transformation to coefficient space. More...
 
virtual int v_GetForceDimension ()=0
 
Array< OneD, NekDoublev_GetMaxStdVelocity (const NekDouble SpeedSoundFactor) override
 
bool v_PreIntegrate (int step) override
 
- Protected Member Functions inherited from Nektar::SolverUtils::AdvectionSystem
SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step) override
 
virtual SOLVER_UTILS_EXPORT Array< OneD, NekDoublev_GetMaxStdVelocity (const NekDouble SpeedSoundFactor=1.0)
 
- 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 void v_InitObject (bool DeclareField=true) override
 Init object for UnsteadySystem class. More...
 
SOLVER_UTILS_EXPORT void v_DoSolve () override
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_PrintStatusInformation (const int step, const NekDouble cpuTime)
 Print Status Information. More...
 
virtual SOLVER_UTILS_EXPORT void v_PrintSummaryStatistics (const NekDouble intTime)
 Print Summary Statistics. More...
 
SOLVER_UTILS_EXPORT void v_DoInitialise (bool dumpInitialConditions=true) override
 Sets up initial conditions. More...
 
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_PreIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT bool v_RequireFwdTrans ()
 
virtual SOLVER_UTILS_EXPORT void v_SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
virtual SOLVER_UTILS_EXPORT bool v_UpdateTimeStepCheck ()
 
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control. More...
 
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...
 
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 void v_InitObject (bool DeclareFeld=true)
 Initialisation object for EquationSystem. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise (bool dumpInitialConditions=true)
 Virtual function for initialisation implementation. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve ()
 Virtual function for solve implementation. 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_GenerateSummary (SummaryList &l)
 Virtual function for generating summary information. 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 MultiRegions::ExpListSharedPtr v_GetPressure (void)
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp (void)
 Virtual function to identify if operator is negated in DoSolve. More...
 
virtual SOLVER_UTILS_EXPORT void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 
virtual SOLVER_UTILS_EXPORT void v_GetVelocity (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &velocity)=0
 
virtual SOLVER_UTILS_EXPORT bool v_HasConstantDensity ()=0
 
virtual SOLVER_UTILS_EXPORT void v_GetDensity (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &density)=0
 
virtual SOLVER_UTILS_EXPORT void v_GetPressure (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &pressure)=0
 
virtual SOLVER_UTILS_EXPORT void v_SetMovingFrameVelocities (const Array< OneD, NekDouble > &vFrameVels, const int step)
 
virtual SOLVER_UTILS_EXPORT bool v_GetMovingFrameVelocities (Array< OneD, NekDouble > &vFrameVels, const int step)
 
virtual SOLVER_UTILS_EXPORT void v_SetMovingFrameDisp (const Array< OneD, NekDouble > &vFrameDisp, const int step)
 
virtual SOLVER_UTILS_EXPORT void v_SetMovingFramePivot (const Array< OneD, NekDouble > &vFramePivot)
 
virtual SOLVER_UTILS_EXPORT bool v_GetMovingFrameDisp (Array< OneD, NekDouble > &vFrameDisp, const int step)
 
virtual SOLVER_UTILS_EXPORT void v_SetAeroForce (Array< OneD, NekDouble > forces)
 
virtual SOLVER_UTILS_EXPORT void v_GetAeroForce (Array< OneD, NekDouble > forces)
 

Protected Attributes

ExtrapolateSharedPtr m_extrapolation
 
IncBoundaryConditionsSharedPtr m_IncNavierStokesBCs
 
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
 pivot point for moving reference frame More...
 
Array< OneD, NekDoublem_aeroForces
 
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
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
Array< OneD, Array< OneD, NekDouble > > m_previousSolution
 Storage for previous solution for steady-state check. More...
 
std::vector< int > m_intVariables
 
NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1). More...
 
NekDouble m_CFLGrowth
 CFL growth rate. More...
 
NekDouble m_CFLEnd
 Maximun cfl in cfl growth. More...
 
int m_abortSteps
 Number of steps between checks for abort conditions. More...
 
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...
 
int m_steadyStateSteps
 Check for steady state at step interval. More...
 
NekDouble m_steadyStateTol
 Tolerance to which steady state should be evaluated at. More...
 
int m_filtersInfosteps
 Number of time steps between outputting filters information. More...
 
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
 
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state. More...
 
std::ofstream m_errFile
 
NekDouble m_epsilon
 Diffusion coefficient. 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_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_iterPIT = 0
 Number of parallel-in-time time iteration. More...
 
int m_windowPIT = 0
 Index of windows for parallel-in-time 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_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_movingFrameData
 Moving reference frame status in the inertial frame X, Y, Z, Theta_x, Theta_y, Theta_z, U, V, W, Omega_x, Omega_y, Omega_z, A_x, A_y, A_z, DOmega_x, DOmega_y, DOmega_z, pivot_x, pivot_y, pivot_z. More...
 
std::vector< std::string > m_strFrameData
 variable name in m_movingFrameData 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...
 
- Protected Attributes inherited from Nektar::SolverUtils::ALEHelper
Array< OneD, MultiRegions::ExpListSharedPtrm_fieldsALE
 
Array< OneD, Array< OneD, NekDouble > > m_gridVelocity
 
Array< OneD, Array< OneD, NekDouble > > m_gridVelocityTrace
 
std::vector< ALEBaseShPtrm_ALEs
 
bool m_ALESolver = false
 
bool m_ImplicitALESolver = false
 
NekDouble m_prevStageTime = 0.0
 
int m_spaceDim
 

Static Protected Attributes

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

Additional Inherited Members

- 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 146 of file IncNavierStokes.h.

Constructor & Destructor Documentation

◆ ~IncNavierStokes()

Nektar::IncNavierStokes::~IncNavierStokes ( void  )
override

Destructor

Definition at line 300 of file IncNavierStokes.cpp.

301{
302}

◆ IncNavierStokes()

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

Constructor.

Constructor. Creates ...

Parameters

param

Definition at line 83 of file IncNavierStokes.cpp.

86 : UnsteadySystem(pSession, pGraph), AdvectionSystem(pSession, pGraph),
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 810 of file IncNavierStokes.cpp.

811{
812 m_forcing.push_back(pForce);
813}
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 994 of file IncNavierStokes.cpp.

995{
996 vector<std::string> vForceList;
997 bool hasForce{false};
998
999 if (!m_session->DefinesElement("Nektar/Forcing"))
1000 {
1001 return hasForce;
1002 }
1003
1004 TiXmlElement *vForcing = m_session->GetElement("Nektar/Forcing");
1005 if (vForcing)
1006 {
1007 TiXmlElement *vForce = vForcing->FirstChildElement("FORCE");
1008 while (vForce)
1009 {
1010 string vType = vForce->Attribute("TYPE");
1011
1012 vForceList.push_back(vType);
1013 vForce = vForce->NextSiblingElement("FORCE");
1014 }
1015 }
1016
1017 for (auto &f : vForceList)
1018 {
1019 if (boost::iequals(f, sForce))
1020 {
1021 hasForce = true;
1022 }
1023 }
1024
1025 return hasForce;
1026}
LibUtilities::SessionReaderSharedPtr m_session
The session reader.

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

Referenced by 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 307 of file IncNavierStokes.cpp.

310{
311 size_t VelDim = m_velocity.size();
312 Array<OneD, Array<OneD, NekDouble>> velocity(VelDim);
313
314 size_t npoints = m_fields[0]->GetNpoints();
315 for (size_t i = 0; i < VelDim; ++i)
316 {
317 velocity[i] = Array<OneD, NekDouble>(npoints);
318 Vmath::Vcopy(npoints, inarray[m_velocity[i]], 1, velocity[i], 1);
319 }
320 for (auto &x : m_forcing)
321 {
322 x->PreApply(m_fields, velocity, velocity, time);
323 }
324
325 m_advObject->Advect(m_nConvectiveFields, m_fields, velocity, inarray,
326 outarray, time);
327}
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.hpp:825

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(), Nektar::VCSMapping::v_EvaluateAdvection_SetPressureBCs(), and Nektar::VCSImplicit::v_EvaluateAdvection_SetPressureBCs().

◆ GetEquationType()

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

Definition at line 245 of file IncNavierStokes.h.

246 {
247 return m_equationType;
248 }
EquationType m_equationType
equation type;

References m_equationType.

◆ GetNConvectiveFields()

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

Definition at line 155 of file IncNavierStokes.h.

156 {
157 return m_nConvectiveFields;
158 }

References m_nConvectiveFields.

◆ SetBoundaryConditions()

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

time dependent boundary conditions updating

Time dependent boundary conditions updating

Definition at line 332 of file IncNavierStokes.cpp.

333{
334 size_t i, n;
335 std::string varName;
336 size_t nvariables = m_fields.size();
337
338 for (i = 0; i < nvariables; ++i)
339 {
340 for (n = 0; n < m_fields[i]->GetBndConditions().size(); ++n)
341 {
342 if (m_fields[i]->GetBndConditions()[n]->IsTimeDependent())
343 {
344 varName = m_session->GetVariable(i);
345 m_fields[i]->EvaluateBoundaryConditions(time, varName);
346 }
347 else if (boost::istarts_with(
348 m_fields[i]->GetBndConditions()[n]->GetUserDefined(),
349 "Womersley"))
350 {
352 }
353 }
354
355 // Set Radiation conditions if required
357 }
358
359 // Enforcing the boundary conditions (Inlet and wall) for the
360 // Moving reference frame
362}
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.

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

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

◆ SetRadiationBoundaryForcing()

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

Set Radiation forcing term.

Probably should be pushed back into ContField?

Definition at line 367 of file IncNavierStokes.cpp.

368{
369 size_t i, n;
370
371 Array<OneD, const SpatialDomains::BoundaryConditionShPtr> BndConds;
372 Array<OneD, MultiRegions::ExpListSharedPtr> BndExp;
373
374 BndConds = m_fields[fieldid]->GetBndConditions();
375 BndExp = m_fields[fieldid]->GetBndCondExpansions();
376
379
380 size_t cnt;
381 size_t elmtid, nq, offset, boundary;
382 Array<OneD, NekDouble> Bvals, U;
383 size_t cnt1 = 0;
384
385 for (cnt = n = 0; n < BndConds.size(); ++n)
386 {
387 std::string type = BndConds[n]->GetUserDefined();
388
389 if ((BndConds[n]->GetBoundaryConditionType() ==
391 (boost::iequals(type, "Radiation")))
392 {
393 size_t nExp = BndExp[n]->GetExpSize();
394 for (i = 0; i < nExp; ++i, cnt++)
395 {
396 elmtid = m_fieldsBCToElmtID[m_velocity[fieldid]][cnt];
397 elmt = m_fields[fieldid]->GetExp(elmtid);
398 offset = m_fields[fieldid]->GetPhys_Offset(elmtid);
399
400 U = m_fields[fieldid]->UpdatePhys() + offset;
401 Bc = BndExp[n]->GetExp(i);
402
403 boundary = m_fieldsBCToTraceID[fieldid][cnt];
404
405 // Get edge values and put into ubc
406 nq = Bc->GetTotPoints();
407 Array<OneD, NekDouble> ubc(nq);
408 elmt->GetTracePhysVals(boundary, Bc, U, ubc);
409
410 Vmath::Vmul(nq,
412 [fieldid][cnt1 + BndExp[n]->GetPhys_Offset(i)],
413 1, &ubc[0], 1, &ubc[0], 1);
414
415 Bvals =
416 BndExp[n]->UpdateCoeffs() + BndExp[n]->GetCoeff_Offset(i);
417
418 Bc->IProductWRTBase(ubc, Bvals);
419 }
420 cnt1 += BndExp[n]->GetTotPoints();
421 }
422 else
423 {
424 cnt += BndExp[n]->GetExpSize();
425 }
426 }
427}
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:66
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
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.hpp:72

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 586 of file IncNavierStokes.cpp.

588{
589 std::string::size_type indxBeg = womStr.find_first_of(':') + 1;
590 string filename = womStr.substr(indxBeg, string::npos);
591
592 TiXmlDocument doc(filename);
593
594 bool loadOkay = doc.LoadFile();
595 ASSERTL0(loadOkay,
596 (std::string("Failed to load file: ") + filename).c_str());
597
598 TiXmlHandle docHandle(&doc);
599
600 int err; /// Error value returned by TinyXML.
601
602 TiXmlElement *nektar = doc.FirstChildElement("NEKTAR");
603 ASSERTL0(nektar, "Unable to find NEKTAR tag in file.");
604
605 TiXmlElement *wombc = nektar->FirstChildElement("WOMERSLEYBC");
606 ASSERTL0(wombc, "Unable to find WOMERSLEYBC tag in file.");
607
608 // read womersley parameters
609 TiXmlElement *womparam = wombc->FirstChildElement("WOMPARAMS");
610 ASSERTL0(womparam, "Unable to find WOMPARAMS tag in file.");
611
612 // Input coefficients
613 TiXmlElement *params = womparam->FirstChildElement("W");
614 map<std::string, std::string> Wparams;
615
616 // read parameter list
617 while (params)
618 {
619
620 std::string propstr;
621 propstr = params->Attribute("PROPERTY");
622
623 ASSERTL0(!propstr.empty(),
624 "Failed to read PROPERTY value Womersley BC Parameter");
625
626 std::string valstr;
627 valstr = params->Attribute("VALUE");
628
629 ASSERTL0(!valstr.empty(),
630 "Failed to read VALUE value Womersley BC Parameter");
631
632 std::transform(propstr.begin(), propstr.end(), propstr.begin(),
633 ::toupper);
634 Wparams[propstr] = valstr;
635
636 params = params->NextSiblingElement("W");
637 }
638 bool parseGood;
639
640 // Read parameters
641
642 ASSERTL0(
643 Wparams.count("RADIUS") == 1,
644 "Failed to find Radius parameter in Womersley boundary conditions");
645 std::vector<NekDouble> rad;
646 ParseUtils::GenerateVector(Wparams["RADIUS"], rad);
647 m_womersleyParams[fldid][bndid]->m_radius = rad[0];
648
649 ASSERTL0(
650 Wparams.count("PERIOD") == 1,
651 "Failed to find period parameter in Womersley boundary conditions");
652 std::vector<NekDouble> period;
653 parseGood = ParseUtils::GenerateVector(Wparams["PERIOD"], period);
654 m_womersleyParams[fldid][bndid]->m_period = period[0];
655
656 ASSERTL0(
657 Wparams.count("AXISNORMAL") == 1,
658 "Failed to find axisnormal parameter in Womersley boundary conditions");
659 std::vector<NekDouble> anorm;
660 parseGood = ParseUtils::GenerateVector(Wparams["AXISNORMAL"], anorm);
661 m_womersleyParams[fldid][bndid]->m_axisnormal[0] = anorm[0];
662 m_womersleyParams[fldid][bndid]->m_axisnormal[1] = anorm[1];
663 m_womersleyParams[fldid][bndid]->m_axisnormal[2] = anorm[2];
664
665 ASSERTL0(
666 Wparams.count("AXISPOINT") == 1,
667 "Failed to find axispoint parameter in Womersley boundary conditions");
668 std::vector<NekDouble> apt;
669 parseGood = ParseUtils::GenerateVector(Wparams["AXISPOINT"], apt);
670 m_womersleyParams[fldid][bndid]->m_axispoint[0] = apt[0];
671 m_womersleyParams[fldid][bndid]->m_axispoint[1] = apt[1];
672 m_womersleyParams[fldid][bndid]->m_axispoint[2] = apt[2];
673
674 // Read Temporal Fourier Coefficients.
675
676 // Find the FourierCoeff tag
677 TiXmlElement *coeff = wombc->FirstChildElement("FOURIERCOEFFS");
678
679 // Input coefficients
680 TiXmlElement *fval = coeff->FirstChildElement("F");
681
682 int indx;
683
684 while (fval)
685 {
686 TiXmlAttribute *fvalAttr = fval->FirstAttribute();
687 std::string attrName(fvalAttr->Name());
688
689 ASSERTL0(attrName == "ID",
690 (std::string("Unknown attribute name: ") + attrName).c_str());
691
692 err = fvalAttr->QueryIntValue(&indx);
693 ASSERTL0(err == TIXML_SUCCESS, "Unable to read attribute ID.");
694
695 std::string coeffStr = fval->FirstChild()->ToText()->ValueStr();
696 vector<NekDouble> coeffvals;
697
698 parseGood = ParseUtils::GenerateVector(coeffStr, coeffvals);
699 ASSERTL0(
700 parseGood,
701 (std::string("Problem reading value of fourier coefficient, ID=") +
702 boost::lexical_cast<string>(indx))
703 .c_str());
704 ASSERTL1(
705 coeffvals.size() == 2,
706 (std::string(
707 "Have not read two entries of Fourier coefficicent from ID=" +
708 boost::lexical_cast<string>(indx))
709 .c_str()));
710
711 m_womersleyParams[fldid][bndid]->m_wom_vel.push_back(
712 NekComplexDouble(coeffvals[0], coeffvals[1]));
713
714 fval = fval->NextSiblingElement("F");
715 }
716
717 // starting point of precalculation
718 size_t i, j, k;
719 // M fourier coefficients
720 size_t M_coeffs = m_womersleyParams[fldid][bndid]->m_wom_vel.size();
721 NekDouble R = m_womersleyParams[fldid][bndid]->m_radius;
722 NekDouble T = m_womersleyParams[fldid][bndid]->m_period;
723 Array<OneD, NekDouble> x0 = m_womersleyParams[fldid][bndid]->m_axispoint;
724
726 // Womersley Number
727 NekComplexDouble omega_c(2.0 * M_PI / T, 0.0);
728 NekComplexDouble alpha_c(R * sqrt(omega_c.real() / m_kinvis), 0.0);
729 NekComplexDouble z1(1.0, 0.0);
730 NekComplexDouble i_pow_3q2(-1.0 / sqrt(2.0), 1.0 / sqrt(2.0));
731
733 BndCondExp = m_fields[fldid]->GetBndCondExpansions()[bndid];
734
736 size_t cnt = 0;
737 size_t nfq;
738 Array<OneD, NekDouble> Bvals;
739
740 size_t exp_npts = BndCondExp->GetExpSize();
741 Array<OneD, NekDouble> wbc(exp_npts, 0.0);
742
743 // allocate time indepedent variables
744 m_womersleyParams[fldid][bndid]->m_poiseuille =
745 Array<OneD, Array<OneD, NekDouble>>(exp_npts);
746 m_womersleyParams[fldid][bndid]->m_zvel =
747 Array<OneD, Array<OneD, Array<OneD, NekComplexDouble>>>(exp_npts);
748 // could use M_coeffs - 1 but need to avoid complicating things
749 Array<OneD, NekComplexDouble> zJ0(M_coeffs);
750 Array<OneD, NekComplexDouble> lamda_n(M_coeffs);
751 Array<OneD, NekComplexDouble> k_c(M_coeffs);
752 NekComplexDouble zJ0r;
753
754 for (k = 1; k < M_coeffs; k++)
755 {
756 k_c[k] = NekComplexDouble((NekDouble)k, 0.0);
757 lamda_n[k] = i_pow_3q2 * alpha_c * sqrt(k_c[k]);
758 zJ0[k] = Polylib::ImagBesselComp(0, lamda_n[k]);
759 }
760
761 // Loop over each element in an expansion
762 for (i = 0; i < exp_npts; ++i, cnt++)
763 {
764 // Get Boundary and trace expansion
765 bc = BndCondExp->GetExp(i);
766 nfq = bc->GetTotPoints();
767
768 Array<OneD, NekDouble> x(nfq, 0.0);
769 Array<OneD, NekDouble> y(nfq, 0.0);
770 Array<OneD, NekDouble> z(nfq, 0.0);
771 bc->GetCoords(x, y, z);
772
773 m_womersleyParams[fldid][bndid]->m_poiseuille[i] =
774 Array<OneD, NekDouble>(nfq);
775 m_womersleyParams[fldid][bndid]->m_zvel[i] =
776 Array<OneD, Array<OneD, NekComplexDouble>>(nfq);
777
778 // Compute coefficients
779 for (j = 0; j < nfq; j++)
780 {
781 rqR = NekComplexDouble(sqrt((x[j] - x0[0]) * (x[j] - x0[0]) +
782 (y[j] - x0[1]) * (y[j] - x0[1]) +
783 (z[j] - x0[2]) * (z[j] - x0[2])) /
784 R,
785 0.0);
786
787 // Compute Poiseulle Flow
788 m_womersleyParams[fldid][bndid]->m_poiseuille[i][j] =
789 m_womersleyParams[fldid][bndid]->m_wom_vel[0].real() *
790 (1. - rqR.real() * rqR.real());
791
792 m_womersleyParams[fldid][bndid]->m_zvel[i][j] =
793 Array<OneD, NekComplexDouble>(M_coeffs);
794
795 // compute the velocity information
796 for (k = 1; k < M_coeffs; k++)
797 {
798 zJ0r = Polylib::ImagBesselComp(0, rqR * lamda_n[k]);
799 m_womersleyParams[fldid][bndid]->m_zvel[i][j][k] =
800 m_womersleyParams[fldid][bndid]->m_wom_vel[k] *
801 (z1 - (zJ0r / zJ0[k]));
802 }
803 }
804 }
805}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
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:130
static NekDouble rad(NekDouble x, NekDouble y)
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::vector< double > z(NPUPPER)
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:1559
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(), tinysimd::sqrt(), and Nektar::UnitTests::z().

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 517 of file IncNavierStokes.cpp.

518{
519 ASSERTL1(m_womersleyParams.count(bndid) == 1,
520 "Womersley parameters for this boundary have not been set up");
521
522 WomersleyParamsSharedPtr WomParam = m_womersleyParams[fldid][bndid];
523 NekComplexDouble zvel;
524 size_t i, j, k;
525
526 size_t M_coeffs = WomParam->m_wom_vel.size();
527
528 NekDouble T = WomParam->m_period;
529 NekDouble axis_normal = WomParam->m_axisnormal[fldid];
530
531 // Womersley Number
532 NekComplexDouble omega_c(2.0 * M_PI / T, 0.0);
533 NekComplexDouble k_c(0.0, 0.0);
534 NekComplexDouble m_time_c(m_time, 0.0);
535 NekComplexDouble zi(0.0, 1.0);
536 NekComplexDouble i_pow_3q2(-1.0 / sqrt(2.0), 1.0 / sqrt(2.0));
537
539 BndCondExp = m_fields[fldid]->GetBndCondExpansions()[bndid];
540
542 size_t cnt = 0;
543 size_t nfq;
544 Array<OneD, NekDouble> Bvals;
545 size_t exp_npts = BndCondExp->GetExpSize();
546 Array<OneD, NekDouble> wbc(exp_npts, 0.0);
547
548 Array<OneD, NekComplexDouble> zt(M_coeffs);
549
550 // preallocate the exponent
551 for (k = 1; k < M_coeffs; k++)
552 {
553 k_c = NekComplexDouble((NekDouble)k, 0.0);
554 zt[k] = std::exp(zi * omega_c * k_c * m_time_c);
555 }
556
557 // Loop over each element in an expansion
558 for (i = 0; i < exp_npts; ++i, cnt++)
559 {
560 // Get Boundary and trace expansion
561 bc = BndCondExp->GetExp(i);
562 nfq = bc->GetTotPoints();
563 Array<OneD, NekDouble> wbc(nfq, 0.0);
564
565 // Compute womersley solution
566 for (j = 0; j < nfq; j++)
567 {
568 wbc[j] = WomParam->m_poiseuille[i][j];
569 for (k = 1; k < M_coeffs; k++)
570 {
571 zvel = WomParam->m_zvel[i][j][k] * zt[k];
572 wbc[j] = wbc[j] + zvel.real();
573 }
574 }
575
576 // Multiply w by normal to get u,v,w component of velocity
577 Vmath::Smul(nfq, axis_normal, wbc, 1, wbc, 1);
578 // get the offset
579 Bvals = BndCondExp->UpdateCoeffs() + BndCondExp->GetCoeff_Offset(i);
580
581 // Push back to Coeff space
582 bc->FwdTrans(wbc, Bvals);
583 }
584}
NekDouble m_time
Current time of simulation.
std::shared_ptr< WomersleyParams > WomersleyParamsSharedPtr
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.hpp:100

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 429 of file IncNavierStokes.cpp.

430{
431 // use static trip since cannot use UserDefinedTag for zero
432 // velocity and have time dependent conditions
433 static bool Setup = false;
434
435 if (Setup == true)
436 {
437 return;
438 }
439 Setup = true;
440
441 size_t i, n;
442
443 Array<OneD, Array<OneD, const SpatialDomains::BoundaryConditionShPtr>>
444 BndConds(m_spacedim);
445 Array<OneD, Array<OneD, MultiRegions::ExpListSharedPtr>> BndExp(m_spacedim);
446
447 for (int i = 0; i < m_spacedim; ++i)
448 {
449 BndConds[i] = m_fields[m_velocity[i]]->GetBndConditions();
450 BndExp[i] = m_fields[m_velocity[i]]->GetBndCondExpansions();
451 }
452
454
455 size_t cnt;
456 size_t elmtid, nq, boundary;
457
458 Array<OneD, Array<OneD, NekDouble>> normals;
459 Array<OneD, NekDouble> Bphys, Bcoeffs;
460
461 size_t fldid = m_velocity[0];
462
463 for (cnt = n = 0; n < BndConds[0].size(); ++n)
464 {
465 if ((BndConds[0][n]->GetBoundaryConditionType() ==
467 (boost::iequals(BndConds[0][n]->GetUserDefined(),
468 "ZeroNormalComponent")))
469 {
470 size_t nExp = BndExp[0][n]->GetExpSize();
471 for (i = 0; i < nExp; ++i, cnt++)
472 {
473 elmtid = m_fieldsBCToElmtID[fldid][cnt];
474 elmt = m_fields[0]->GetExp(elmtid);
475 boundary = m_fieldsBCToTraceID[fldid][cnt];
476
477 normals = elmt->GetTraceNormal(boundary);
478
479 nq = BndExp[0][n]->GetExp(i)->GetTotPoints();
480 Array<OneD, NekDouble> normvel(nq, 0.0);
481
482 for (int k = 0; k < m_spacedim; ++k)
483 {
484 Bphys = BndExp[k][n]->UpdatePhys() +
485 BndExp[k][n]->GetPhys_Offset(i);
486 Bc = BndExp[k][n]->GetExp(i);
487 Vmath::Vvtvp(nq, normals[k], 1, Bphys, 1, normvel, 1,
488 normvel, 1);
489 }
490
491 // negate normvel for next step
492 Vmath::Neg(nq, normvel, 1);
493
494 for (int k = 0; k < m_spacedim; ++k)
495 {
496 Bphys = BndExp[k][n]->UpdatePhys() +
497 BndExp[k][n]->GetPhys_Offset(i);
498 Bcoeffs = BndExp[k][n]->UpdateCoeffs() +
499 BndExp[k][n]->GetCoeff_Offset(i);
500 Bc = BndExp[k][n]->GetExp(i);
501 Vmath::Vvtvp(nq, normvel, 1, normals[k], 1, Bphys, 1, Bphys,
502 1);
503 Bc->FwdTransBndConstrained(Bphys, Bcoeffs);
504 }
505 }
506 }
507 else
508 {
509 cnt += BndExp[0][n]->GetExpSize();
510 }
511 }
512}
int m_spacedim
Spatial dimension (>= expansion dim).
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.hpp:292
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.hpp:366

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_GetAeroForce()

void Nektar::IncNavierStokes::v_GetAeroForce ( Array< OneD, NekDouble forces)
overridevirtual

Reimplemented from Nektar::SolverUtils::FluidInterface.

Definition at line 983 of file IncNavierStokes.cpp.

984{
985 if (m_aeroForces.size() >= 6)
986 {
987 Vmath::Vcopy(6, m_aeroForces, 1, forces, 1);
988 }
989}
Array< OneD, NekDouble > m_aeroForces

References m_aeroForces, and Vmath::Vcopy().

◆ 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 867 of file IncNavierStokes.cpp.

870{
871 int nPts = physfield[0].size();
872 Vmath::Fill(nPts, 1.0, density, 1);
873}
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.hpp:54

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 818 of file IncNavierStokes.cpp.

820{
821 size_t nvel = m_velocity.size();
822 size_t nelmt = m_fields[0]->GetExpSize();
823
824 Array<OneD, NekDouble> stdVelocity(nelmt, 0.0);
825 Array<OneD, Array<OneD, NekDouble>> velfields;
826
827 if (m_HomogeneousType == eHomogeneous1D) // just do check on 2D info
828 {
829 velfields = Array<OneD, Array<OneD, NekDouble>>(2);
830
831 for (size_t i = 0; i < 2; ++i)
832 {
833 velfields[i] = m_fields[m_velocity[i]]->UpdatePhys();
834 }
835 }
836 else
837 {
838 velfields = Array<OneD, Array<OneD, NekDouble>>(nvel);
839
840 for (size_t i = 0; i < nvel; ++i)
841 {
842 velfields[i] = m_fields[m_velocity[i]]->UpdatePhys();
843 }
844 }
845
846 stdVelocity = m_extrapolation->GetMaxStdVelocity(velfields);
847
848 return stdVelocity;
849}
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_GetMovingFrameDisp()

bool Nektar::IncNavierStokes::v_GetMovingFrameDisp ( Array< OneD, NekDouble > &  vFrameDisp,
const int  step 
)
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 946 of file IncNavierStokes.cpp.

948{
949 if (m_movingFrameData.size())
950 {
951 ASSERTL0(
952 vFrameDisp.size() == 6,
953 "Arrays have different size, cannot get moving frame displacement");
954 Vmath::Vcopy(vFrameDisp.size(), m_movingFrameData + 21 * step, 1,
955 vFrameDisp, 1);
956 return true;
957 }
958 else
959 {
960 return false;
961 }
962}
Array< OneD, NekDouble > m_movingFrameData
Moving reference frame status in the inertial frame X, Y, Z, Theta_x, Theta_y, Theta_z,...

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

◆ v_GetMovingFrameVelocities()

bool Nektar::IncNavierStokes::v_GetMovingFrameVelocities ( Array< OneD, NekDouble > &  vFrameVels,
const int  step 
)
overridevirtual

Reimplemented from Nektar::SolverUtils::FluidInterface.

Definition at line 907 of file IncNavierStokes.cpp.

909{
910 if (m_movingFrameData.size())
911 {
912 ASSERTL0(vFrameVels.size() == 12,
913 "Arrays have different dimensions, cannot get moving frame "
914 "velocities");
915 Vmath::Vcopy(vFrameVels.size(), m_movingFrameData + 6 + 21 * step, 1,
916 vFrameVels, 1);
917 return true;
918 }
919 else
920 {
921 return false;
922 }
923}

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

◆ v_GetPressure() [1/2]

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

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 275 of file IncNavierStokes.h.

276 {
277 return m_pressure;
278 }
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 854 of file IncNavierStokes.cpp.

857{
858 if (physfield.size())
859 {
860 pressure = physfield[physfield.size() - 1];
861 }
862}

References 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 878 of file IncNavierStokes.cpp.

881{
882 for (int i = 0; i < m_spacedim; ++i)
883 {
884 velocity[i] = physfield[i];
885 }
886}

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

◆ v_HasConstantDensity()

bool Nektar::IncNavierStokes::v_HasConstantDensity ( )
inlineoverridevirtual

Implements Nektar::SolverUtils::FluidInterface.

Definition at line 170 of file IncNavierStokes.h.

171 {
172 return true;
173 }

◆ v_InitObject()

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

Initialisation object for EquationSystem.

Continuous field

Setting up the normals

Setting up the normals

Reimplemented from Nektar::SolverUtils::AdvectionSystem.

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

Definition at line 91 of file IncNavierStokes.cpp.

92{
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:
138 vConvectiveType = "NoAdvection";
139 break;
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 if (DefinedForcing("MovingReferenceFrame"))
167 {
168 // 0-5(inertial disp), 6-11(inertial vel), 12-17(inertial acce) current
169 // 18-21(body pivot)
170 // 21-26(inertial disp), 27-32(body vel), 33-38(body acce) next step
171 // 39-41(body pivot)
173 "X", "Y", "Z", "Theta_x", "Theta_y", "Theta_z",
174 "U", "V", "W", "Omega_x", "Omega_y", "Omega_z",
175 "A_x", "A_y", "A_z", "DOmega_x", "DOmega_y", "DOmega_z",
176 "X0", "Y0", "Z0"};
177 m_movingFrameData = Array<OneD, NekDouble>(42, 0.0);
178 m_aeroForces = Array<OneD, NekDouble>(6, 0.0);
179 }
180
181 // Forcing terms
182 m_forcing = SolverUtils::Forcing::Load(m_session, shared_from_this(),
184
185 // check to see if any Robin boundary conditions and if so set
186 // up m_field to boundary condition maps;
187 m_fieldsBCToElmtID = Array<OneD, Array<OneD, int>>(numfields);
188 m_fieldsBCToTraceID = Array<OneD, Array<OneD, int>>(numfields);
189 m_fieldsRadiationFactor = Array<OneD, Array<OneD, NekDouble>>(numfields);
190
191 for (size_t i = 0; i < m_fields.size(); ++i)
192 {
193 bool Set = false;
194
195 Array<OneD, const SpatialDomains::BoundaryConditionShPtr> BndConds;
196 Array<OneD, MultiRegions::ExpListSharedPtr> BndExp;
197 int radpts = 0;
198
199 BndConds = m_fields[i]->GetBndConditions();
200 BndExp = m_fields[i]->GetBndCondExpansions();
201 for (size_t n = 0; n < BndConds.size(); ++n)
202 {
203 if (boost::iequals(BndConds[n]->GetUserDefined(), "Radiation"))
204 {
205 ASSERTL0(
206 BndConds[n]->GetBoundaryConditionType() ==
208 "Radiation boundary condition must be of type Robin <R>");
209
210 if (Set == false)
211 {
212 m_fields[i]->GetBoundaryToElmtMap(m_fieldsBCToElmtID[i],
214 Set = true;
215 }
216 radpts += BndExp[n]->GetTotPoints();
217 }
218 if (boost::iequals(BndConds[n]->GetUserDefined(),
219 "ZeroNormalComponent"))
220 {
221 ASSERTL0(BndConds[n]->GetBoundaryConditionType() ==
223 "Zero Normal Component boundary condition option must "
224 "be of type Dirichlet <D>");
225
226 if (Set == false)
227 {
228 m_fields[i]->GetBoundaryToElmtMap(m_fieldsBCToElmtID[i],
230 Set = true;
231 }
232 }
233 }
234
235 m_fieldsRadiationFactor[i] = Array<OneD, NekDouble>(radpts);
236
237 radpts = 0; // reset to use as a counter
238
239 for (size_t n = 0; n < BndConds.size(); ++n)
240 {
241 if (boost::iequals(BndConds[n]->GetUserDefined(), "Radiation"))
242 {
243
244 int npoints = BndExp[n]->GetNpoints();
245 Array<OneD, NekDouble> x0(npoints, 0.0);
246 Array<OneD, NekDouble> x1(npoints, 0.0);
247 Array<OneD, NekDouble> x2(npoints, 0.0);
248 Array<OneD, NekDouble> tmpArray;
249
250 BndExp[n]->GetCoords(x0, x1, x2);
251
252 LibUtilities::Equation coeff =
253 std::static_pointer_cast<
254 SpatialDomains::RobinBoundaryCondition>(BndConds[n])
255 ->m_robinPrimitiveCoeff;
256
257 coeff.Evaluate(x0, x1, x2, m_time,
258 tmpArray = m_fieldsRadiationFactor[i] + radpts);
259 // Vmath::Neg(npoints,tmpArray = m_fieldsRadiationFactor[i]+
260 // radpts,1);
261 radpts += npoints;
262 }
263 }
264 }
265
266 // Set up maping for womersley BC - and load variables
267 for (size_t i = 0; i < m_fields.size(); ++i)
268 {
269 for (size_t n = 0; n < m_fields[i]->GetBndConditions().size(); ++n)
270 {
271 if (boost::istarts_with(
272 m_fields[i]->GetBndConditions()[n]->GetUserDefined(),
273 "Womersley"))
274 {
275 // assumes that boundary condition is applied in normal
276 // direction and is decomposed for each direction. There could
277 // be a unique file for each direction
278 m_womersleyParams[i][n] =
280 m_spacedim);
281 // Read in fourier coeffs and precompute coefficients
283 i, n, m_fields[i]->GetBndConditions()[n]->GetUserDefined());
284
285 m_fields[i]->GetBoundaryToElmtMap(m_fieldsBCToElmtID[i],
287 }
288 }
289 }
290
291 // Set up Field Meta Data for output files
292 m_fieldMetaDataMap["Kinvis"] = boost::lexical_cast<std::string>(m_kinvis);
293 m_fieldMetaDataMap["TimeStep"] =
294 boost::lexical_cast<std::string>(m_timestep);
295}
virtual int v_GetForceDimension()=0
bool DefinedForcing(const std::string &sForce)
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.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareField=true) override
Initialisation object for EquationSystem.
NekDouble m_timestep
Time step size.
std::vector< std::string > m_strFrameData
variable name in m_movingFrameData
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:118
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:43
@ 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(), Nektar::kEquationTypeStr, Nektar::SolverUtils::Forcing::Load(), Nektar::SolverUtils::AdvectionSystem::m_advObject, m_aeroForces, 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_movingFrameData, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_spacedim, Nektar::SolverUtils::EquationSystem::m_strFrameData, Nektar::SolverUtils::EquationSystem::m_time, Nektar::SolverUtils::EquationSystem::m_timestep, m_velocity, m_womersleyParams, 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 1031 of file IncNavierStokes.cpp.

1032{
1033 m_extrapolation->SubStepSaveFields(step);
1034 m_extrapolation->SubStepAdvance(step, m_time);
1036 return false;
1037}
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_SetAeroForce()

void Nektar::IncNavierStokes::v_SetAeroForce ( Array< OneD, NekDouble forces)
overridevirtual

Reimplemented from Nektar::SolverUtils::FluidInterface.

Definition at line 975 of file IncNavierStokes.cpp.

976{
977 if (m_aeroForces.size() >= 6)
978 {
979 Vmath::Vcopy(6, forces, 1, m_aeroForces, 1);
980 }
981}

References m_aeroForces, and Vmath::Vcopy().

◆ v_SetMovingFrameDisp()

void Nektar::IncNavierStokes::v_SetMovingFrameDisp ( const Array< OneD, NekDouble > &  vFrameDisp,
const int  step 
)
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 929 of file IncNavierStokes.cpp.

931{
932 if (m_movingFrameData.size())
933 {
934 ASSERTL0(
935 vFrameDisp.size() == 6,
936 "Arrays have different size, cannot set moving frame displacement");
937 Array<OneD, NekDouble> temp = m_movingFrameData + 21 * step;
938 Vmath::Vcopy(vFrameDisp.size(), vFrameDisp, 1, temp, 1);
939 }
940}

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

◆ v_SetMovingFramePivot()

void Nektar::IncNavierStokes::v_SetMovingFramePivot ( const Array< OneD, NekDouble > &  vFramePivot)
overridevirtual

Reimplemented from Nektar::SolverUtils::FluidInterface.

Definition at line 964 of file IncNavierStokes.cpp.

966{
967 ASSERTL0(vFramePivot.size() == 3,
968 "Arrays have different size, cannot set moving frame pivot");
969 Array<OneD, NekDouble> temp = m_movingFrameData + 18;
970 Vmath::Vcopy(vFramePivot.size(), vFramePivot, 1, temp, 1);
971 temp = m_movingFrameData + 39;
972 Vmath::Vcopy(vFramePivot.size(), vFramePivot, 1, temp, 1);
973}

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

◆ v_SetMovingFrameVelocities()

void Nektar::IncNavierStokes::v_SetMovingFrameVelocities ( const Array< OneD, NekDouble > &  vFrameVels,
const int  step 
)
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 894 of file IncNavierStokes.cpp.

896{
897 if (m_movingFrameData.size())
898 {
899 ASSERTL0(vFrameVels.size() == 12,
900 "Arrays have different dimensions, cannot set moving frame "
901 "velocities");
902 Array<OneD, NekDouble> temp = m_movingFrameData + 6 + 21 * step;
903 Vmath::Vcopy(vFrameVels.size(), vFrameVels, 1, temp, 1);
904 }
905}

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

◆ v_TransCoeffToPhys()

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

Virtual function for transformation to physical space.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Reimplemented in Nektar::VelocityCorrectionScheme.

Definition at line 280 of file IncNavierStokes.h.

281 {
282 ASSERTL0(false, "This method is not defined in this class");
283 }

References ASSERTL0.

◆ v_TransPhysToCoeff()

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

Virtual function for transformation to coefficient space.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Reimplemented in Nektar::VelocityCorrectionScheme.

Definition at line 285 of file IncNavierStokes.h.

286 {
287 ASSERTL0(false, "This method is not defined in this class");
288 }

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 249 of file IncNavierStokes.h.

◆ m_aeroForces

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

Definition at line 239 of file IncNavierStokes.h.

Referenced by v_GetAeroForce(), v_InitObject(), and v_SetAeroForce().

◆ m_energysteps

int Nektar::IncNavierStokes::m_energysteps
protected

dump energy to file at steps time

Definition at line 221 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 227 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 229 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 231 of file IncNavierStokes.h.

Referenced by SetRadiationBoundaryForcing(), and v_InitObject().

◆ m_forcing

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

◆ m_IncNavierStokesBCs

IncBoundaryConditionsSharedPtr Nektar::IncNavierStokes::m_IncNavierStokesBCs
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 235 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 201 of file IncNavierStokes.h.

◆ m_nConvectiveFields

int Nektar::IncNavierStokes::m_nConvectiveFields
protected

◆ m_pivotPoint

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

pivot point for moving reference frame

Definition at line 238 of file IncNavierStokes.h.

◆ 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

int which identifies which components of m_fields contains the velocity (u,v,w);

Definition at line 214 of file IncNavierStokes.h.

Referenced by Nektar::VCSImplicit::AddImplicitSkewSymAdvection(), Nektar::CoupledLinearNS::Continuation(), Nektar::CoupledLinearNS::DefineForcingTerm(), EvaluateAdvectionTerms(), Nektar::CoupledLinearNS::EvaluateNewtonRHS(), Nektar::CoupledLinearNS::InfNorm(), Nektar::CoupledLinearNS::L2Norm(), Nektar::SmoothedProfileMethod::ReadPhi(), Nektar::SmoothedProfileMethod::SetCorrectionPressureBCs(), SetRadiationBoundaryForcing(), Nektar::SmoothedProfileMethod::SetUpCorrectionPressure(), Nektar::CoupledLinearNS::SetUpCoupledMatrix(), Nektar::SmoothedProfileMethod::SetUpExpansions(), Nektar::VelocityCorrectionScheme::SetUpExtrapolation(), Nektar::VelocityCorrectionScheme::SetUpSVV(), SetZeroNormalVelocity(), Nektar::CoupledLinearNS::Solve(), Nektar::SmoothedProfileMethod::SolveCorrectedVelocity(), Nektar::CoupledLinearNS::SolveLinearNS(), Nektar::CoupledLinearNS::SolveSteadyNavierStokes(), Nektar::CoupledLinearNS::SolveUnsteadyStokesSystem(), Nektar::SmoothedProfileMethod::UpdateForcing(), Nektar::VCSImplicit::v_DoInitialise(), Nektar::CoupledLinearNS::v_DoInitialise(), Nektar::VCSImplicit::v_EvaluateAdvection_SetPressureBCs(), v_GetMaxStdVelocity(), Nektar::CoupledLinearNS::v_InitObject(), v_InitObject(), Nektar::SmoothedProfileMethod::v_InitObject(), Nektar::VCSMapping::v_InitObject(), Nektar::VelocityCorrectionScheme::v_SetUpPressureForcing(), Nektar::VCSImplicit::v_SetUpPressureForcing(), Nektar::VelocityCorrectionScheme::v_SetUpViscousForcing(), Nektar::VCSMapping::v_SetUpViscousForcing(), Nektar::VCSImplicit::v_SetUpViscousForcing(), and Nektar::VCSImplicit::v_SolveViscous().

◆ m_womersleyParams

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

Womersley parameters if required.

Definition at line 273 of file IncNavierStokes.h.

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