Nektar++
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 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 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, 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 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::FieldMetaDataMapUpdateFieldMetaDataMap ()
 Get hold of FieldInfoMap so it can be updated. More...
 
SOLVER_UTILS_EXPORT NekDouble GetFinalTime ()
 Return final time. More...
 
SOLVER_UTILS_EXPORT int GetNcoeffs ()
 
SOLVER_UTILS_EXPORT int GetNcoeffs (const int eid)
 
SOLVER_UTILS_EXPORT int GetNumExpModes ()
 
SOLVER_UTILS_EXPORT const Array< OneD, int > GetNumExpModesPerExp ()
 
SOLVER_UTILS_EXPORT int GetNvariables ()
 
SOLVER_UTILS_EXPORT const std::string GetVariable (unsigned int i)
 
SOLVER_UTILS_EXPORT int GetTraceTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTraceNpoints ()
 
SOLVER_UTILS_EXPORT int GetExpSize ()
 
SOLVER_UTILS_EXPORT int GetPhys_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetCoeff_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTotPoints (int n)
 
SOLVER_UTILS_EXPORT int GetNpoints ()
 
SOLVER_UTILS_EXPORT int 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 SetStepsToOne ()
 
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 string &s1, const string &s2)
 Perform a case-insensitive string comparison. 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...
 
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)
 
- 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 string &s1, const string &s2)
 
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Virtual function for the L_inf error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false)
 Virtual function for the L_2 error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 
virtual SOLVER_UTILS_EXPORT void v_EvaluateExactSolution (unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 
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_subSteppingScheme
 bool to identify if using a substepping scheme More...
 
bool m_SmoothAdvection
 bool to identify if advection term smoothing is requested More...
 
LibUtilities::TimeIntegrationWrapperSharedPtr m_subStepIntegrationScheme
 
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...
 
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...
 
map< std::string, Array< OneD, Array< OneD, float > > > m_interpWeights
 Map of the interpolation weights for a specific filename. More...
 
map< std::string, Array< OneD, Array< OneD, unsigned int > > > m_interpInds
 Map of the interpolation indices for a specific filename. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 Array holding all dependent variables. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_base
 Base fields. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_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...
 
NekDouble m_fintime
 Finish time of the simulation. More...
 
NekDouble m_timestep
 Time step size. More...
 
NekDouble m_lambda
 Lambda constant in real system if one required. More...
 
NekDouble m_checktime
 Time between checkpoints. More...
 
int m_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...
 
int m_NumMode
 Mode to use in case of single mode analysis. 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 251 of file IncNavierStokes.cpp.

252  {
253  }
Nektar::IncNavierStokes::IncNavierStokes ( const LibUtilities::SessionReaderSharedPtr pSession)
protected

Constructor.

Constructor. Creates ...

Parameters

Definition at line 56 of file IncNavierStokes.cpp.

56  :
57  UnsteadySystem(pSession),
58  AdvectionSystem(pSession),
59  m_subSteppingScheme(false),
60  m_SmoothAdvection(false),
62  {
63  }
bool m_subSteppingScheme
bool to identify if using a substepping scheme
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 450 of file IncNavierStokes.cpp.

References m_forcing.

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

451  {
452  m_forcing.push_back(pForce);
453  }
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 461 of file IncNavierStokes.cpp.

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

Referenced by v_PostIntegrate().

462  {
463  static NekDouble previousL2 = 0.0;
464  bool returnval = false;
465 
466  NekDouble L2 = 0.0;
467 
468  // calculate L2 discrete summation
469  int ncoeffs = m_fields[0]->GetNcoeffs();
470 
471  for(int i = 0; i < m_fields.num_elements(); ++i)
472  {
473  L2 += Vmath::Dot(ncoeffs,m_fields[i]->GetCoeffs(),1,m_fields[i]->GetCoeffs(),1);
474  }
475 
476  if(fabs(L2-previousL2) < ncoeffs*m_steadyStateTol)
477  {
478  returnval = true;
479  }
480 
481  previousL2 = L2;
482 
483  return returnval;
484  }
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:891
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 320 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::VelocityCorrectionScheme::EvaluateAdvection_SetPressureBCs(), and Nektar::CoupledLinearNS::EvaluateNewtonRHS().

323  {
324  int i;
325  int nqtot = m_fields[0]->GetTotPoints();
326  int VelDim = m_velocity.num_elements();
327  Array<OneD, Array<OneD, NekDouble> > velocity(VelDim);
329 
330  for(i = 0; i < VelDim; ++i)
331  {
332  if(m_fields[i]->GetWaveSpace() && !m_singleMode && !m_halfMode)
333  {
334  velocity[i] = Array<OneD, NekDouble>(nqtot,0.0);
335  m_fields[i]->HomogeneousBwdTrans(inarray[m_velocity[i]],velocity[i]);
336  }
337  else
338  {
339  velocity[i] = inarray[m_velocity[i]];
340  }
341  }
342 
343  // Set up Derivative work space;
344  if(wk.num_elements())
345  {
346  ASSERTL0(wk.num_elements() >= nqtot*VelDim,
347  "Workspace is not sufficient");
348  Deriv = wk;
349  }
350  else
351  {
352  Deriv = Array<OneD, NekDouble> (nqtot*VelDim);
353  }
354 
356  velocity, inarray, outarray, m_time);
357  }
bool m_singleMode
Flag to determine if single homogeneous mode is used.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
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 536 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().

537  {
538  int n_element = m_fields[0]->GetExpSize();
539 
541 
542  elmtid = Vmath::Imax(n_element,cfl,1);
543  NekDouble CFL,CFL_loc;
544 
545  CFL = CFL_loc = cfl[elmtid];
546  m_comm->AllReduce(CFL,LibUtilities::ReduceMax);
547 
548  // unshuffle elmt id if data is not stored in consecutive order.
549  elmtid = m_fields[0]->GetExp(elmtid)->GetGeom()->GetGlobalID();
550  if(CFL != CFL_loc)
551  {
552  elmtid = -1;
553  }
554 
555  m_comm->AllReduce(elmtid,LibUtilities::ReduceMax);
556 
557  // express element id with respect to plane
559  {
560  elmtid = elmtid%m_fields[0]->GetPlane(0)->GetExpSize();
561  }
562  return CFL;
563  }
int Imax(int n, const T *x, const int incx)
Return the index of the maximum element in x.
Definition: Vmath.cpp:732
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 489 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().

490  {
491  int n_vel = m_velocity.num_elements();
492  int n_element = m_fields[0]->GetExpSize();
493 
494  const Array<OneD, int> ExpOrder = GetNumExpModesPerExp();
495  Array<OneD, int> ExpOrderList (n_element, ExpOrder);
496 
497  const NekDouble cLambda = 0.2; // Spencer book pag. 317
498 
499  Array<OneD, NekDouble> cfl (n_element, 0.0);
500  Array<OneD, NekDouble> stdVelocity(n_element, 0.0);
502 
503  if(m_HomogeneousType == eHomogeneous1D) // just do check on 2D info
504  {
505  velfields = Array<OneD, Array<OneD, NekDouble> >(2);
506 
507  for(int i = 0; i < 2; ++i)
508  {
509  velfields[i] = m_fields[m_velocity[i]]->UpdatePhys();
510  }
511  }
512  else
513  {
514  velfields = Array<OneD, Array<OneD, NekDouble> >(n_vel);
515 
516  for(int i = 0; i < n_vel; ++i)
517  {
518  velfields[i] = m_fields[m_velocity[i]]->UpdatePhys();
519  }
520  }
521 
522  stdVelocity = m_extrapolation->GetMaxStdVelocity(velfields);
523 
524  for(int el = 0; el < n_element; ++el)
525  {
526  cfl[el] = m_timestep*(stdVelocity[el] * cLambda *
527  (ExpOrder[el]-1) * (ExpOrder[el]-1));
528  }
529 
530  return cfl;
531  }
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 195 of file IncNavierStokes.h.

References m_equationType.

196  {
197  return m_equationType;
198  }
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 362 of file IncNavierStokes.cpp.

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

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

363  {
364  int i, n;
365  std::string varName;
366  int nvariables = m_fields.num_elements();
367 
368  for (i = 0; i < nvariables; ++i)
369  {
370  for(n = 0; n < m_fields[i]->GetBndConditions().num_elements(); ++n)
371  {
372  if(m_fields[i]->GetBndConditions()[n]->IsTimeDependent() ||
373  m_fields[i]->GetBndConditions()[n]->GetUserDefined() ==
374  "MovingBody")
375  {
376  varName = m_session->GetVariable(i);
377  m_fields[i]->EvaluateBoundaryConditions(time, varName);
378  }
379 
380  }
381 
382  // Set Radiation conditions if required
384  }
385  }
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 390 of file IncNavierStokes.cpp.

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

Referenced by SetBoundaryConditions().

391  {
392  int i,n;
393 
396 
397  BndConds = m_fields[fieldid]->GetBndConditions();
398  BndExp = m_fields[fieldid]->GetBndCondExpansions();
399 
402 
403  int cnt;
404  int elmtid,nq,offset, boundary;
405  Array<OneD, NekDouble> Bvals, U;
406  int cnt1 = 0;
407 
408  for(cnt = n = 0; n < BndConds.num_elements(); ++n)
409  {
410  std::string type = BndConds[n]->GetUserDefined();
411 
412  if((BndConds[n]->GetBoundaryConditionType() == SpatialDomains::eRobin)&&(boost::iequals(type,"Radiation")))
413  {
414  for(i = 0; i < BndExp[n]->GetExpSize(); ++i,cnt++)
415  {
416  elmtid = m_fieldsBCToElmtID[fieldid][cnt];
417  elmt = m_fields[fieldid]->GetExp(elmtid);
418  offset = m_fields[fieldid]->GetPhys_Offset(elmtid);
419 
420  U = m_fields[fieldid]->UpdatePhys() + offset;
421  Bc = BndExp[n]->GetExp(i);
422 
423  boundary = m_fieldsBCToTraceID[fieldid][cnt];
424 
425  // Get edge values and put into ubc
426  nq = Bc->GetTotPoints();
427  Array<OneD, NekDouble> ubc(nq);
428  elmt->GetTracePhysVals(boundary,Bc,U,ubc);
429 
430  Vmath::Vmul(nq,&m_fieldsRadiationFactor[fieldid][cnt1 +
431  BndExp[n]->GetPhys_Offset(i)],1,&ubc[0],1,&ubc[0],1);
432 
433  Bvals = BndExp[n]->UpdateCoeffs()+BndExp[n]->GetCoeff_Offset(i);
434 
435  Bc->IProductWRTBase(ubc,Bvals);
436  }
437  cnt1 += BndExp[n]->GetTotPoints();
438  }
439  else
440  {
441  cnt += BndExp[n]->GetExpSize();
442  }
443  }
444  }
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.
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::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 259 of file IncNavierStokes.cpp.

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

262  {
263  ASSERTL1(flux.num_elements() == m_velocity.num_elements(),"Dimension of flux array and velocity array do not match");
264 
265  for(int j = 0; j < flux.num_elements(); ++j)
266  {
267  Vmath::Vmul(GetNpoints(), physfield[i], 1, m_fields[m_velocity[j]]->GetPhys(), 1, flux[j], 1);
268  }
269  }
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:191
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 215 of file IncNavierStokes.h.

References m_pressure.

216  {
217  return m_pressure;
218  }
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, and Nektar::VelocityCorrectionScheme.

Definition at line 65 of file IncNavierStokes.cpp.

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, >::CreateInstance(), 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, and v_GetForceDimension().

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

66  {
67  AdvectionSystem::v_InitObject();
68 
69  int i,j;
70  int numfields = m_fields.num_elements();
71  std::string velids[] = {"u","v","w"};
72 
73  // Set up Velocity field to point to the first m_expdim of m_fields;
75 
76  for(i = 0; i < m_spacedim; ++i)
77  {
78  for(j = 0; j < numfields; ++j)
79  {
80  std::string var = m_boundaryConditions->GetVariable(j);
81  if(boost::iequals(velids[i], var))
82  {
83  m_velocity[i] = j;
84  break;
85  }
86 
87  ASSERTL0(j != numfields, "Failed to find field: " + var);
88  }
89  }
90 
91  // Set up equation type enum using kEquationTypeStr
92  for(i = 0; i < (int) eEquationTypeSize; ++i)
93  {
94  bool match;
95  m_session->MatchSolverInfo("EQTYPE",kEquationTypeStr[i],match,false);
96  if(match)
97  {
99  break;
100  }
101  }
102  ASSERTL0(i != eEquationTypeSize,"EQTYPE not found in SOLVERINFO section");
103 
104  // This probably should to into specific implementations
105  // Equation specific Setups
106  switch(m_equationType)
107  {
108  case eSteadyStokes:
109  case eSteadyOseen:
110  case eSteadyNavierStokes:
111  case eSteadyLinearisedNS:
112  break;
114  case eUnsteadyStokes:
115  {
116  m_session->LoadParameter("IO_InfoSteps", m_infosteps, 0);
117  m_session->LoadParameter("IO_CFLSteps", m_cflsteps, 0);
118  m_session->LoadParameter("SteadyStateSteps", m_steadyStateSteps, 0);
119  m_session->LoadParameter("SteadyStateTol", m_steadyStateTol, 1e-6);
120 
121  // check to see if any user defined boundary condition is
122  // indeed implemented
123 
124  for(int n = 0; n < m_fields[0]->GetBndConditions().num_elements(); ++n)
125  {
126  std::string type =m_fields[0]->GetBndConditions()[n]->GetUserDefined();
127  if(!type.empty())
128  // Time Dependent Boundary Condition (if no user
129  // defined then this is empty)
130  ASSERTL0 (boost::iequals(type,"Wall_Forces") ||
131  boost::iequals(type,"TimeDependent") ||
132  boost::iequals(type,"MovingBody") ||
133  boost::iequals(type,"Radiation") ||
134  boost::iequals(type,"I") ||
135  boost::iequals(type,"HOutflow"),
136  "Unknown USERDEFINEDTYPE boundary condition");
137  }
138  }
139  break;
140  case eNoEquationType:
141  default:
142  ASSERTL0(false,"Unknown or undefined equation type");
143  }
144 
145  m_session->LoadParameter("Kinvis", m_kinvis);
146 
147  // Default advection type per solver
148  std::string vConvectiveType;
149  switch(m_equationType)
150  {
151  case eUnsteadyStokes:
152  case eSteadyLinearisedNS:
153  vConvectiveType = "NoAdvection";
154  break;
156  case eSteadyNavierStokes:
157  vConvectiveType = "Convective";
158  break;
160  vConvectiveType = "Linearised";
161  break;
162  default:
163  break;
164  }
165 
166  // Check if advection type overridden
167  if (m_session->DefinesTag("AdvectiveType") && m_equationType != eUnsteadyStokes &&
169  {
170  vConvectiveType = m_session->GetTag("AdvectiveType");
171  }
172 
173  // Initialise advection
174  m_advObject = SolverUtils::GetAdvectionFactory().CreateInstance(vConvectiveType, vConvectiveType);
175  m_advObject->InitObject( m_session, m_fields);
176 
177  // Forcing terms
180 
181  // check to see if any Robin boundary conditions and if so set
182  // up m_field to boundary condition maps;
186 
187  for (i = 0; i < m_fields.num_elements(); ++i)
188  {
189  bool Set = false;
190 
193  int radpts = 0;
194 
195  BndConds = m_fields[i]->GetBndConditions();
196  BndExp = m_fields[i]->GetBndCondExpansions();
197  for(int n = 0; n < BndConds.num_elements(); ++n)
198  {
199  if(boost::iequals(BndConds[n]->GetUserDefined(),"Radiation"))
200  {
201  ASSERTL0(BndConds[n]->GetBoundaryConditionType() == SpatialDomains::eRobin,
202  "Radiation boundary condition must be of type Robin <R>");
203 
204  if(Set == false)
205  {
206  m_fields[i]->GetBoundaryToElmtMap(m_fieldsBCToElmtID[i],m_fieldsBCToTraceID[i]);
207  Set = true;
208  }
209  radpts += BndExp[n]->GetTotPoints();
210  }
211  }
212 
213  m_fieldsRadiationFactor[i] = Array<OneD, NekDouble>(radpts);
214 
215  radpts = 0; // reset to use as a counter
216 
217  for(int n = 0; n < BndConds.num_elements(); ++n)
218  {
219  if(boost::iequals(BndConds[n]->GetUserDefined(),"Radiation"))
220  {
221 
222  int npoints = BndExp[n]->GetNpoints();
223  Array<OneD, NekDouble> x0(npoints,0.0);
224  Array<OneD, NekDouble> x1(npoints,0.0);
225  Array<OneD, NekDouble> x2(npoints,0.0);
226  Array<OneD, NekDouble> tmpArray;
227 
228  BndExp[n]->GetCoords(x0,x1,x2);
229 
230  LibUtilities::Equation coeff =
231  boost::static_pointer_cast<
233  >(BndConds[n])->m_robinPrimitiveCoeff;
234 
235  coeff.Evaluate(x0,x1,x2,m_time,
236  tmpArray = m_fieldsRadiationFactor[i]+ radpts);
237  //Vmath::Neg(npoints,tmpArray = m_fieldsRadiationFactor[i]+ radpts,1);
238  radpts += npoints;
239  }
240  }
241  }
242 
243  // Set up Field Meta Data for output files
244  m_fieldMetaDataMap["Kinvis"] = boost::lexical_cast<std::string>(m_kinvis);
245  m_fieldMetaDataMap["TimeStep"] = boost::lexical_cast<std::string>(m_timestep);
246  }
EquationType m_equationType
equation type;
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
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:84
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.
NekDouble Evaluate() const
Definition: Equation.h:102
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 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 274 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().

276  {
277  /// Counter variable
278  int i;
279 
280  /// Number of trace points
281  int nTracePts = GetTraceNpoints();
282 
283  /// Number of spatial dimensions
284  int nDimensions = m_spacedim;
285 
286  /// Forward state array
287  Array<OneD, NekDouble> Fwd(2*nTracePts);
288 
289  /// Backward state array
290  Array<OneD, NekDouble> Bwd = Fwd + nTracePts;
291 
292  /// Normal velocity array
293  Array<OneD, NekDouble> Vn (nTracePts, 0.0);
294 
295  // Extract velocity field along the trace space and multiply by trace normals
296  for(i = 0; i < nDimensions; ++i)
297  {
298  m_fields[0]->ExtractTracePhys(m_fields[m_velocity[i]]->GetPhys(), Fwd);
299  Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, Fwd, 1, Vn, 1, Vn, 1);
300  }
301 
302  /// Compute the numerical fluxes at the trace points
303  for(i = 0; i < numflux.num_elements(); ++i)
304  {
305  /// Extract forwards/backwards trace spaces
306  m_fields[i]->GetFwdBwdTracePhys(physfield[i], Fwd, Bwd);
307 
308  /// Upwind between elements
309  m_fields[i]->GetTrace()->Upwind(Vn, Fwd, Bwd, numflux[i]);
310 
311  /// Calculate the numerical fluxes multipling Fwd or Bwd
312  /// by the normal advection velocity
313  Vmath::Vmul(nTracePts, numflux[i], 1, Vn, 1, numflux[i], 1);
314  }
315  }
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 580 of file IncNavierStokes.cpp.

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

581  {
582  if(m_cflsteps && !((step+1)%m_cflsteps))
583  {
584  int elmtid;
585  NekDouble cfl = GetCFLEstimate(elmtid);
586 
587  if(m_comm->GetRank() == 0)
588  {
589  cout << "CFL (zero plane): "<< cfl << " (in elmt "
590  << elmtid << ")" << endl;
591  }
592  }
593 
594  if(m_steadyStateSteps && step && (!((step+1)%m_steadyStateSteps)))
595  {
596  if(CalcSteadyState() == true)
597  {
598  cout << "Reached Steady State to tolerance "
599  << m_steadyStateTol << endl;
600  return true;
601  }
602  }
603 
604  return false;
605  }
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 569 of file IncNavierStokes.cpp.

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

570  {
571  m_extrapolation->SubStepSaveFields(step);
572  m_extrapolation->SubStepAdvance(m_intSoln,step,m_time);
573  return false;
574  }
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 220 of file IncNavierStokes.h.

References ASSERTL0.

221  {
222  ASSERTL0(false,"This method is not defined in this class");
223  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
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 225 of file IncNavierStokes.h.

References ASSERTL0.

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

Member Data Documentation

int Nektar::IncNavierStokes::m_cflsteps
protected

dump cfl estimate

Definition at line 172 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 170 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 182 of file IncNavierStokes.h.

Referenced by SetRadiationBoundaryForcing(), and v_InitObject().

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

Mapping from BCs to Elmt Edge IDs.

Definition at line 184 of file IncNavierStokes.h.

Referenced by SetRadiationBoundaryForcing(), and v_InitObject().

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

RHS Factor for Radiation Condition.

Definition at line 186 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 190 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

bool to identify if advection term smoothing is requested

Definition at line 151 of file IncNavierStokes.h.

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

int Nektar::IncNavierStokes::m_steadyStateSteps
protected

Check for steady state at step interval.

Definition at line 174 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 176 of file IncNavierStokes.h.

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

LibUtilities::TimeIntegrationWrapperSharedPtr Nektar::IncNavierStokes::m_subStepIntegrationScheme
protected
bool Nektar::IncNavierStokes::m_subSteppingScheme
protected

bool to identify if using a substepping scheme

Definition at line 149 of file IncNavierStokes.h.

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

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