Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Protected Member Functions | 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:
Inheritance graph
[legend]
Collaboration diagram for Nektar::IncNavierStokes:
Collaboration graph
[legend]

Public Member Functions

virtual ~IncNavierStokes ()
 
virtual void v_InitObject ()
 Init object for UnsteadySystem class. More...
 
virtual void v_GetFluxVector (const int i, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
 
virtual void v_NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numflux)
 
int GetNConvectiveFields (void)
 
Array< OneD, int > & GetVelocity (void)
 
Array< OneD, NekDoubleGetElmtCFLVals (void)
 
NekDouble GetCFLEstimate (int &elmtid)
 
void AddForcing (const SolverUtils::ForcingSharedPtr &pForce)
 
- Public Member Functions inherited from Nektar::SolverUtils::AdvectionSystem
SOLVER_UTILS_EXPORT AdvectionSystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 
virtual SOLVER_UTILS_EXPORT ~AdvectionSystem ()
 
AdvectionSharedPtr GetAdvObject ()
 Returns the advection object held by this instance. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
virtual SOLVER_UTILS_EXPORT ~UnsteadySystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT NekDouble GetTimeStep (const Array< OneD, const Array< OneD, NekDouble > > &inarray)
 Calculate the larger time-step mantaining the problem stable. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT void SetUpTraceNormals (void)
 
SOLVER_UTILS_EXPORT void InitObject ()
 Initialises the members of this object. More...
 
SOLVER_UTILS_EXPORT void DoInitialise ()
 Perform any initialisation necessary before solving the problem. More...
 
SOLVER_UTILS_EXPORT void DoSolve ()
 Solve the problem. More...
 
SOLVER_UTILS_EXPORT void TransCoeffToPhys ()
 Transform from coefficient to physical space. More...
 
SOLVER_UTILS_EXPORT void TransPhysToCoeff ()
 Transform from physical to coefficient space. More...
 
SOLVER_UTILS_EXPORT void Output ()
 Perform output operations after solve. More...
 
SOLVER_UTILS_EXPORT NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation. More...
 
SOLVER_UTILS_EXPORT std::string GetSessionName ()
 Get Session name. More...
 
template<class T >
boost::shared_ptr< T > as ()
 
SOLVER_UTILS_EXPORT void ResetSessionName (std::string newname)
 Reset Session name. More...
 
SOLVER_UTILS_EXPORT
LibUtilities::SessionReaderSharedPtr 
GetSession ()
 Get Session name. More...
 
SOLVER_UTILS_EXPORT
MultiRegions::ExpListSharedPtr 
GetPressure ()
 Get pressure field if available. More...
 
SOLVER_UTILS_EXPORT void PrintSummary (std::ostream &out)
 Print a summary of parameters and solver characteristics. More...
 
SOLVER_UTILS_EXPORT void SetLambda (NekDouble lambda)
 Set parameter m_lambda. More...
 
SOLVER_UTILS_EXPORT void EvaluateFunction (Array< OneD, Array< OneD, NekDouble > > &pArray, std::string pFunctionName, const NekDouble pTime=0.0, const int domain=0)
 Evaluates a function as specified in the session file. More...
 
SOLVER_UTILS_EXPORT void EvaluateFunction (std::vector< std::string > pFieldNames, Array< OneD, Array< OneD, NekDouble > > &pFields, const std::string &pName, const NekDouble &pTime=0.0, const int domain=0)
 Populate given fields with the function from session. More...
 
SOLVER_UTILS_EXPORT void EvaluateFunction (std::vector< std::string > pFieldNames, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const std::string &pName, const NekDouble &pTime=0.0, const int domain=0)
 Populate given fields with the function from session. More...
 
SOLVER_UTILS_EXPORT void EvaluateFunction (std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, const NekDouble &pTime=0.0, const int domain=0)
 
SOLVER_UTILS_EXPORT std::string DescribeFunction (std::string pFieldName, const std::string &pFunctionName, const int domain)
 Provide a description of a function for a given field name. More...
 
SOLVER_UTILS_EXPORT void InitialiseBaseFlow (Array< OneD, Array< OneD, NekDouble > > &base)
 Perform initialisation of the base flow. More...
 
SOLVER_UTILS_EXPORT void SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 Initialise the data in the dependent fields. More...
 
SOLVER_UTILS_EXPORT void EvaluateExactSolution (int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 Evaluates an exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln, bool Normalised=false)
 Compute the L2 error between fields and a given exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, bool Normalised=false)
 Compute the L2 error of the fields. More...
 
SOLVER_UTILS_EXPORT Array
< OneD, NekDouble
ErrorExtraPoints (unsigned int field)
 Compute error (L2 and L_inf) over an larger set of quadrature points return [L2 Linf]. More...
 
SOLVER_UTILS_EXPORT void WeakAdvectionGreensDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\nabla \phi \cdot F) $. More...
 
SOLVER_UTILS_EXPORT void WeakAdvectionDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\phi, \nabla \cdot F) $. More...
 
SOLVER_UTILS_EXPORT void WeakAdvectionNonConservativeForm (const Array< OneD, Array< OneD, NekDouble > > &V, const Array< OneD, const NekDouble > &u, Array< OneD, NekDouble > &outarray, bool UseContCoeffs=false)
 Compute the inner product $ (\phi, V\cdot \nabla u) $. More...
 
f SOLVER_UTILS_EXPORT void AdvectionNonConservativeForm (const Array< OneD, Array< OneD, NekDouble > > &V, const Array< OneD, const NekDouble > &u, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wk=NullNekDouble1DArray)
 Compute the non-conservative advection. More...
 
SOLVER_UTILS_EXPORT void WeakDGAdvection (const Array< OneD, Array< OneD, NekDouble > > &InField, Array< OneD, Array< OneD, NekDouble > > &OutField, bool NumericalFluxIncludesNormal=true, bool InFieldIsInPhysSpace=false, int nvariables=0)
 Calculate the weak discontinuous Galerkin advection. More...
 
SOLVER_UTILS_EXPORT void WeakDGDiffusion (const Array< OneD, Array< OneD, NekDouble > > &InField, Array< OneD, Array< OneD, NekDouble > > &OutField, bool NumericalFluxIncludesNormal=true, bool InFieldIsInPhysSpace=false)
 Calculate weak DG Diffusion in the LDG form. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n)
 Write checkpoint file of m_fields. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 Write checkpoint file of custom data fields. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_BaseFlow (const int n)
 Write base flow file of m_fields. More...
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname)
 Write field data to the given filename. More...
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 Write input fields to the given filename. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 Input field data from the given file. More...
 
SOLVER_UTILS_EXPORT void ImportFldToMultiDomains (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int ndomains)
 Input field data from the given file to multiple domains. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, std::vector< std::string > &fieldStr, Array< OneD, Array< OneD, NekDouble > > &coeffs)
 Output a field. Input field data into array from the given file. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, MultiRegions::ExpListSharedPtr &pField, std::string &pFieldName)
 Output a field. Input field data into ExpList from the given file. More...
 
SOLVER_UTILS_EXPORT void ScanForHistoryPoints ()
 Builds map of which element holds each history point. More...
 
SOLVER_UTILS_EXPORT void WriteHistoryData (std::ostream &out)
 Probe each history point and write to file. More...
 
SOLVER_UTILS_EXPORT void SessionSummary (SummaryList &vSummary)
 Write out a session summary. More...
 
SOLVER_UTILS_EXPORT Array
< OneD,
MultiRegions::ExpListSharedPtr > & 
UpdateFields ()
 
SOLVER_UTILS_EXPORT
LibUtilities::FieldMetaDataMap
UpdateFieldMetaDataMap ()
 Get hold of FieldInfoMap so it can be updated. More...
 
SOLVER_UTILS_EXPORT NekDouble GetFinalTime ()
 Return final time. More...
 
SOLVER_UTILS_EXPORT int GetNcoeffs ()
 
SOLVER_UTILS_EXPORT int GetNcoeffs (const int eid)
 
SOLVER_UTILS_EXPORT int GetNumExpModes ()
 
SOLVER_UTILS_EXPORT const
Array< OneD, int > 
GetNumExpModesPerExp ()
 
SOLVER_UTILS_EXPORT int GetNvariables ()
 
SOLVER_UTILS_EXPORT const
std::string 
GetVariable (unsigned int i)
 
SOLVER_UTILS_EXPORT int GetTraceTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTraceNpoints ()
 
SOLVER_UTILS_EXPORT int GetExpSize ()
 
SOLVER_UTILS_EXPORT int GetPhys_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetCoeff_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTotPoints (int n)
 
SOLVER_UTILS_EXPORT int GetNpoints ()
 
SOLVER_UTILS_EXPORT int GetNumElmVelocity ()
 
SOLVER_UTILS_EXPORT int GetSteps ()
 
SOLVER_UTILS_EXPORT NekDouble GetTimeStep ()
 
SOLVER_UTILS_EXPORT void CopyFromPhysField (const int i, Array< OneD, NekDouble > &output)
 
SOLVER_UTILS_EXPORT void CopyToPhysField (const int i, Array< OneD, NekDouble > &output)
 
SOLVER_UTILS_EXPORT void SetSteps (const int steps)
 
SOLVER_UTILS_EXPORT void ZeroPhysFields ()
 
SOLVER_UTILS_EXPORT void FwdTransFields ()
 
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
 
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &fluxX, Array< OneD, Array< OneD, NekDouble > > &fluxY)
 
SOLVER_UTILS_EXPORT void GetFluxVector (const int i, const int j, Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &flux)
 
SOLVER_UTILS_EXPORT void NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numflux)
 
SOLVER_UTILS_EXPORT void NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxX, Array< OneD, Array< OneD, NekDouble > > &numfluxY)
 
SOLVER_UTILS_EXPORT void NumFluxforScalar (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &uflux)
 
SOLVER_UTILS_EXPORT void NumFluxforVector (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
 
SOLVER_UTILS_EXPORT void SetModifiedBasis (const bool modbasis)
 
SOLVER_UTILS_EXPORT int NoCaseStringCompare (const std::string &s1, const std::string &s2)
 Perform a case-insensitive string comparison. More...
 
SOLVER_UTILS_EXPORT int GetCheckpointNumber ()
 
SOLVER_UTILS_EXPORT void SetCheckpointNumber (int num)
 
SOLVER_UTILS_EXPORT int GetCheckpointSteps ()
 
SOLVER_UTILS_EXPORT void SetCheckpointSteps (int num)
 
SOLVER_UTILS_EXPORT void SetTime (const NekDouble time)
 
SOLVER_UTILS_EXPORT void SetInitialStep (const int step)
 
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time. More...
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp ()
 Virtual function to identify if operator is negated in DoSolve. More...
 

Protected Member Functions

 IncNavierStokes (const LibUtilities::SessionReaderSharedPtr &pSession)
 Constructor. More...
 
EquationType GetEquationType (void)
 
void EvaluateAdvectionTerms (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, Array< OneD, NekDouble > &wk=NullNekDouble1DArray)
 
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...
 
bool CalcSteadyState (void)
 evaluate steady state More...
 
virtual
MultiRegions::ExpListSharedPtr 
v_GetPressure ()
 
virtual void v_TransCoeffToPhys (void)
 Virtual function for transformation to physical space. More...
 
virtual void v_TransPhysToCoeff (void)
 Virtual function for transformation to coefficient space. More...
 
virtual int v_GetForceDimension ()=0
 
virtual bool v_PreIntegrate (int step)
 
virtual bool v_PostIntegrate (int step)
 
- Protected Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
SOLVER_UTILS_EXPORT UnsteadySystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 Initialises UnsteadySystem class members. More...
 
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve ()
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise ()
 Sets up initial conditions. More...
 
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary (SummaryList &s)
 Print a summary of time stepping parameters. More...
 
virtual SOLVER_UTILS_EXPORT void v_AppendOutput1D (Array< OneD, Array< OneD, NekDouble > > &solution1D)
 Print the solution at each solution point in a txt file. More...
 
virtual SOLVER_UTILS_EXPORT void v_NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxX, Array< OneD, Array< OneD, NekDouble > > &numfluxY)
 
virtual SOLVER_UTILS_EXPORT void v_NumFluxforScalar (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &uflux)
 
virtual SOLVER_UTILS_EXPORT void v_NumFluxforVector (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
 
virtual SOLVER_UTILS_EXPORT
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_SteadyStateCheck (int step)
 
SOLVER_UTILS_EXPORT void CheckForRestartTime (NekDouble &time)
 
SOLVER_UTILS_EXPORT void SVVVarDiffCoeff (const Array< OneD, Array< OneD, NekDouble > > vel, StdRegions::VarCoeffMap &varCoeffMap)
 Evaluate the SVV diffusion coefficient according to Moura's paper where it should proportional to h time velocity. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession)
 Initialises EquationSystem class members. More...
 
int nocase_cmp (const std::string &s1, const std::string &s2)
 
virtual SOLVER_UTILS_EXPORT
NekDouble 
v_LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Virtual function for the L_inf error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT
NekDouble 
v_L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false)
 Virtual function for the L_2 error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT void v_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)
 
SOLVER_UTILS_EXPORT void SetUpBaseFields (SpatialDomains::MeshGraphSharedPtr &mesh)
 
SOLVER_UTILS_EXPORT void ImportFldBase (std::string pInfile, SpatialDomains::MeshGraphSharedPtr pGraph)
 
virtual SOLVER_UTILS_EXPORT void v_Output (void)
 
virtual SOLVER_UTILS_EXPORT void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 

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::ForcingSharedPtr
m_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...
 
int m_cflsteps
 dump cfl estimate 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...
 
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...
 
- Protected Attributes inherited from Nektar::SolverUtils::AdvectionSystem
SolverUtils::AdvectionSharedPtr m_advObject
 Advection term. More...
 
- Protected Attributes inherited from Nektar::SolverUtils::UnsteadySystem
int m_infosteps
 Number of time steps between outputting status information. More...
 
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
LibUtilities::TimeIntegrationSolutionSharedPtr m_intSoln
 
NekDouble m_epsilon
 
bool m_explicitDiffusion
 Indicates if explicit or implicit treatment of diffusion is used. More...
 
bool m_explicitAdvection
 Indicates if explicit or implicit treatment of advection is used. More...
 
bool m_explicitReaction
 Indicates if explicit or implicit treatment of reaction is used. More...
 
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state. More...
 
std::vector< int > m_intVariables
 
std::vector< FilterSharedPtrm_filters
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output. More...
 
std::map< std::string, Array
< OneD, Array< OneD, float > > > 
m_interpWeights
 Map of the interpolation weights for a specific filename. More...
 
std::map< std::string, Array
< OneD, Array< OneD, unsigned
int > > > 
m_interpInds
 Map of the interpolation indices for a specific filename. More...
 
Array< OneD,
MultiRegions::ExpListSharedPtr
m_fields
 Array holding all dependent variables. More...
 
Array< OneD,
MultiRegions::ExpListSharedPtr
m_base
 Base fields. More...
 
Array< OneD,
MultiRegions::ExpListSharedPtr
m_derivedfields
 Array holding all dependent variables. More...
 
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
 Pointer to boundary conditions object. More...
 
SpatialDomains::MeshGraphSharedPtr m_graph
 Pointer to graph defining mesh. More...
 
std::string m_sessionName
 Name of the session. More...
 
NekDouble m_time
 Current time of simulation. More...
 
int m_initialStep
 Number of the step where the simulation should begin. More...
 
NekDouble m_fintime
 Finish time of the simulation. More...
 
NekDouble m_timestep
 Time step size. More...
 
NekDouble m_lambda
 Lambda constant in real system if one required. More...
 
std::set< std::string > m_loadedFields
 
NekDouble m_checktime
 Time between checkpoints. More...
 
int m_nchk
 Number of checkpoints written so far. More...
 
int m_steps
 Number of steps to take. More...
 
int m_checksteps
 Number of steps between checkpoints. More...
 
int m_spacedim
 Spatial dimension (>= expansion dim). More...
 
int m_expdim
 Expansion dimension. More...
 
bool m_singleMode
 Flag to determine if single homogeneous mode is used. More...
 
bool m_halfMode
 Flag to determine if half homogeneous mode is used. More...
 
bool m_multipleModes
 Flag to determine if use multiple homogenenous modes are used. More...
 
bool m_useFFT
 Flag to determine if FFT is used for homogeneous transform. More...
 
bool m_homogen_dealiasing
 Flag to determine if dealiasing is used for homogeneous simulations. More...
 
bool m_specHP_dealiasing
 Flag to determine if dealisising is usde for the Spectral/hp element discretisation. More...
 
enum MultiRegions::ProjectionType m_projectionType
 Type of projection; e.g continuous or discontinuous. More...
 
Array< OneD, Array< OneD,
NekDouble > > 
m_traceNormals
 Array holding trace normals for DG simulations in the forwards direction. More...
 
Array< OneD, Array< OneD,
Array< OneD, NekDouble > > > 
m_gradtan
 1 x nvariable x nq More...
 
Array< OneD, Array< OneD,
Array< OneD, NekDouble > > > 
m_tanbasis
 2 x m_spacedim x nq More...
 
Array< OneD, bool > m_checkIfSystemSingular
 Flag to indicate if the fields should be checked for singularity. More...
 
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
 Map to identify relevant solver info to dump in output fields. More...
 
int m_NumQuadPointsError
 Number of Quadrature points used to work out the error. More...
 
enum HomogeneousType m_HomogeneousType
 
NekDouble m_LhomX
 physical length in X direction (if homogeneous) More...
 
NekDouble m_LhomY
 physical length in Y direction (if homogeneous) More...
 
NekDouble m_LhomZ
 physical length in Z direction (if homogeneous) More...
 
int m_npointsX
 number of points in X direction (if homogeneous) More...
 
int m_npointsY
 number of points in Y direction (if homogeneous) More...
 
int m_npointsZ
 number of points in Z direction (if homogeneous) More...
 
int m_HomoDirec
 number of homogenous directions More...
 

Additional Inherited Members

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

Detailed Description

This class is the base class for Navier Stokes problems.

Definition at line 104 of file IncNavierStokes.h.

Constructor & Destructor Documentation

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

Destructor

Definition at line 245 of file IncNavierStokes.cpp.

246  {
247  }
Nektar::IncNavierStokes::IncNavierStokes ( const LibUtilities::SessionReaderSharedPtr pSession)
protected

Constructor.

Constructor. Creates ...

Parameters
param

Definition at line 58 of file IncNavierStokes.cpp.

58  :
59  UnsteadySystem(pSession),
60  AdvectionSystem(pSession),
61  m_SmoothAdvection(false),
63  {
64  }
bool m_SmoothAdvection
bool to identify if advection term smoothing is requested
int m_steadyStateSteps
Check for steady state at step interval.
SOLVER_UTILS_EXPORT AdvectionSystem(const LibUtilities::SessionReaderSharedPtr &pSession)
SOLVER_UTILS_EXPORT UnsteadySystem(const LibUtilities::SessionReaderSharedPtr &pSession)
Initialises UnsteadySystem class members.

Member Function Documentation

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

Add an additional forcing term programmatically.

Definition at line 528 of file IncNavierStokes.cpp.

References m_forcing.

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

529  {
530  m_forcing.push_back(pForce);
531  }
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Forcing terms.
bool Nektar::IncNavierStokes::CalcSteadyState ( void  )
protected

evaluate steady state

Decide if at a steady state if the discrerte L2 sum of the coefficients is the same as the previous step to within the tolerance m_steadyStateTol;

Definition at line 539 of file IncNavierStokes.cpp.

References Vmath::Dot(), Nektar::SolverUtils::EquationSystem::m_fields, and m_steadyStateTol.

Referenced by v_PostIntegrate().

540  {
541  static NekDouble previousL2 = 0.0;
542  bool returnval = false;
543 
544  NekDouble L2 = 0.0;
545 
546  // calculate L2 discrete summation
547  int ncoeffs = m_fields[0]->GetNcoeffs();
548 
549  for(int i = 0; i < m_fields.num_elements(); ++i)
550  {
551  L2 += Vmath::Dot(ncoeffs,m_fields[i]->GetCoeffs(),1,m_fields[i]->GetCoeffs(),1);
552  }
553 
554  if(fabs(L2-previousL2) < ncoeffs*m_steadyStateTol)
555  {
556  returnval = true;
557  }
558 
559  previousL2 = L2;
560 
561  return returnval;
562  }
double NekDouble
T Dot(int n, const T *w, const T *x)
vvtvp (vector times vector times vector): z = w*x*y
Definition: Vmath.cpp:900
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
NekDouble m_steadyStateTol
Tolerance to which steady state should be evaluated at.
void Nektar::IncNavierStokes::EvaluateAdvectionTerms ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
Array< OneD, NekDouble > &  wk = NullNekDouble1DArray 
)
protected

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

Definition at line 314 of file IncNavierStokes.cpp.

References ASSERTL0, Nektar::SolverUtils::AdvectionSystem::m_advObject, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_halfMode, m_nConvectiveFields, Nektar::SolverUtils::EquationSystem::m_singleMode, Nektar::SolverUtils::EquationSystem::m_time, and m_velocity.

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

317  {
318  int i;
319  int nqtot = m_fields[0]->GetTotPoints();
320  int VelDim = m_velocity.num_elements();
321  Array<OneD, Array<OneD, NekDouble> > velocity(VelDim);
322  Array<OneD, NekDouble > Deriv;
323 
324  for(i = 0; i < VelDim; ++i)
325  {
326  if(m_fields[i]->GetWaveSpace() && !m_singleMode && !m_halfMode)
327  {
328  velocity[i] = Array<OneD, NekDouble>(nqtot,0.0);
329  m_fields[i]->HomogeneousBwdTrans(inarray[m_velocity[i]],velocity[i]);
330  }
331  else
332  {
333  velocity[i] = inarray[m_velocity[i]];
334  }
335  }
336 
337  // Set up Derivative work space;
338  if(wk.num_elements())
339  {
340  ASSERTL0(wk.num_elements() >= nqtot*VelDim,
341  "Workspace is not sufficient");
342  Deriv = wk;
343  }
344  else
345  {
346  Deriv = Array<OneD, NekDouble> (nqtot*VelDim);
347  }
348 
350  velocity, inarray, outarray, m_time);
351  }
bool m_singleMode
Flag to determine if single homogeneous mode is used.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
NekDouble m_time
Current time of simulation.
bool m_halfMode
Flag to determine if half homogeneous mode is used.
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;.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
NekDouble Nektar::IncNavierStokes::GetCFLEstimate ( int &  elmtid)

Definition at line 614 of file IncNavierStokes.cpp.

References Nektar::SolverUtils::EquationSystem::eHomogeneous1D, GetElmtCFLVals(), Vmath::Imax(), Nektar::SolverUtils::EquationSystem::m_comm, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_HomogeneousType, and Nektar::LibUtilities::ReduceMax.

Referenced by v_PostIntegrate().

615  {
616  int n_element = m_fields[0]->GetExpSize();
617 
618  Array<OneD, NekDouble> cfl = GetElmtCFLVals();
619 
620  elmtid = Vmath::Imax(n_element,cfl,1);
621  NekDouble CFL,CFL_loc;
622 
623  CFL = CFL_loc = cfl[elmtid];
624  m_comm->AllReduce(CFL,LibUtilities::ReduceMax);
625 
626  // unshuffle elmt id if data is not stored in consecutive order.
627  elmtid = m_fields[0]->GetExp(elmtid)->GetGeom()->GetGlobalID();
628  if(CFL != CFL_loc)
629  {
630  elmtid = -1;
631  }
632 
633  m_comm->AllReduce(elmtid,LibUtilities::ReduceMax);
634 
635  // express element id with respect to plane
637  {
638  elmtid = elmtid%m_fields[0]->GetPlane(0)->GetExpSize();
639  }
640  return CFL;
641  }
int Imax(int n, const T *x, const int incx)
Return the index of the maximum element in x.
Definition: Vmath.cpp:741
LibUtilities::CommSharedPtr m_comm
Communicator.
double NekDouble
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
enum HomogeneousType m_HomogeneousType
Array< OneD, NekDouble > GetElmtCFLVals(void)
Array< OneD, NekDouble > Nektar::IncNavierStokes::GetElmtCFLVals ( void  )

Definition at line 567 of file IncNavierStokes.cpp.

References Nektar::SolverUtils::EquationSystem::eHomogeneous1D, Nektar::SolverUtils::EquationSystem::GetNumExpModesPerExp(), m_extrapolation, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_HomogeneousType, Nektar::SolverUtils::EquationSystem::m_timestep, and m_velocity.

Referenced by GetCFLEstimate().

568  {
569  int n_vel = m_velocity.num_elements();
570  int n_element = m_fields[0]->GetExpSize();
571 
572  const Array<OneD, int> ExpOrder = GetNumExpModesPerExp();
573  Array<OneD, int> ExpOrderList (n_element, ExpOrder);
574 
575  const NekDouble cLambda = 0.2; // Spencer book pag. 317
576 
577  Array<OneD, NekDouble> cfl (n_element, 0.0);
578  Array<OneD, NekDouble> stdVelocity(n_element, 0.0);
579  Array<OneD, Array<OneD, NekDouble> > velfields;
580 
581  if(m_HomogeneousType == eHomogeneous1D) // just do check on 2D info
582  {
583  velfields = Array<OneD, Array<OneD, NekDouble> >(2);
584 
585  for(int i = 0; i < 2; ++i)
586  {
587  velfields[i] = m_fields[m_velocity[i]]->UpdatePhys();
588  }
589  }
590  else
591  {
592  velfields = Array<OneD, Array<OneD, NekDouble> >(n_vel);
593 
594  for(int i = 0; i < n_vel; ++i)
595  {
596  velfields[i] = m_fields[m_velocity[i]]->UpdatePhys();
597  }
598  }
599 
600  stdVelocity = m_extrapolation->GetMaxStdVelocity(velfields);
601 
602  for(int el = 0; el < n_element; ++el)
603  {
604  cfl[el] = m_timestep*(stdVelocity[el] * cLambda *
605  (ExpOrder[el]-1) * (ExpOrder[el]-1));
606  }
607 
608  return cfl;
609  }
NekDouble m_timestep
Time step size.
SOLVER_UTILS_EXPORT const Array< OneD, int > GetNumExpModesPerExp()
Array< OneD, int > m_velocity
int which identifies which components of m_fields contains the velocity (u,v,w);
ExtrapolateSharedPtr m_extrapolation
double NekDouble
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
enum HomogeneousType m_HomogeneousType
EquationType Nektar::IncNavierStokes::GetEquationType ( void  )
inlineprotected

Definition at line 191 of file IncNavierStokes.h.

References m_equationType.

192  {
193  return m_equationType;
194  }
EquationType m_equationType
equation type;
int Nektar::IncNavierStokes::GetNConvectiveFields ( void  )
inline

Definition at line 122 of file IncNavierStokes.h.

References m_nConvectiveFields.

123  {
124  return m_nConvectiveFields;
125  }
int m_nConvectiveFields
Number of fields to be convected;.
Array<OneD, int>& Nektar::IncNavierStokes::GetVelocity ( void  )
inline

Definition at line 127 of file IncNavierStokes.h.

References m_velocity.

128  {
129  return m_velocity;
130  }
Array< OneD, int > m_velocity
int which identifies which components of m_fields contains the velocity (u,v,w);
void Nektar::IncNavierStokes::SetBoundaryConditions ( NekDouble  time)
protected

time dependent boundary conditions updating

Time dependent boundary conditions updating

Definition at line 356 of file IncNavierStokes.cpp.

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

Referenced by Nektar::VelocityCorrectionScheme::SolveUnsteadyStokesSystem(), Nektar::CoupledLinearNS::SolveUnsteadyStokesSystem(), and Nektar::VelocityCorrectionScheme::v_DoInitialise().

357  {
358  int i, n;
359  std::string varName;
360  int nvariables = m_fields.num_elements();
361 
362  for (i = 0; i < nvariables; ++i)
363  {
364  for(n = 0; n < m_fields[i]->GetBndConditions().num_elements(); ++n)
365  {
366  if(m_fields[i]->GetBndConditions()[n]->IsTimeDependent())
367  {
368  varName = m_session->GetVariable(i);
369  m_fields[i]->EvaluateBoundaryConditions(time, varName);
370  }
371 
372  }
373 
374  // Set Radiation conditions if required
376  }
377 
379  }
void SetZeroNormalVelocity()
Set Normal Velocity Component to Zero.
void SetRadiationBoundaryForcing(int fieldid)
Set Radiation forcing term.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
void Nektar::IncNavierStokes::SetRadiationBoundaryForcing ( int  fieldid)
protected

Set Radiation forcing term.

Probably should be pushed back into ContField?

Definition at line 384 of file IncNavierStokes.cpp.

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

385  {
386  int i,n;
387 
388  Array<OneD, const SpatialDomains::BoundaryConditionShPtr > BndConds;
389  Array<OneD, MultiRegions::ExpListSharedPtr> BndExp;
390 
391  BndConds = m_fields[fieldid]->GetBndConditions();
392  BndExp = m_fields[fieldid]->GetBndCondExpansions();
393 
396 
397  int cnt;
398  int elmtid,nq,offset, boundary;
399  Array<OneD, NekDouble> Bvals, U;
400  int cnt1 = 0;
401 
402  for(cnt = n = 0; n < BndConds.num_elements(); ++n)
403  {
404  std::string type = BndConds[n]->GetUserDefined();
405 
406  if((BndConds[n]->GetBoundaryConditionType() == SpatialDomains::eRobin)&&(boost::iequals(type,"Radiation")))
407  {
408  for(i = 0; i < BndExp[n]->GetExpSize(); ++i,cnt++)
409  {
410  elmtid = m_fieldsBCToElmtID[m_velocity[fieldid]][cnt];
411  elmt = m_fields[fieldid]->GetExp(elmtid);
412  offset = m_fields[fieldid]->GetPhys_Offset(elmtid);
413 
414  U = m_fields[fieldid]->UpdatePhys() + offset;
415  Bc = BndExp[n]->GetExp(i);
416 
417  boundary = m_fieldsBCToTraceID[fieldid][cnt];
418 
419  // Get edge values and put into ubc
420  nq = Bc->GetTotPoints();
421  Array<OneD, NekDouble> ubc(nq);
422  elmt->GetTracePhysVals(boundary,Bc,U,ubc);
423 
424  Vmath::Vmul(nq,&m_fieldsRadiationFactor[fieldid][cnt1 +
425  BndExp[n]->GetPhys_Offset(i)],1,&ubc[0],1,&ubc[0],1);
426 
427  Bvals = BndExp[n]->UpdateCoeffs()+BndExp[n]->GetCoeff_Offset(i);
428 
429  Bc->IProductWRTBase(ubc,Bvals);
430  }
431  cnt1 += BndExp[n]->GetTotPoints();
432  }
433  else
434  {
435  cnt += BndExp[n]->GetExpSize();
436  }
437  }
438  }
Array< OneD, Array< OneD, int > > m_fieldsBCToTraceID
Mapping from BCs to Elmt Edge IDs.
Array< OneD, int > m_velocity
int which identifies which components of m_fields contains the velocity (u,v,w);
Array< OneD, Array< OneD, NekDouble > > m_fieldsRadiationFactor
RHS Factor for Radiation Condition.
SOLVER_UTILS_EXPORT int GetPhys_Offset(int n)
Array< OneD, Array< OneD, int > > m_fieldsBCToElmtID
Mapping from BCs to Elmt IDs.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
boost::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.cpp:169
void Nektar::IncNavierStokes::SetZeroNormalVelocity ( )
protected

Set Normal Velocity Component to Zero.

Definition at line 441 of file IncNavierStokes.cpp.

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

442  {
443  // use static trip since cannot use UserDefinedTag for zero
444  // velocity and have time dependent conditions
445  static bool Setup = false;
446 
447  if(Setup == true)
448  {
449  return;
450  }
451  Setup = true;
452 
453  int i,n;
454 
455  Array<OneD, Array<OneD, const SpatialDomains::BoundaryConditionShPtr > >
456  BndConds(m_spacedim);
457  Array<OneD, Array<OneD, MultiRegions::ExpListSharedPtr> >
458  BndExp(m_spacedim);
459 
460 
461  for(i = 0; i < m_spacedim; ++i)
462  {
463  BndConds[i] = m_fields[m_velocity[i]]->GetBndConditions();
464  BndExp[i] = m_fields[m_velocity[i]]->GetBndCondExpansions();
465  }
466 
468 
469  int cnt;
470  int elmtid,nq, boundary;
471 
472  Array<OneD, Array<OneD, NekDouble> > normals;
473  Array<OneD, NekDouble> Bphys,Bcoeffs;
474 
475  int fieldid = m_velocity[0];
476 
477  for(cnt = n = 0; n < BndConds[0].num_elements(); ++n)
478  {
479  if((BndConds[0][n]->GetBoundaryConditionType() == SpatialDomains::eDirichlet)&& (boost::iequals(BndConds[0][n]->GetUserDefined(),"ZeroNormalComponent")))
480  {
481  for(i = 0; i < BndExp[0][n]->GetExpSize(); ++i,cnt++)
482  {
483  elmtid = m_fieldsBCToElmtID[fieldid][cnt];
484  elmt = m_fields[0]->GetExp(elmtid);
485  boundary = m_fieldsBCToTraceID[fieldid][cnt];
486 
487  normals = elmt->GetSurfaceNormal(boundary);
488 
489  nq = BndExp[0][n]->GetExp(i)->GetTotPoints();
490  Array<OneD, NekDouble> normvel(nq,0.0);
491 
492  for(int k = 0; k < m_spacedim; ++k)
493  {
494  Bphys = BndExp[k][n]->UpdatePhys()+
495  BndExp[k][n]->GetPhys_Offset(i);
496  Bc = BndExp[k][n]->GetExp(i);
497  Vmath::Vvtvp(nq,normals[k],1,Bphys,1,normvel,1,
498  normvel,1);
499  }
500 
501  // negate normvel for next step
502  Vmath::Neg(nq,normvel,1);
503 
504  for(int k = 0; k < m_spacedim; ++k)
505  {
506  Bphys = BndExp[k][n]->UpdatePhys()+
507  BndExp[k][n]->GetPhys_Offset(i);
508  Bcoeffs = BndExp[k][n]->UpdateCoeffs()+
509  BndExp[k][n]->GetCoeff_Offset(i);
510  Bc = BndExp[k][n]->GetExp(i);
511  Vmath::Vvtvp(nq,normvel,1,normals[k],1,Bphys,1,
512  Bphys,1);
513  Bc->FwdTrans_BndConstrained(Bphys,Bcoeffs);
514  }
515  }
516  }
517  else
518  {
519  cnt += BndExp[0][n]->GetExpSize();
520  }
521  }
522  }
Array< OneD, Array< OneD, int > > m_fieldsBCToTraceID
Mapping from BCs to Elmt Edge IDs.
Array< OneD, int > m_velocity
int which identifies which components of m_fields contains the velocity (u,v,w);
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:428
int m_spacedim
Spatial dimension (>= expansion dim).
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:382
Array< OneD, Array< OneD, int > > m_fieldsBCToElmtID
Mapping from BCs to Elmt IDs.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
boost::shared_ptr< StdExpansion > StdExpansionSharedPtr
void Nektar::IncNavierStokes::v_GetFluxVector ( const int  i,
Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, NekDouble > > &  flux 
)
virtual

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 253 of file IncNavierStokes.cpp.

References ASSERTL1, Nektar::SolverUtils::EquationSystem::GetNpoints(), Nektar::SolverUtils::EquationSystem::m_fields, m_velocity, and Vmath::Vmul().

256  {
257  ASSERTL1(flux.num_elements() == m_velocity.num_elements(),"Dimension of flux array and velocity array do not match");
258 
259  for(int j = 0; j < flux.num_elements(); ++j)
260  {
261  Vmath::Vmul(GetNpoints(), physfield[i], 1, m_fields[m_velocity[j]]->GetPhys(), 1, flux[j], 1);
262  }
263  }
Array< OneD, int > m_velocity
int which identifies which components of m_fields contains the velocity (u,v,w);
SOLVER_UTILS_EXPORT int GetNpoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218
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:169
virtual int Nektar::IncNavierStokes::v_GetForceDimension ( )
protectedpure virtual
virtual MultiRegions::ExpListSharedPtr Nektar::IncNavierStokes::v_GetPressure ( void  )
inlineprotectedvirtual

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 214 of file IncNavierStokes.h.

References m_pressure.

215  {
216  return m_pressure;
217  }
MultiRegions::ExpListSharedPtr m_pressure
Pointer to field holding pressure field.
void Nektar::IncNavierStokes::v_InitObject ( )
virtual

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::SolverUtils::AdvectionSystem.

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

Definition at line 66 of file IncNavierStokes.cpp.

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, >::CreateInstance(), Nektar::SpatialDomains::eDirichlet, Nektar::eEquationTypeSize, Nektar::eNoEquationType, Nektar::SpatialDomains::eRobin, Nektar::eSteadyLinearisedNS, Nektar::eSteadyNavierStokes, Nektar::eSteadyOseen, Nektar::eSteadyStokes, Nektar::eUnsteadyLinearisedNS, Nektar::eUnsteadyNavierStokes, Nektar::eUnsteadyStokes, Nektar::LibUtilities::Equation::Evaluate(), Nektar::SolverUtils::GetAdvectionFactory(), Nektar::kEquationTypeStr, Nektar::SolverUtils::Forcing::Load(), Nektar::SolverUtils::AdvectionSystem::m_advObject, Nektar::SolverUtils::EquationSystem::m_boundaryConditions, m_cflsteps, m_equationType, Nektar::SolverUtils::EquationSystem::m_fieldMetaDataMap, Nektar::SolverUtils::EquationSystem::m_fields, m_fieldsBCToElmtID, m_fieldsBCToTraceID, m_fieldsRadiationFactor, m_forcing, Nektar::SolverUtils::UnsteadySystem::m_infosteps, m_kinvis, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_spacedim, m_steadyStateSteps, m_steadyStateTol, Nektar::SolverUtils::EquationSystem::m_time, Nektar::SolverUtils::EquationSystem::m_timestep, m_velocity, v_GetForceDimension(), and Nektar::SolverUtils::AdvectionSystem::v_InitObject().

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

67  {
69 
70  int i,j;
71  int numfields = m_fields.num_elements();
72  std::string velids[] = {"u","v","w"};
73 
74  // Set up Velocity field to point to the first m_expdim of m_fields;
75  m_velocity = Array<OneD,int>(m_spacedim);
76 
77  for(i = 0; i < m_spacedim; ++i)
78  {
79  for(j = 0; j < numfields; ++j)
80  {
81  std::string var = m_boundaryConditions->GetVariable(j);
82  if(boost::iequals(velids[i], var))
83  {
84  m_velocity[i] = j;
85  break;
86  }
87 
88  ASSERTL0(j != numfields, "Failed to find field: " + var);
89  }
90  }
91 
92  // Set up equation type enum using kEquationTypeStr
93  for(i = 0; i < (int) eEquationTypeSize; ++i)
94  {
95  bool match;
96  m_session->MatchSolverInfo("EQTYPE",kEquationTypeStr[i],match,false);
97  if(match)
98  {
100  break;
101  }
102  }
103  ASSERTL0(i != eEquationTypeSize,"EQTYPE not found in SOLVERINFO section");
104 
105  // This probably should to into specific implementations
106  // Equation specific Setups
107  switch(m_equationType)
108  {
109  case eSteadyStokes:
110  case eSteadyOseen:
111  case eSteadyNavierStokes:
112  case eSteadyLinearisedNS:
113  break;
115  case eUnsteadyStokes:
116  {
117  m_session->LoadParameter("IO_InfoSteps", m_infosteps, 0);
118  m_session->LoadParameter("IO_CFLSteps", m_cflsteps, 0);
119  m_session->LoadParameter("SteadyStateSteps", m_steadyStateSteps, 0);
120  m_session->LoadParameter("SteadyStateTol", m_steadyStateTol, 1e-6);
121  }
122  break;
123  case eNoEquationType:
124  default:
125  ASSERTL0(false,"Unknown or undefined equation type");
126  }
127 
128  m_session->LoadParameter("Kinvis", m_kinvis);
129 
130  // Default advection type per solver
131  std::string vConvectiveType;
132  switch(m_equationType)
133  {
134  case eUnsteadyStokes:
135  case eSteadyLinearisedNS:
136  vConvectiveType = "NoAdvection";
137  break;
139  case eSteadyNavierStokes:
140  vConvectiveType = "Convective";
141  break;
143  vConvectiveType = "Linearised";
144  break;
145  default:
146  break;
147  }
148 
149  // Check if advection type overridden
150  if (m_session->DefinesTag("AdvectiveType") && m_equationType != eUnsteadyStokes &&
152  {
153  vConvectiveType = m_session->GetTag("AdvectiveType");
154  }
155 
156  // Initialise advection
157  m_advObject = SolverUtils::GetAdvectionFactory().CreateInstance(vConvectiveType, vConvectiveType);
158  m_advObject->InitObject( m_session, m_fields);
159 
160  // Forcing terms
163 
164  // check to see if any Robin boundary conditions and if so set
165  // up m_field to boundary condition maps;
166  m_fieldsBCToElmtID = Array<OneD, Array<OneD, int> >(numfields);
167  m_fieldsBCToTraceID = Array<OneD, Array<OneD, int> >(numfields);
168  m_fieldsRadiationFactor = Array<OneD, Array<OneD, NekDouble> > (numfields);
169 
170  for (i = 0; i < m_fields.num_elements(); ++i)
171  {
172  bool Set = false;
173 
174  Array<OneD, const SpatialDomains::BoundaryConditionShPtr > BndConds;
175  Array<OneD, MultiRegions::ExpListSharedPtr> BndExp;
176  int radpts = 0;
177 
178  BndConds = m_fields[i]->GetBndConditions();
179  BndExp = m_fields[i]->GetBndCondExpansions();
180  for(int n = 0; n < BndConds.num_elements(); ++n)
181  {
182  if(boost::iequals(BndConds[n]->GetUserDefined(),"Radiation"))
183  {
184  ASSERTL0(BndConds[n]->GetBoundaryConditionType() == SpatialDomains::eRobin,
185  "Radiation boundary condition must be of type Robin <R>");
186 
187  if(Set == false)
188  {
189  m_fields[i]->GetBoundaryToElmtMap(m_fieldsBCToElmtID[i],m_fieldsBCToTraceID[i]);
190  Set = true;
191  }
192  radpts += BndExp[n]->GetTotPoints();
193  }
194  if(boost::iequals(BndConds[n]->GetUserDefined(),"ZeroNormalComponent"))
195  {
196  ASSERTL0(BndConds[n]->GetBoundaryConditionType() == SpatialDomains::eDirichlet,
197  "Zero Normal Component boundary condition option must be of type Dirichlet <D>");
198 
199  if(Set == false)
200  {
201  m_fields[i]->GetBoundaryToElmtMap(m_fieldsBCToElmtID[i],m_fieldsBCToTraceID[i]);
202  Set = true;
203  }
204  }
205  }
206 
207  m_fieldsRadiationFactor[i] = Array<OneD, NekDouble>(radpts);
208 
209  radpts = 0; // reset to use as a counter
210 
211  for(int n = 0; n < BndConds.num_elements(); ++n)
212  {
213  if(boost::iequals(BndConds[n]->GetUserDefined(),"Radiation"))
214  {
215 
216  int npoints = BndExp[n]->GetNpoints();
217  Array<OneD, NekDouble> x0(npoints,0.0);
218  Array<OneD, NekDouble> x1(npoints,0.0);
219  Array<OneD, NekDouble> x2(npoints,0.0);
220  Array<OneD, NekDouble> tmpArray;
221 
222  BndExp[n]->GetCoords(x0,x1,x2);
223 
224  LibUtilities::Equation coeff =
225  boost::static_pointer_cast<
226  SpatialDomains::RobinBoundaryCondition
227  >(BndConds[n])->m_robinPrimitiveCoeff;
228 
229  coeff.Evaluate(x0,x1,x2,m_time,
230  tmpArray = m_fieldsRadiationFactor[i]+ radpts);
231  //Vmath::Neg(npoints,tmpArray = m_fieldsRadiationFactor[i]+ radpts,1);
232  radpts += npoints;
233  }
234  }
235  }
236 
237  // Set up Field Meta Data for output files
238  m_fieldMetaDataMap["Kinvis"] = boost::lexical_cast<std::string>(m_kinvis);
239  m_fieldMetaDataMap["TimeStep"] = boost::lexical_cast<std::string>(m_timestep);
240  }
EquationType m_equationType
equation type;
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
Array< OneD, Array< OneD, int > > m_fieldsBCToTraceID
Mapping from BCs to Elmt Edge IDs.
tBaseSharedPtr CreateInstance(tKey idKey BOOST_PP_COMMA_IF(MAX_PARAM) BOOST_PP_ENUM_BINARY_PARAMS(MAX_PARAM, tParam, x))
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:162
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
NekDouble m_time
Current time of simulation.
NekDouble m_kinvis
Kinematic viscosity.
NekDouble m_timestep
Time step size.
Array< OneD, int > m_velocity
int which identifies which components of m_fields contains the velocity (u,v,w);
static SOLVER_UTILS_EXPORT std::vector< ForcingSharedPtr > Load(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields=0)
Definition: Forcing.cpp:86
Array< OneD, Array< OneD, NekDouble > > m_fieldsRadiationFactor
RHS Factor for Radiation Condition.
int m_cflsteps
dump cfl estimate
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Forcing terms.
int m_steadyStateSteps
Check for steady state at step interval.
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
Pointer to boundary conditions object.
int m_spacedim
Spatial dimension (>= expansion dim).
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
Map to identify relevant solver info to dump in output fields.
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:46
Array< OneD, Array< OneD, int > > m_fieldsBCToElmtID
Mapping from BCs to Elmt IDs.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
const std::string kEquationTypeStr[]
int m_infosteps
Number of time steps between outputting status information.
NekDouble m_steadyStateTol
Tolerance to which steady state should be evaluated at.
virtual SOLVER_UTILS_EXPORT void v_InitObject()
Init object for UnsteadySystem class.
virtual int v_GetForceDimension()=0
void Nektar::IncNavierStokes::v_NumericalFlux ( Array< OneD, Array< OneD, NekDouble > > &  physfield,
Array< OneD, Array< OneD, NekDouble > > &  numflux 
)
virtual

Calcualate numerical fluxes

Counter variable

Number of trace points

Number of spatial dimensions

Forward state array

Backward state array

Normal velocity array

Compute the numerical fluxes at the trace points

Extract forwards/backwards trace spaces

Upwind between elements

Calculate the numerical fluxes multipling Fwd or Bwd by the normal advection velocity

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 268 of file IncNavierStokes.cpp.

References Nektar::SolverUtils::EquationSystem::GetTraceNpoints(), Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_spacedim, Nektar::SolverUtils::EquationSystem::m_traceNormals, m_velocity, Vmath::Vmul(), and Vmath::Vvtvp().

270  {
271  /// Counter variable
272  int i;
273 
274  /// Number of trace points
275  int nTracePts = GetTraceNpoints();
276 
277  /// Number of spatial dimensions
278  int nDimensions = m_spacedim;
279 
280  /// Forward state array
281  Array<OneD, NekDouble> Fwd(2*nTracePts);
282 
283  /// Backward state array
284  Array<OneD, NekDouble> Bwd = Fwd + nTracePts;
285 
286  /// Normal velocity array
287  Array<OneD, NekDouble> Vn (nTracePts, 0.0);
288 
289  // Extract velocity field along the trace space and multiply by trace normals
290  for(i = 0; i < nDimensions; ++i)
291  {
292  m_fields[0]->ExtractTracePhys(m_fields[m_velocity[i]]->GetPhys(), Fwd);
293  Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, Fwd, 1, Vn, 1, Vn, 1);
294  }
295 
296  /// Compute the numerical fluxes at the trace points
297  for(i = 0; i < numflux.num_elements(); ++i)
298  {
299  /// Extract forwards/backwards trace spaces
300  m_fields[i]->GetFwdBwdTracePhys(physfield[i], Fwd, Bwd);
301 
302  /// Upwind between elements
303  m_fields[i]->GetTrace()->Upwind(Vn, Fwd, Bwd, numflux[i]);
304 
305  /// Calculate the numerical fluxes multipling Fwd or Bwd
306  /// by the normal advection velocity
307  Vmath::Vmul(nTracePts, numflux[i], 1, Vn, 1, numflux[i], 1);
308  }
309  }
Array< OneD, int > m_velocity
int which identifies which components of m_fields contains the velocity (u,v,w);
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:428
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
int m_spacedim
Spatial dimension (>= expansion dim).
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetTraceNpoints()
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:169
bool Nektar::IncNavierStokes::v_PostIntegrate ( int  step)
protectedvirtual

Estimate CFL and perform steady-state check

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 658 of file IncNavierStokes.cpp.

References CalcSteadyState(), GetCFLEstimate(), m_cflsteps, Nektar::SolverUtils::EquationSystem::m_comm, m_steadyStateSteps, and m_steadyStateTol.

659  {
660  if(m_cflsteps && !((step+1)%m_cflsteps))
661  {
662  int elmtid;
663  NekDouble cfl = GetCFLEstimate(elmtid);
664 
665  if(m_comm->GetRank() == 0)
666  {
667  cout << "CFL (zero plane): "<< cfl << " (in elmt "
668  << elmtid << ")" << endl;
669  }
670  }
671 
672  if(m_steadyStateSteps && step && (!((step+1)%m_steadyStateSteps)))
673  {
674  if(CalcSteadyState() == true)
675  {
676  cout << "Reached Steady State to tolerance "
677  << m_steadyStateTol << endl;
678  return true;
679  }
680  }
681 
682  return false;
683  }
bool CalcSteadyState(void)
evaluate steady state
LibUtilities::CommSharedPtr m_comm
Communicator.
int m_cflsteps
dump cfl estimate
int m_steadyStateSteps
Check for steady state at step interval.
double NekDouble
NekDouble m_steadyStateTol
Tolerance to which steady state should be evaluated at.
NekDouble GetCFLEstimate(int &elmtid)
bool Nektar::IncNavierStokes::v_PreIntegrate ( int  step)
protectedvirtual

Perform the extrapolation.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 647 of file IncNavierStokes.cpp.

References m_extrapolation, Nektar::SolverUtils::UnsteadySystem::m_intSoln, and Nektar::SolverUtils::EquationSystem::m_time.

648  {
649  m_extrapolation->SubStepSaveFields(step);
650  m_extrapolation->SubStepAdvance(m_intSoln,step,m_time);
651  return false;
652  }
NekDouble m_time
Current time of simulation.
ExtrapolateSharedPtr m_extrapolation
LibUtilities::TimeIntegrationSolutionSharedPtr m_intSoln
virtual void Nektar::IncNavierStokes::v_TransCoeffToPhys ( void  )
inlineprotectedvirtual

Virtual function for transformation to physical space.

Reimplemented from Nektar::SolverUtils::EquationSystem.

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

Definition at line 219 of file IncNavierStokes.h.

References ASSERTL0.

220  {
221  ASSERTL0(false,"This method is not defined in this class");
222  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
virtual void Nektar::IncNavierStokes::v_TransPhysToCoeff ( void  )
inlineprotectedvirtual

Virtual function for transformation to coefficient space.

Reimplemented from Nektar::SolverUtils::EquationSystem.

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

Definition at line 224 of file IncNavierStokes.h.

References ASSERTL0.

225  {
226  ASSERTL0(false,"This method is not defined in this class");
227  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
void Nektar::IncNavierStokes::WriteModalEnergy ( void  )
protected

Member Data Documentation

int Nektar::IncNavierStokes::m_cflsteps
protected

dump cfl estimate

Definition at line 168 of file IncNavierStokes.h.

Referenced by v_InitObject(), and v_PostIntegrate().

int Nektar::IncNavierStokes::m_energysteps
protected

dump energy to file at steps time

Definition at line 166 of file IncNavierStokes.h.

EquationType Nektar::IncNavierStokes::m_equationType
protected
ExtrapolateSharedPtr Nektar::IncNavierStokes::m_extrapolation
protected
Array<OneD, Array<OneD, int> > Nektar::IncNavierStokes::m_fieldsBCToElmtID
protected

Mapping from BCs to Elmt IDs.

Definition at line 178 of file IncNavierStokes.h.

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

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

Mapping from BCs to Elmt Edge IDs.

Definition at line 180 of file IncNavierStokes.h.

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

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

RHS Factor for Radiation Condition.

Definition at line 182 of file IncNavierStokes.h.

Referenced by SetRadiationBoundaryForcing(), and v_InitObject().

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

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

Definition at line 186 of file IncNavierStokes.h.

NekDouble Nektar::IncNavierStokes::m_kinvis
protected
std::ofstream Nektar::IncNavierStokes::m_mdlFile
protected

modal energy file

Definition at line 146 of file IncNavierStokes.h.

int Nektar::IncNavierStokes::m_nConvectiveFields
protected
MultiRegions::ExpListSharedPtr Nektar::IncNavierStokes::m_pressure
protected
bool Nektar::IncNavierStokes::m_SmoothAdvection
protected
int Nektar::IncNavierStokes::m_steadyStateSteps
protected

Check for steady state at step interval.

Definition at line 170 of file IncNavierStokes.h.

Referenced by v_InitObject(), and v_PostIntegrate().

NekDouble Nektar::IncNavierStokes::m_steadyStateTol
protected

Tolerance to which steady state should be evaluated at.

Definition at line 172 of file IncNavierStokes.h.

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

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