Nektar++
|
This class is the base class for Navier Stokes problems. More...
#include <IncNavierStokes.h>
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) |
int | GetNConvectiveFields (void) |
Array< OneD, int > & | GetVelocity (void) |
Array< OneD, NekDouble > | GetElmtCFLVals (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. | |
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. | |
template<class T > | |
boost::shared_ptr< T > | as () |
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 . | |
SOLVER_UTILS_EXPORT void | WeakAdvectionDivergenceForm (const Array< OneD, Array< OneD, NekDouble > > &F, Array< OneD, NekDouble > &outarray) |
Compute the inner product . | |
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 . | |
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 | Checkpoint_BaseFlow (const int n) |
Write base flow file of m_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 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 |
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::AdvectionSystem | |
SolverUtils::AdvectionSharedPtr | m_advObject |
Advection term. | |
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< FilterSharedPtr > | m_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_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... |
This class is the base class for Navier Stokes problems.
Definition at line 104 of file IncNavierStokes.h.
|
virtual |
|
protected |
Constructor.
Constructor. Creates ...
\param |
Definition at line 56 of file IncNavierStokes.cpp.
void Nektar::IncNavierStokes::AddForcing | ( | const SolverUtils::ForcingSharedPtr & | pForce | ) |
Add an additional forcing term programmatically.
Definition at line 464 of file IncNavierStokes.cpp.
References m_forcing.
Referenced by Nektar::VortexWaveInteraction::ExecuteRoll().
|
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 475 of file IncNavierStokes.cpp.
References Nektar::Dot(), Nektar::SolverUtils::EquationSystem::m_fields, and m_steadyStateTol.
Referenced by v_PostIntegrate().
|
protected |
Evaluation -N(V) for all fields except pressure using m_velocity
Definition at line 325 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().
NekDouble Nektar::IncNavierStokes::GetCFLEstimate | ( | int & | elmtid | ) |
Definition at line 550 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().
Definition at line 503 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().
|
inlineprotected |
Definition at line 195 of file IncNavierStokes.h.
References m_equationType.
|
inline |
Definition at line 122 of file IncNavierStokes.h.
References m_nConvectiveFields.
time dependent boundary conditions updating
Time dependent boundary conditions updating
Reimplemented from Nektar::SolverUtils::EquationSystem.
Definition at line 367 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().
|
protected |
Set Radiation forcing term.
Probably should be pushed back into ContField?
Definition at line 394 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().
|
virtual |
Reimplemented from Nektar::SolverUtils::EquationSystem.
Definition at line 264 of file IncNavierStokes.cpp.
References ASSERTL1, Nektar::SolverUtils::EquationSystem::GetNpoints(), Nektar::SolverUtils::EquationSystem::m_fields, m_velocity, and Vmath::Vmul().
|
protectedpure virtual |
Referenced by v_InitObject().
|
inlineprotectedvirtual |
Reimplemented from Nektar::SolverUtils::EquationSystem.
Definition at line 215 of file IncNavierStokes.h.
References m_pressure.
|
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::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::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().
|
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 279 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().
|
protectedvirtual |
Estimate CFL and perform steady-state check
Reimplemented from Nektar::SolverUtils::UnsteadySystem.
Definition at line 594 of file IncNavierStokes.cpp.
References CalcSteadyState(), GetCFLEstimate(), m_cflsteps, Nektar::SolverUtils::EquationSystem::m_comm, m_steadyStateSteps, and m_steadyStateTol.
|
protectedvirtual |
Perform the extrapolation.
Reimplemented from Nektar::SolverUtils::UnsteadySystem.
Definition at line 583 of file IncNavierStokes.cpp.
References m_extrapolation, Nektar::SolverUtils::UnsteadySystem::m_intSoln, and Nektar::SolverUtils::EquationSystem::m_time.
Virtual function for transformation to physical space.
Reimplemented from Nektar::SolverUtils::EquationSystem.
Definition at line 220 of file IncNavierStokes.h.
References ASSERTL0.
Virtual function for transformation to coefficient space.
Reimplemented from Nektar::SolverUtils::EquationSystem.
Definition at line 225 of file IncNavierStokes.h.
References ASSERTL0.
|
protected |
dump cfl estimate
Definition at line 172 of file IncNavierStokes.h.
Referenced by v_InitObject(), and v_PostIntegrate().
|
protected |
dump energy to file at steps time
Definition at line 170 of file IncNavierStokes.h.
|
protected |
equation type;
Definition at line 179 of file IncNavierStokes.h.
Referenced by GetEquationType(), Nektar::CoupledLinearNS::SolveLinearNS(), Nektar::CoupledLinearNS::v_DoInitialise(), Nektar::CoupledLinearNS::v_DoSolve(), Nektar::VelocityCorrectionScheme::v_InitObject(), v_InitObject(), and Nektar::CoupledLinearNS::v_InitObject().
|
protected |
Definition at line 143 of file IncNavierStokes.h.
Referenced by Nektar::VelocityCorrectionScheme::EvaluateAdvection_SetPressureBCs(), GetElmtCFLVals(), Nektar::VelocityCorrectionScheme::SolveUnsteadyStokesSystem(), Nektar::VelocityCorrectionScheme::v_InitObject(), Nektar::CoupledLinearNS::v_InitObject(), and v_PreIntegrate().
Mapping from BCs to Elmt IDs.
Definition at line 182 of file IncNavierStokes.h.
Referenced by SetRadiationBoundaryForcing(), and v_InitObject().
Mapping from BCs to Elmt Edge IDs.
Definition at line 184 of file IncNavierStokes.h.
Referenced by SetRadiationBoundaryForcing(), and v_InitObject().
RHS Factor for Radiation Condition.
Definition at line 186 of file IncNavierStokes.h.
Referenced by SetRadiationBoundaryForcing(), and v_InitObject().
|
protected |
Forcing terms.
Definition at line 156 of file IncNavierStokes.h.
Referenced by AddForcing(), Nektar::CoupledLinearNS::EvaluateAdvection(), Nektar::VelocityCorrectionScheme::EvaluateAdvection_SetPressureBCs(), Nektar::CoupledLinearNS::Solve(), and v_InitObject().
|
protected |
Number of time integration steps AND Order of extrapolation for pressure boundary conditions.
Definition at line 190 of file IncNavierStokes.h.
|
protected |
Kinematic viscosity.
Definition at line 168 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().
|
protected |
modal energy file
Definition at line 146 of file IncNavierStokes.h.
|
protected |
Number of fields to be convected;.
Definition at line 159 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().
|
protected |
Pointer to field holding pressure field.
Definition at line 166 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(), Nektar::CoupledLinearNS::v_InitObject(), and Nektar::CoupledLinearNS::v_Output().
|
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().
|
protected |
Check for steady state at step interval.
Definition at line 174 of file IncNavierStokes.h.
Referenced by v_InitObject(), and v_PostIntegrate().
|
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().
|
protected |
Definition at line 153 of file IncNavierStokes.h.
Referenced by Nektar::VelocityCorrectionScheme::v_GenerateSummary().
|
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().
|
protected |
int which identifies which components of m_fields contains the velocity (u,v,w);
Definition at line 163 of file IncNavierStokes.h.
Referenced by Nektar::CoupledLinearNS::Continuation(), Nektar::CoupledLinearNS::DefineForcingTerm(), 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().