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 Reference

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.
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)
AdvectionTermSharedPtr GetAdvObject (void)
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::UnsteadySystem
virtual SOLVER_UTILS_EXPORT ~UnsteadySystem ()
 Destructor.
SOLVER_UTILS_EXPORT NekDouble GetTimeStep (const Array< OneD, const Array< OneD, NekDouble > > &inarray)
 Calculate the larger time-step mantaining the problem stable.
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor.
SOLVER_UTILS_EXPORT void SetUpTraceNormals (void)
SOLVER_UTILS_EXPORT void InitObject ()
 Initialises the members of this object.
SOLVER_UTILS_EXPORT void DoInitialise ()
 Perform any initialisation necessary before solving the problem.
SOLVER_UTILS_EXPORT void DoSolve ()
 Solve the problem.
SOLVER_UTILS_EXPORT void TransCoeffToPhys ()
 Transform from coefficient to physical space.
SOLVER_UTILS_EXPORT void TransPhysToCoeff ()
 Transform from physical to coefficient space.
SOLVER_UTILS_EXPORT void Output ()
 Perform output operations after solve.
SOLVER_UTILS_EXPORT NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation.
SOLVER_UTILS_EXPORT std::string GetSessionName ()
 Get Session name.
SOLVER_UTILS_EXPORT void ResetSessionName (std::string newname)
 Reset Session name.
SOLVER_UTILS_EXPORT
LibUtilities::SessionReaderSharedPtr 
GetSession ()
 Get Session name.
SOLVER_UTILS_EXPORT
MultiRegions::ExpListSharedPtr 
GetPressure ()
 Get pressure field if available.
SOLVER_UTILS_EXPORT void PrintSummary (std::ostream &out)
 Print a summary of parameters and solver characteristics.
SOLVER_UTILS_EXPORT void SetLambda (NekDouble lambda)
 Set parameter m_lambda.
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.
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.
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.
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.
SOLVER_UTILS_EXPORT void InitialiseBaseFlow (Array< OneD, Array< OneD, NekDouble > > &base)
 Perform initialisation of the base flow.
SOLVER_UTILS_EXPORT void SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 Initialise the data in the dependent fields.
SOLVER_UTILS_EXPORT void EvaluateExactSolution (int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 Evaluates an exact solution.
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.
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, bool Normalised=false)
 Compute the L2 error of the fields.
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].
SOLVER_UTILS_EXPORT void WeakAdvectionGreensDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\nabla \phi \cdot F) $.
SOLVER_UTILS_EXPORT void WeakAdvectionDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray)
 Compute the inner product $ (\phi, \nabla \cdot F) $.
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) $.
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.
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.
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.
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n)
 Write checkpoint file of m_fields.
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.
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname)
 Write field data to the given filename.
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.
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 Input field data from the given file.
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.
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.
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.
SOLVER_UTILS_EXPORT void ScanForHistoryPoints ()
 Builds map of which element holds each history point.
SOLVER_UTILS_EXPORT void WriteHistoryData (std::ostream &out)
 Probe each history point and write to file.
SOLVER_UTILS_EXPORT void SessionSummary (SummaryList &vSummary)
 Write out a session summary.
SOLVER_UTILS_EXPORT Array
< OneD,
MultiRegions::ExpListSharedPtr > & 
UpdateFields ()
SOLVER_UTILS_EXPORT
LibUtilities::FieldMetaDataMap
UpdateFieldMetaDataMap ()
 Get hold of FieldInfoMap so it can be updated.
SOLVER_UTILS_EXPORT NekDouble GetFinalTime ()
 Return final time.
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.

Protected Member Functions

 IncNavierStokes (const LibUtilities::SessionReaderSharedPtr &pSession)
 Constructor.
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
void SetRadiationBoundaryForcing (int fieldid)
 Set Radiation forcing term.
bool CalcSteadyState (void)
 evaluate steady state
virtual
MultiRegions::ExpListSharedPtr 
v_GetPressure ()
virtual void v_TransCoeffToPhys (void)
 Virtual function for transformation to physical space.
virtual void v_TransPhysToCoeff (void)
 Virtual function for transformation to coefficient space.
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.
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control.
virtual SOLVER_UTILS_EXPORT void v_DoSolve ()
 Solves an unsteady problem.
virtual SOLVER_UTILS_EXPORT void v_DoInitialise ()
 Sets up initial conditions.
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary (SummaryList &s)
 Print a summary of time stepping parameters.
virtual SOLVER_UTILS_EXPORT void v_AppendOutput1D (Array< OneD, Array< OneD, NekDouble > > &solution1D)
 Print the solution at each solution point in a txt file.
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.
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.
int nocase_cmp (const string &s1, const 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.
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.
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
bool m_subSteppingScheme
 bool to identify if using a substepping scheme
bool m_SmoothAdvection
 bool to identify if advection term smoothing is requested
LibUtilities::TimeIntegrationWrapperSharedPtr m_subStepIntegrationScheme
AdvectionTermSharedPtr m_advObject
 Advection term.
std::vector
< SolverUtils::ForcingSharedPtr
m_forcing
 Forcing terms.
int m_nConvectiveFields
 Number of fields to be convected;.
Array< OneD, int > m_velocity
 int which identifies which components of m_fields contains the velocity (u,v,w);
MultiRegions::ExpListSharedPtr m_pressure
 Pointer to field holding pressure field.
NekDouble m_kinvis
 Kinematic viscosity.
int m_energysteps
 dump energy to file at steps time
int m_cflsteps
 dump cfl estimate
int m_steadyStateSteps
 Check for steady state at step interval.
NekDouble m_steadyStateTol
 Tolerance to which steady state should be evaluated at.
EquationType m_equationType
 equation type;
Array< OneD, Array< OneD, int > > m_fieldsBCToElmtID
 Mapping from BCs to Elmt IDs.
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.
int m_intSteps
 Number of time integration steps AND Order of extrapolation for pressure boundary conditions.
- Protected Attributes inherited from Nektar::SolverUtils::UnsteadySystem
int m_infosteps
 Number of time steps between outputting status information.
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
 Wrapper to the time integration scheme.
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use.
LibUtilities::TimeIntegrationSolutionSharedPtr m_intSoln
NekDouble m_epsilon
bool m_explicitDiffusion
 Indicates if explicit or implicit treatment of diffusion is used.
bool m_explicitAdvection
 Indicates if explicit or implicit treatment of advection is used.
bool m_explicitReaction
 Indicates if explicit or implicit treatment of reaction is used.
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state.
std::vector< int > m_intVariables
std::vector< FilterSharedPtrm_filters
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator.
LibUtilities::SessionReaderSharedPtr m_session
 The session reader.
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output.
Array< OneD,
MultiRegions::ExpListSharedPtr
m_fields
 Array holding all dependent variables.
Array< OneD,
MultiRegions::ExpListSharedPtr
m_base
 Base fields.
Array< OneD,
MultiRegions::ExpListSharedPtr
m_derivedfields
 Array holding all dependent variables.
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
 Pointer to boundary conditions object.
SpatialDomains::MeshGraphSharedPtr m_graph
 Pointer to graph defining mesh.
std::string m_filename
 Filename.
std::string m_sessionName
 Name of the session.
NekDouble m_time
 Current time of simulation.
NekDouble m_fintime
 Finish time of the simulation.
NekDouble m_timestep
 Time step size.
NekDouble m_lambda
 Lambda constant in real system if one required.
NekDouble m_checktime
 Time between checkpoints.
int m_steps
 Number of steps to take.
int m_checksteps
 Number of steps between checkpoints.
int m_spacedim
 Spatial dimension (>= expansion dim).
int m_expdim
 Expansion dimension.
bool m_SingleMode
 Flag to determine if single homogeneous mode is used.
bool m_HalfMode
 Flag to determine if half homogeneous mode is used.
bool m_MultipleModes
 Flag to determine if use multiple homogenenous modes are used.
bool m_useFFT
 Flag to determine if FFT is used for homogeneous transform.
bool m_homogen_dealiasing
 Flag to determine if dealiasing is used for homogeneous simulations.
bool m_specHP_dealiasing
 Flag to determine if dealisising is usde for the Spectral/hp element discretisation.
enum MultiRegions::ProjectionType m_projectionType
 Type of projection; e.g continuous or discontinuous.
Array< OneD, Array< OneD,
NekDouble > > 
m_traceNormals
 Array holding trace normals for DG simulations in the forwards direction.
Array< OneD, Array< OneD,
Array< OneD, NekDouble > > > 
m_gradtan
 1 x nvariable x nq
Array< OneD, Array< OneD,
Array< OneD, NekDouble > > > 
m_tanbasis
 2 x m_spacedim x nq
Array< OneD, bool > m_checkIfSystemSingular
 Flag to indicate if the fields should be checked for singularity.
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
 Map to identify relevant solver info to dump in output fields.
int m_NumQuadPointsError
 Number of Quadrature points used to work out the error.
enum HomogeneousType m_HomogeneousType
NekDouble m_LhomX
 physical length in X direction (if homogeneous)
NekDouble m_LhomY
 physical length in Y direction (if homogeneous)
NekDouble m_LhomZ
 physical length in Z direction (if homogeneous)
int m_npointsX
 number of points in X direction (if homogeneous)
int m_npointsY
 number of points in Y direction (if homogeneous)
int m_npointsZ
 number of points in Z direction (if homogeneous)
int m_HomoDirec
 number of homogenous directions
int m_NumMode
 Mode to use in case of single mode analysis.

Additional Inherited Members

- Public Attributes inherited from Nektar::SolverUtils::UnsteadySystem
NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1).
- 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 271 of file IncNavierStokes.cpp.

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

Constructor.

Constructor. Creates ...

Parameters
\param

Definition at line 56 of file IncNavierStokes.cpp.

Member Function Documentation

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

Add an additional forcing term programmatically.

Definition at line 469 of file IncNavierStokes.cpp.

References m_forcing.

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

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

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

Referenced by v_PostIntegrate().

{
static NekDouble previousL2 = 0.0;
bool returnval = false;
NekDouble L2 = 0.0;
// calculate L2 discrete summation
int ncoeffs = m_fields[0]->GetNcoeffs();
for(int i = 0; i < m_fields.num_elements(); ++i)
{
L2 += Vmath::Dot(ncoeffs,m_fields[i]->GetCoeffs(),1,m_fields[i]->GetCoeffs(),1);
}
if(fabs(L2-previousL2) < ncoeffs*m_steadyStateTol)
{
returnval = true;
}
previousL2 = L2;
return returnval;
}
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 340 of file IncNavierStokes.cpp.

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

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

{
int i;
int nqtot = m_fields[0]->GetTotPoints();
int VelDim = m_velocity.num_elements();
Array<OneD, Array<OneD, NekDouble> > velocity(VelDim);
Array<OneD, NekDouble > Deriv;
for(i = 0; i < VelDim; ++i)
{
velocity[i] = inarray[m_velocity[i]];
}
// Set up Derivative work space;
if(wk.num_elements())
{
ASSERTL0(wk.num_elements() >= nqtot*VelDim,
"Workspace is not sufficient");
Deriv = wk;
}
else
{
Deriv = Array<OneD, NekDouble> (nqtot*VelDim);
}
m_velocity,inarray,outarray,m_time,Deriv);
}
AdvectionTermSharedPtr Nektar::IncNavierStokes::GetAdvObject ( void  )
inline

Definition at line 122 of file IncNavierStokes.h.

References m_advObject.

{
return m_advObject;
}
NekDouble Nektar::IncNavierStokes::GetCFLEstimate ( int &  elmtid)

Definition at line 555 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().

{
int n_element = m_fields[0]->GetExpSize();
Array<OneD, NekDouble> cfl = GetElmtCFLVals();
elmtid = Vmath::Imax(n_element,cfl,1);
NekDouble CFL,CFL_loc;
CFL = CFL_loc = cfl[elmtid];
// unshuffle elmt id if data is not stored in consecutive order.
elmtid = m_fields[0]->GetExp(elmtid)->GetGeom()->GetGlobalID();
if(CFL != CFL_loc)
{
elmtid = -1;
}
m_comm->AllReduce(elmtid,LibUtilities::ReduceMax);
// express element id with respect to plane
{
elmtid = elmtid%m_fields[0]->GetPlane(0)->GetExpSize();
}
return CFL;
}
Array< OneD, NekDouble > Nektar::IncNavierStokes::GetElmtCFLVals ( void  )

Definition at line 508 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().

{
int n_vel = m_velocity.num_elements();
int n_element = m_fields[0]->GetExpSize();
const Array<OneD, int> ExpOrder = GetNumExpModesPerExp();
Array<OneD, int> ExpOrderList (n_element, ExpOrder);
const NekDouble cLambda = 0.2; // Spencer book pag. 317
Array<OneD, NekDouble> cfl (n_element, 0.0);
Array<OneD, NekDouble> stdVelocity(n_element, 0.0);
Array<OneD, Array<OneD, NekDouble> > velfields;
if(m_HomogeneousType == eHomogeneous1D) // just do check on 2D info
{
velfields = Array<OneD, Array<OneD, NekDouble> >(2);
for(int i = 0; i < 2; ++i)
{
velfields[i] = m_fields[m_velocity[i]]->UpdatePhys();
}
}
else
{
velfields = Array<OneD, Array<OneD, NekDouble> >(n_vel);
for(int i = 0; i < n_vel; ++i)
{
velfields[i] = m_fields[m_velocity[i]]->UpdatePhys();
}
}
stdVelocity = m_extrapolation->GetMaxStdVelocity(velfields);
for(int el = 0; el < n_element; ++el)
{
cfl[el] = m_timestep*(stdVelocity[el] * cLambda *
(ExpOrder[el]-1) * (ExpOrder[el]-1));
}
return cfl;
}
EquationType Nektar::IncNavierStokes::GetEquationType ( void  )
inlineprotected

Definition at line 204 of file IncNavierStokes.h.

References m_equationType.

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

Definition at line 128 of file IncNavierStokes.h.

References m_nConvectiveFields.

{
}
Array<OneD, int>& Nektar::IncNavierStokes::GetVelocity ( void  )
inline

Definition at line 133 of file IncNavierStokes.h.

References m_velocity.

{
return m_velocity;
}
void Nektar::IncNavierStokes::SetBoundaryConditions ( NekDouble  time)
protected

time dependent boundary conditions updating

Time dependent boundary conditions updating

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 372 of file IncNavierStokes.cpp.

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

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

{
int i, n;
std::string varName;
int nvariables = m_fields.num_elements();
for (i = 0; i < nvariables; ++i)
{
for(n = 0; n < m_fields[i]->GetBndConditions().num_elements(); ++n)
{
if(m_fields[i]->GetBndConditions()[n]->GetUserDefined() ==
{
varName = m_session->GetVariable(i);
m_fields[i]->EvaluateBoundaryConditions(time, varName);
}
}
// Set Radiation conditions if required
}
}
void Nektar::IncNavierStokes::SetRadiationBoundaryForcing ( int  fieldid)
protected

Set Radiation forcing term.

Probably should be pushed back into ContField?

Definition at line 399 of file IncNavierStokes.cpp.

References ASSERTL0, Nektar::SpatialDomains::eHigh, Nektar::SpatialDomains::eHighOutflow, Nektar::SpatialDomains::eNoUserDefined, Nektar::SpatialDomains::eRadiation, Nektar::SpatialDomains::eRobin, Nektar::SpatialDomains::eTimeDependent, Nektar::SpatialDomains::eWall_Forces, Nektar::SolverUtils::EquationSystem::GetPhys_Offset(), Nektar::SolverUtils::EquationSystem::m_fields, m_fieldsBCToElmtID, m_fieldsBCToTraceID, m_fieldsRadiationFactor, and Vmath::Vmul().

Referenced by SetBoundaryConditions().

{
int i,n;
Array<OneD, const SpatialDomains::BoundaryConditionShPtr > BndConds;
Array<OneD, MultiRegions::ExpListSharedPtr> BndExp;
BndConds = m_fields[fieldid]->GetBndConditions();
BndExp = m_fields[fieldid]->GetBndCondExpansions();
int cnt;
int elmtid,nq,offset, boundary;
Array<OneD, NekDouble> Bvals, U;
int cnt1 = 0;
for(cnt = n = 0; n < BndConds.num_elements(); ++n)
{
SpatialDomains::BndUserDefinedType type = BndConds[n]->GetUserDefined();
if((BndConds[n]->GetBoundaryConditionType() == SpatialDomains::eRobin)&&(type == SpatialDomains::eRadiation))
{
for(i = 0; i < BndExp[n]->GetExpSize(); ++i,cnt++)
{
elmtid = m_fieldsBCToElmtID[fieldid][cnt];
elmt = m_fields[fieldid]->GetExp(elmtid);
offset = m_fields[fieldid]->GetPhys_Offset(elmtid);
U = m_fields[fieldid]->UpdatePhys() + offset;
Bc = BndExp[n]->GetExp(i);
boundary = m_fieldsBCToTraceID[fieldid][cnt];
// Get edge values and put into ubc
nq = Bc->GetTotPoints();
Array<OneD, NekDouble> ubc(nq);
elmt->GetTracePhysVals(boundary,Bc,U,ubc);
Vmath::Vmul(nq,&m_fieldsRadiationFactor[fieldid][cnt1 +
BndExp[n]->GetPhys_Offset(i)],1,&ubc[0],1,&ubc[0],1);
Bvals = BndExp[n]->UpdateCoeffs()+BndExp[n]->GetCoeff_Offset(i);
Bc->IProductWRTBase(ubc,Bvals);
}
cnt1 += BndExp[n]->GetTotPoints();
}
else if(type == SpatialDomains::eNoUserDefined ||
{
cnt += BndExp[n]->GetExpSize();
}
else
{
ASSERTL0(false,"Unknown USERDEFINEDTYPE in pressure boundary condition");
}
}
}
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 279 of file IncNavierStokes.cpp.

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

{
ASSERTL1(flux.num_elements() == m_velocity.num_elements(),"Dimension of flux array and velocity array do not match");
for(int j = 0; j < flux.num_elements(); ++j)
{
Vmath::Vmul(GetNpoints(), physfield[i], 1, m_fields[m_velocity[j]]->GetPhys(), 1, flux[j], 1);
}
}
virtual int Nektar::IncNavierStokes::v_GetForceDimension ( )
protectedpure virtual

Referenced by v_InitObject().

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

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 224 of file IncNavierStokes.h.

References m_pressure.

{
return m_pressure;
}
void Nektar::IncNavierStokes::v_InitObject ( )
virtual

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

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

Definition at line 64 of file IncNavierStokes.cpp.

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, >::CreateInstance(), Nektar::eEquationTypeSize, Nektar::SpatialDomains::eHighOutflow, Nektar::SpatialDomains::eI, Nektar::eNoEquationType, Nektar::SpatialDomains::eNoUserDefined, Nektar::SpatialDomains::eRadiation, Nektar::SpatialDomains::eRobin, Nektar::eSteadyLinearisedNS, Nektar::eSteadyNavierStokes, Nektar::eSteadyOseen, Nektar::eSteadyStokes, Nektar::SpatialDomains::eTimeDependent, Nektar::eUnsteadyLinearisedNS, Nektar::eUnsteadyNavierStokes, Nektar::eUnsteadyStokes, Nektar::LibUtilities::Equation::Evaluate(), Nektar::SpatialDomains::eWall_Forces, Nektar::GetAdvectionTermFactory(), Nektar::GetExtrapolateFactory(), Nektar::kEquationTypeStr, Nektar::SolverUtils::Forcing::Load(), m_advObject, Nektar::SolverUtils::EquationSystem::m_boundaryConditions, m_cflsteps, m_equationType, m_extrapolation, Nektar::SolverUtils::EquationSystem::m_fieldMetaDataMap, Nektar::SolverUtils::EquationSystem::m_fields, m_fieldsBCToElmtID, m_fieldsBCToTraceID, m_fieldsRadiationFactor, m_forcing, Nektar::SolverUtils::EquationSystem::m_graph, Nektar::SolverUtils::UnsteadySystem::m_infosteps, m_kinvis, m_pressure, 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().

{
int i,j;
int numfields = m_fields.num_elements();
std::string velids[] = {"u","v","w"};
// Set up Velocity field to point to the first m_expdim of m_fields;
m_velocity = Array<OneD,int>(m_spacedim);
for(i = 0; i < m_spacedim; ++i)
{
for(j = 0; j < numfields; ++j)
{
std::string var = m_boundaryConditions->GetVariable(j);
if(boost::iequals(velids[i], var))
{
m_velocity[i] = j;
break;
}
ASSERTL0(j != numfields, "Failed to find field: " + var);
}
}
// Set up equation type enum using kEquationTypeStr
for(i = 0; i < (int) eEquationTypeSize; ++i)
{
bool match;
m_session->MatchSolverInfo("EQTYPE",kEquationTypeStr[i],match,false);
if(match)
{
break;
}
}
ASSERTL0(i != eEquationTypeSize,"EQTYPE not found in SOLVERINFO section");
// This probably should to into specific implementations
// Equation specific Setups
{
case eSteadyOseen:
break;
{
m_session->LoadParameter("IO_InfoSteps", m_infosteps, 0);
m_session->LoadParameter("IO_CFLSteps", m_cflsteps, 0);
m_session->LoadParameter("SteadyStateSteps", m_steadyStateSteps, 0);
m_session->LoadParameter("SteadyStateTol", m_steadyStateTol, 1e-6);
// check to see if any user defined boundary condition is
// indeed implemented
for(int n = 0; n < m_fields[0]->GetBndConditions().num_elements(); ++n)
{
// Time Dependent Boundary Condition (if no user
// defined then this is empty)
m_fields[0]->GetBndConditions()[n]->GetUserDefined() ==
m_fields[0]->GetBndConditions()[n]->GetUserDefined() ==
m_fields[0]->GetBndConditions()[n]->GetUserDefined() ==
m_fields[0]->GetBndConditions()[n]->GetUserDefined() ==
m_fields[0]->GetBndConditions()[n]->GetUserDefined() ==
m_fields[0]->GetBndConditions()[n]->GetUserDefined() ==
"Unknown USERDEFINEDTYPE boundary condition");
}
}
break;
default:
ASSERTL0(false,"Unknown or undefined equation type");
}
m_session->LoadParameter("Kinvis", m_kinvis);
// Default advection type per solver
std::string vConvectiveType;
{
vConvectiveType = "NoAdvection";
break;
vConvectiveType = "Convective";
break;
vConvectiveType = "Linearised";
break;
default:
break;
}
// Check if advection type overridden
if (m_session->DefinesTag("AdvectiveType") && m_equationType != eUnsteadyStokes)
{
vConvectiveType = m_session->GetTag("AdvectiveType");
}
// Initialise advection
// Forcing terms
// check to see if any Robin boundary conditions and if so set
// up m_field to boundary condition maps;
m_fieldsBCToElmtID = Array<OneD, Array<OneD, int> >(numfields);
m_fieldsBCToTraceID = Array<OneD, Array<OneD, int> >(numfields);
m_fieldsRadiationFactor = Array<OneD, Array<OneD, NekDouble> > (numfields);
for (i = 0; i < m_fields.num_elements(); ++i)
{
bool Set = false;
Array<OneD, const SpatialDomains::BoundaryConditionShPtr > BndConds;
Array<OneD, MultiRegions::ExpListSharedPtr> BndExp;
int radpts = 0;
BndConds = m_fields[i]->GetBndConditions();
BndExp = m_fields[i]->GetBndCondExpansions();
for(int n = 0; n < BndConds.num_elements(); ++n)
{
if(BndConds[n]->GetUserDefined() == SpatialDomains::eRadiation)
{
ASSERTL0(BndConds[n]->GetBoundaryConditionType() == SpatialDomains::eRobin,
"Radiation boundary condition must be of type Robin <R>");
if(Set == false)
{
m_fields[i]->GetBoundaryToElmtMap(m_fieldsBCToElmtID[i],m_fieldsBCToTraceID[i]);
Set = true;
}
radpts += BndExp[n]->GetTotPoints();
}
}
m_fieldsRadiationFactor[i] = Array<OneD, NekDouble>(radpts);
radpts = 0; // reset to use as a counter
for(int n = 0; n < BndConds.num_elements(); ++n)
{
if(BndConds[n]->GetUserDefined() == SpatialDomains::eRadiation)
{
int npoints = BndExp[n]->GetNpoints();
Array<OneD, NekDouble> x0(npoints,0.0);
Array<OneD, NekDouble> x1(npoints,0.0);
Array<OneD, NekDouble> x2(npoints,0.0);
Array<OneD, NekDouble> tmpArray;
BndExp[n]->GetCoords(x0,x1,x2);
boost::static_pointer_cast<
>(BndConds[n])->m_robinPrimitiveCoeff;
coeff.Evaluate(x0,x1,x2,m_time,
tmpArray = m_fieldsRadiationFactor[i]+ radpts);
//Vmath::Neg(npoints,tmpArray = m_fieldsRadiationFactor[i]+ radpts,1);
radpts += npoints;
}
}
}
// Set up Field Meta Data for output files
m_fieldMetaDataMap["Kinvis"] = boost::lexical_cast<std::string>(m_kinvis);
m_fieldMetaDataMap["TimeStep"] = boost::lexical_cast<std::string>(m_timestep);
// creation of the extrapolation object
{
std::string vExtrapolation = "Standard";
if (m_session->DefinesSolverInfo("Extrapolation"))
{
vExtrapolation = m_session->GetSolverInfo("Extrapolation");
}
vExtrapolation,
}
}
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 294 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().

{
/// Counter variable
int i;
/// Number of trace points
int nTracePts = GetTraceNpoints();
/// Number of spatial dimensions
int nDimensions = m_spacedim;
/// Forward state array
Array<OneD, NekDouble> Fwd(2*nTracePts);
/// Backward state array
Array<OneD, NekDouble> Bwd = Fwd + nTracePts;
/// Normal velocity array
Array<OneD, NekDouble> Vn (nTracePts, 0.0);
// Extract velocity field along the trace space and multiply by trace normals
for(i = 0; i < nDimensions; ++i)
{
m_fields[0]->ExtractTracePhys(m_fields[m_velocity[i]]->GetPhys(), Fwd);
Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, Fwd, 1, Vn, 1, Vn, 1);
}
/// Compute the numerical fluxes at the trace points
for(i = 0; i < numflux.num_elements(); ++i)
{
/// Extract forwards/backwards trace spaces
m_fields[i]->GetFwdBwdTracePhys(physfield[i], Fwd, Bwd);
/// Upwind between elements
m_fields[i]->GetTrace()->Upwind(Vn, Fwd, Bwd, numflux[i]);
/// Calculate the numerical fluxes multipling Fwd or Bwd
/// by the normal advection velocity
Vmath::Vmul(nTracePts, numflux[i], 1, Vn, 1, numflux[i], 1);
}
}
bool Nektar::IncNavierStokes::v_PostIntegrate ( int  step)
protectedvirtual

Estimate CFL and perform steady-state check

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 599 of file IncNavierStokes.cpp.

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

{
if(m_cflsteps && !((step+1)%m_cflsteps))
{
int elmtid;
NekDouble cfl = GetCFLEstimate(elmtid);
if(m_comm->GetRank() == 0)
{
cout << "CFL (zero plane): "<< cfl << " (in elmt "
<< elmtid << ")" << endl;
}
}
if(m_steadyStateSteps && step && (!((step+1)%m_steadyStateSteps)))
{
if(CalcSteadyState() == true)
{
cout << "Reached Steady State to tolerance "
<< m_steadyStateTol << endl;
return true;
}
}
return false;
}
bool Nektar::IncNavierStokes::v_PreIntegrate ( int  step)
protectedvirtual

Perform the extrapolation.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 588 of file IncNavierStokes.cpp.

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

{
m_extrapolation->SubStepSaveFields(step);
m_extrapolation->SubStepAdvance(m_intSoln,step,m_time);
return false;
}
virtual void Nektar::IncNavierStokes::v_TransCoeffToPhys ( void  )
inlineprotectedvirtual

Virtual function for transformation to physical space.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 229 of file IncNavierStokes.h.

References ASSERTL0.

{
ASSERTL0(false,"This method is not defined in this class");
}
virtual void Nektar::IncNavierStokes::v_TransPhysToCoeff ( void  )
inlineprotectedvirtual

Virtual function for transformation to coefficient space.

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 234 of file IncNavierStokes.h.

References ASSERTL0.

{
ASSERTL0(false,"This method is not defined in this class");
}
void Nektar::IncNavierStokes::WriteModalEnergy ( void  )
protected

Member Data Documentation

AdvectionTermSharedPtr Nektar::IncNavierStokes::m_advObject
protected

Advection term.

Definition at line 162 of file IncNavierStokes.h.

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

int Nektar::IncNavierStokes::m_cflsteps
protected

dump cfl estimate

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

EquationType Nektar::IncNavierStokes::m_equationType
protected

equation type;

Definition at line 188 of file IncNavierStokes.h.

Referenced by GetEquationType(), Nektar::CoupledLinearNS::SolveLinearNS(), Nektar::CoupledLinearNS::v_DoInitialise(), Nektar::CoupledLinearNS::v_DoSolve(), and v_InitObject().

ExtrapolateSharedPtr Nektar::IncNavierStokes::m_extrapolation
protected

Definition at line 149 of file IncNavierStokes.h.

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

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

Mapping from BCs to Elmt IDs.

Definition at line 191 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 193 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 195 of file IncNavierStokes.h.

Referenced by SetRadiationBoundaryForcing(), and v_InitObject().

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

Forcing terms.

Definition at line 165 of file IncNavierStokes.h.

Referenced by AddForcing(), Nektar::CoupledLinearNS::EvaluateAdvection(), Nektar::VelocityCorrectionScheme::EvaluateAdvection_SetPressureBCs(), Nektar::CoupledLinearNS::Solve(), and v_InitObject().

int Nektar::IncNavierStokes::m_intSteps
protected

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

Definition at line 199 of file IncNavierStokes.h.

NekDouble Nektar::IncNavierStokes::m_kinvis
protected

Kinematic viscosity.

Definition at line 177 of file IncNavierStokes.h.

Referenced by Nektar::CoupledLinearNS::Continuation(), Nektar::VelocityCorrectionScheme::EvaluateAdvection_SetPressureBCs(), Nektar::CoupledLinearNS::EvaluateNewtonRHS(), Nektar::CoupledLinearNS::SetUpCoupledMatrix(), Nektar::VelocityCorrectionScheme::SetUpViscousForcing(), Nektar::CoupledLinearNS::SolveSteadyNavierStokes(), Nektar::VelocityCorrectionScheme::SolveUnsteadyStokesSystem(), Nektar::VelocityCorrectionScheme::v_DoInitialise(), Nektar::CoupledLinearNS::v_DoInitialise(), Nektar::CoupledLinearNS::v_DoSolve(), Nektar::VelocityCorrectionScheme::v_InitObject(), and v_InitObject().

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

modal energy file

Definition at line 152 of file IncNavierStokes.h.

int Nektar::IncNavierStokes::m_nConvectiveFields
protected

Number of fields to be convected;.

Definition at line 168 of file IncNavierStokes.h.

Referenced by Nektar::VelocityCorrectionScheme::EvaluateAdvection_SetPressureBCs(), EvaluateAdvectionTerms(), GetNConvectiveFields(), Nektar::VelocityCorrectionScheme::SolveUnsteadyStokesSystem(), Nektar::CoupledLinearNS::SolveUnsteadyStokesSystem(), Nektar::VelocityCorrectionScheme::v_DoInitialise(), Nektar::VelocityCorrectionScheme::v_InitObject(), and Nektar::CoupledLinearNS::v_InitObject().

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

Pointer to field holding pressure field.

Definition at line 175 of file IncNavierStokes.h.

Referenced by Nektar::VelocityCorrectionScheme::EvaluateAdvection_SetPressureBCs(), Nektar::CoupledLinearNS::SetUpCoupledMatrix(), Nektar::VelocityCorrectionScheme::SetUpViscousForcing(), Nektar::CoupledLinearNS::SolveLinearNS(), Nektar::VelocityCorrectionScheme::SolveUnsteadyStokesSystem(), v_GetPressure(), Nektar::VelocityCorrectionScheme::v_InitObject(), v_InitObject(), Nektar::CoupledLinearNS::v_InitObject(), and Nektar::CoupledLinearNS::v_Output().

bool Nektar::IncNavierStokes::m_SmoothAdvection
protected

bool to identify if advection term smoothing is requested

Definition at line 157 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 183 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 185 of file IncNavierStokes.h.

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

LibUtilities::TimeIntegrationWrapperSharedPtr Nektar::IncNavierStokes::m_subStepIntegrationScheme
protected

Definition at line 159 of file IncNavierStokes.h.

Referenced by Nektar::VelocityCorrectionScheme::v_GenerateSummary().

bool Nektar::IncNavierStokes::m_subSteppingScheme
protected

bool to identify if using a substepping scheme

Definition at line 155 of file IncNavierStokes.h.

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

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

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

Definition at line 172 of file IncNavierStokes.h.

Referenced by Nektar::CoupledLinearNS::Continuation(), Nektar::CoupledLinearNS::DefineForcingTerm(), Nektar::VelocityCorrectionScheme::EvaluateAdvection_SetPressureBCs(), EvaluateAdvectionTerms(), Nektar::CoupledLinearNS::EvaluateNewtonRHS(), GetElmtCFLVals(), GetVelocity(), Nektar::CoupledLinearNS::InfNorm(), Nektar::CoupledLinearNS::L2Norm(), Nektar::CoupledLinearNS::SetUpCoupledMatrix(), Nektar::VelocityCorrectionScheme::SetUpPressureForcing(), Nektar::VelocityCorrectionScheme::SetUpViscousForcing(), Nektar::CoupledLinearNS::Solve(), Nektar::CoupledLinearNS::SolveLinearNS(), Nektar::CoupledLinearNS::SolveSteadyNavierStokes(), Nektar::CoupledLinearNS::SolveUnsteadyStokesSystem(), Nektar::CoupledLinearNS::v_DoInitialise(), v_GetFluxVector(), Nektar::VelocityCorrectionScheme::v_InitObject(), v_InitObject(), Nektar::CoupledLinearNS::v_InitObject(), and v_NumericalFlux().