| 
    Nektar++
    
   | 
 
#include <VelocityCorrectionScheme.h>


Static Public Member Functions | |
| static  SolverUtils::EquationSystemSharedPtr  | create (const LibUtilities::SessionReaderSharedPtr &pSession) | 
| Creates an instance of this class.  More... | |
Static Public Attributes | |
| static std::string | className | 
| Name of class.  More... | |
Protected Member Functions | |
| virtual void | v_GenerateSummary (SolverUtils::SummaryList &s) | 
| Print a summary of time stepping parameters.  More... | |
| 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 void | v_DoInitialise (void) | 
| Sets up initial conditions.  More... | |
| virtual Array< OneD, bool > | v_GetSystemSingularChecks () | 
| virtual int | v_GetForceDimension () | 
| virtual void | v_SetUpPressureForcing (const Array< OneD, const Array< OneD, NekDouble > > &fields, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt) | 
| virtual void | v_SetUpViscousForcing (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt) | 
| virtual void | v_SolvePressure (const Array< OneD, NekDouble > &Forcing) | 
| virtual void | v_SolveViscous (const Array< OneD, const Array< OneD, NekDouble > > &Forcing, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble aii_Dt) | 
| virtual void | v_EvaluateAdvection_SetPressureBCs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time) | 
  Protected Member Functions inherited from Nektar::IncNavierStokes | |
| IncNavierStokes (const LibUtilities::SessionReaderSharedPtr &pSession) | |
| Constructor.  More... | |
| EquationType | GetEquationType (void) | 
| void | EvaluateAdvectionTerms (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, Array< OneD, NekDouble > &wk=NullNekDouble1DArray) | 
| void | WriteModalEnergy (void) | 
| void | SetBoundaryConditions (NekDouble time) | 
| time dependent boundary conditions updating  More... | |
| void | SetRadiationBoundaryForcing (int fieldid) | 
| Set Radiation forcing term.  More... | |
| void | SetZeroNormalVelocity () | 
| Set Normal Velocity Component to Zero.  More... | |
| bool | CalcSteadyState (void) | 
| evaluate steady state  More... | |
| virtual  MultiRegions::ExpListSharedPtr  | v_GetPressure () | 
| virtual 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_AppendOutput1D (Array< OneD, Array< OneD, NekDouble > > &solution1D) | 
| Print the solution at each solution point in a txt file.  More... | |
| virtual SOLVER_UTILS_EXPORT void | v_NumericalFlux (Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &numfluxX, Array< OneD, Array< OneD, NekDouble > > &numfluxY) | 
| virtual SOLVER_UTILS_EXPORT void | v_NumFluxforScalar (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &uflux) | 
| virtual SOLVER_UTILS_EXPORT void | v_NumFluxforVector (const Array< OneD, Array< OneD, NekDouble > > &ufield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux) | 
| virtual SOLVER_UTILS_EXPORT  NekDouble  | v_GetTimeStep (const Array< OneD, const Array< OneD, NekDouble > > &inarray) | 
| Return the timestep to be used for the next step in the time-marching loop.  More... | |
| virtual SOLVER_UTILS_EXPORT bool | v_SteadyStateCheck (int step) | 
| SOLVER_UTILS_EXPORT void | CheckForRestartTime (NekDouble &time) | 
| SOLVER_UTILS_EXPORT void | SVVVarDiffCoeff (const Array< OneD, Array< OneD, NekDouble > > vel, StdRegions::VarCoeffMap &varCoeffMap) | 
| Evaluate the SVV diffusion coefficient according to Moura's paper where it should proportional to h time velocity.  More... | |
  Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem | |
| SOLVER_UTILS_EXPORT | EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession) | 
| Initialises EquationSystem class members.  More... | |
| int | nocase_cmp (const std::string &s1, const std::string &s2) | 
| virtual SOLVER_UTILS_EXPORT  NekDouble  | v_LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray) | 
| Virtual function for the L_inf error computation between fields and a given exact solution.  More... | |
| virtual SOLVER_UTILS_EXPORT  NekDouble  | v_L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false) | 
| Virtual function for the L_2 error computation between fields and a given exact solution.  More... | |
| virtual SOLVER_UTILS_EXPORT void | v_SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0) | 
| virtual SOLVER_UTILS_EXPORT void | v_EvaluateExactSolution (unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time) | 
| SOLVER_UTILS_EXPORT void | SetUpBaseFields (SpatialDomains::MeshGraphSharedPtr &mesh) | 
| SOLVER_UTILS_EXPORT void | ImportFldBase (std::string pInfile, SpatialDomains::MeshGraphSharedPtr pGraph) | 
| virtual SOLVER_UTILS_EXPORT void | v_Output (void) | 
| virtual SOLVER_UTILS_EXPORT void | v_ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables) | 
Protected Attributes | |
| bool | m_useHomo1DSpecVanVisc | 
| bool to identify if spectral vanishing viscosity is active.  More... | |
| bool | m_useSpecVanVisc | 
| bool to identify if spectral vanishing viscosity is active.  More... | |
| NekDouble | m_sVVCutoffRatio | 
| cutt off ratio from which to start decayhing modes  More... | |
| NekDouble | m_sVVDiffCoeff | 
| Diffusion coefficient of SVV modes.  More... | |
| Array< OneD, NekDouble > | m_saved_aii_Dt | 
| Save aiiDt value to use as a trip to reset global matrix setup.  More... | |
| StdRegions::VarCoeffMap | m_varCoeffLap | 
| Variable Coefficient map for the Laplacian which can be activated as part of SVV or otherwise.  More... | |
  Protected Attributes inherited from Nektar::IncNavierStokes | |
| ExtrapolateSharedPtr | m_extrapolation | 
| std::ofstream | m_mdlFile | 
| modal energy file  More... | |
| bool | m_SmoothAdvection | 
| bool to identify if advection term smoothing is requested  More... | |
| std::vector < SolverUtils::ForcingSharedPtr >  | m_forcing | 
| Forcing terms.  More... | |
| int | m_nConvectiveFields | 
| Number of fields to be convected;.  More... | |
| Array< OneD, int > | m_velocity | 
| int which identifies which components of m_fields contains the velocity (u,v,w);  More... | |
| MultiRegions::ExpListSharedPtr | m_pressure | 
| Pointer to field holding pressure field.  More... | |
| NekDouble | m_kinvis | 
| Kinematic viscosity.  More... | |
| int | m_energysteps | 
| dump energy to file at steps time  More... | |
| int | m_cflsteps | 
| dump cfl estimate  More... | |
| int | m_steadyStateSteps | 
| Check for steady state at step interval.  More... | |
| NekDouble | m_steadyStateTol | 
| Tolerance to which steady state should be evaluated at.  More... | |
| EquationType | m_equationType | 
| equation type;  More... | |
| Array< OneD, Array< OneD, int > > | m_fieldsBCToElmtID | 
| Mapping from BCs to Elmt IDs.  More... | |
| Array< OneD, Array< OneD, int > > | m_fieldsBCToTraceID | 
| Mapping from BCs to Elmt Edge IDs.  More... | |
| Array< OneD, Array< OneD,  NekDouble > >  | m_fieldsRadiationFactor | 
| RHS Factor for Radiation Condition.  More... | |
| int | m_intSteps | 
| Number of time integration steps AND Order of extrapolation for pressure boundary conditions.  More... | |
  Protected Attributes inherited from Nektar::SolverUtils::AdvectionSystem | |
| SolverUtils::AdvectionSharedPtr | m_advObject | 
| Advection term.  More... | |
  Protected Attributes inherited from Nektar::SolverUtils::UnsteadySystem | |
| int | m_infosteps | 
| Number of time steps between outputting status information.  More... | |
| LibUtilities::TimeIntegrationWrapperSharedPtr | m_intScheme | 
| Wrapper to the time integration scheme.  More... | |
| LibUtilities::TimeIntegrationSchemeOperators | m_ode | 
| The time integration scheme operators to use.  More... | |
| LibUtilities::TimeIntegrationSolutionSharedPtr | m_intSoln | 
| NekDouble | m_epsilon | 
| bool | m_explicitDiffusion | 
| Indicates if explicit or implicit treatment of diffusion is used.  More... | |
| bool | m_explicitAdvection | 
| Indicates if explicit or implicit treatment of advection is used.  More... | |
| bool | m_explicitReaction | 
| Indicates if explicit or implicit treatment of reaction is used.  More... | |
| bool | m_homoInitialFwd | 
| Flag to determine if simulation should start in homogeneous forward transformed state.  More... | |
| std::vector< int > | m_intVariables | 
| std::vector< FilterSharedPtr > | m_filters | 
  Protected Attributes inherited from Nektar::SolverUtils::EquationSystem | |
| LibUtilities::CommSharedPtr | m_comm | 
| Communicator.  More... | |
| LibUtilities::SessionReaderSharedPtr | m_session | 
| The session reader.  More... | |
| LibUtilities::FieldIOSharedPtr | m_fld | 
| Field input/output.  More... | |
| std::map< std::string, Array < OneD, Array< OneD, float > > >  | m_interpWeights | 
| Map of the interpolation weights for a specific filename.  More... | |
| std::map< std::string, Array < OneD, Array< OneD, unsigned int > > >  | m_interpInds | 
| Map of the interpolation indices for a specific filename.  More... | |
| Array< OneD,  MultiRegions::ExpListSharedPtr >  | m_fields | 
| Array holding all dependent variables.  More... | |
| Array< OneD,  MultiRegions::ExpListSharedPtr >  | m_base | 
| Base fields.  More... | |
| Array< OneD,  MultiRegions::ExpListSharedPtr >  | m_derivedfields | 
| Array holding all dependent variables.  More... | |
| SpatialDomains::BoundaryConditionsSharedPtr | m_boundaryConditions | 
| Pointer to boundary conditions object.  More... | |
| SpatialDomains::MeshGraphSharedPtr | m_graph | 
| Pointer to graph defining mesh.  More... | |
| std::string | m_sessionName | 
| Name of the session.  More... | |
| NekDouble | m_time | 
| Current time of simulation.  More... | |
| int | m_initialStep | 
| Number of the step where the simulation should begin.  More... | |
| NekDouble | m_fintime | 
| Finish time of the simulation.  More... | |
| NekDouble | m_timestep | 
| Time step size.  More... | |
| NekDouble | m_lambda | 
| Lambda constant in real system if one required.  More... | |
| std::set< std::string > | m_loadedFields | 
| NekDouble | m_checktime | 
| Time between checkpoints.  More... | |
| int | m_nchk | 
| Number of checkpoints written so far.  More... | |
| int | m_steps | 
| Number of steps to take.  More... | |
| int | m_checksteps | 
| Number of steps between checkpoints.  More... | |
| int | m_spacedim | 
| Spatial dimension (>= expansion dim).  More... | |
| int | m_expdim | 
| Expansion dimension.  More... | |
| bool | m_singleMode | 
| Flag to determine if single homogeneous mode is used.  More... | |
| bool | m_halfMode | 
| Flag to determine if half homogeneous mode is used.  More... | |
| bool | m_multipleModes | 
| Flag to determine if use multiple homogenenous modes are used.  More... | |
| bool | m_useFFT | 
| Flag to determine if FFT is used for homogeneous transform.  More... | |
| bool | m_homogen_dealiasing | 
| Flag to determine if dealiasing is used for homogeneous simulations.  More... | |
| bool | m_specHP_dealiasing | 
| Flag to determine if dealisising is usde for the Spectral/hp element discretisation.  More... | |
| enum MultiRegions::ProjectionType | m_projectionType | 
| Type of projection; e.g continuous or discontinuous.  More... | |
| Array< OneD, Array< OneD,  NekDouble > >  | m_traceNormals | 
| Array holding trace normals for DG simulations in the forwards direction.  More... | |
| Array< OneD, Array< OneD,  Array< OneD, NekDouble > > >  | m_gradtan | 
| 1 x nvariable x nq  More... | |
| Array< OneD, Array< OneD,  Array< OneD, NekDouble > > >  | m_tanbasis | 
| 2 x m_spacedim x nq  More... | |
| Array< OneD, bool > | m_checkIfSystemSingular | 
| Flag to indicate if the fields should be checked for singularity.  More... | |
| LibUtilities::FieldMetaDataMap | m_fieldMetaDataMap | 
| Map to identify relevant solver info to dump in output fields.  More... | |
| int | m_NumQuadPointsError | 
| Number of Quadrature points used to work out the error.  More... | |
| enum HomogeneousType | m_HomogeneousType | 
| NekDouble | m_LhomX | 
| physical length in X direction (if homogeneous)  More... | |
| NekDouble | m_LhomY | 
| physical length in Y direction (if homogeneous)  More... | |
| NekDouble | m_LhomZ | 
| physical length in Z direction (if homogeneous)  More... | |
| int | m_npointsX | 
| number of points in X direction (if homogeneous)  More... | |
| int | m_npointsY | 
| number of points in Y direction (if homogeneous)  More... | |
| int | m_npointsZ | 
| number of points in Z direction (if homogeneous)  More... | |
| int | m_HomoDirec | 
| number of homogenous directions  More... | |
Additional Inherited Members | |
  Public Attributes inherited from Nektar::SolverUtils::UnsteadySystem | |
| NekDouble | m_cflSafetyFactor | 
| CFL safety factor (comprise between 0 to 1).  More... | |
  Protected Types inherited from Nektar::SolverUtils::EquationSystem | |
| enum | HomogeneousType { eHomogeneous1D, eHomogeneous2D, eHomogeneous3D, eNotHomogeneous } | 
| Parameter for homogeneous expansions.  More... | |
Definition at line 43 of file VelocityCorrectionScheme.h.
| Nektar::VelocityCorrectionScheme::VelocityCorrectionScheme | ( | const LibUtilities::SessionReaderSharedPtr & | pSession | ) | 
Constructor.
Constructor. Creates ...
| param | 
Definition at line 57 of file VelocityCorrectionScheme.cpp.
      
  | 
  virtual | 
      
  | 
  inlinestatic | 
Creates an instance of this class.
Definition at line 48 of file VelocityCorrectionScheme.h.
References Nektar::MemoryManager< DataType >::AllocateSharedPtr().
      
  | 
  inline | 
Definition at line 103 of file VelocityCorrectionScheme.h.
References v_EvaluateAdvection_SetPressureBCs().
Referenced by v_InitObject().
      
  | 
  inline | 
Definition at line 68 of file VelocityCorrectionScheme.h.
References v_SetUpPressureForcing().
Referenced by SolveUnsteadyStokesSystem().
      
  | 
  inline | 
Definition at line 76 of file VelocityCorrectionScheme.h.
References v_SetUpViscousForcing().
Referenced by SolveUnsteadyStokesSystem().
      
  | 
  inline | 
Definition at line 84 of file VelocityCorrectionScheme.h.
References v_SolvePressure().
Referenced by SolveUnsteadyStokesSystem().
| void Nektar::VelocityCorrectionScheme::SolveUnsteadyStokesSystem | ( | const Array< OneD, const Array< OneD, NekDouble > > & | inarray, | 
| Array< OneD, Array< OneD, NekDouble > > & | outarray, | ||
| const NekDouble | time, | ||
| const NekDouble | aii_Dt | ||
| ) | 
Implicit part of the method - Poisson + nConv*Helmholtz
Definition at line 333 of file VelocityCorrectionScheme.cpp.
References Nektar::IncNavierStokes::m_extrapolation, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::IncNavierStokes::m_kinvis, Nektar::IncNavierStokes::m_nConvectiveFields, Nektar::IncNavierStokes::SetBoundaryConditions(), SetUpPressureForcing(), SetUpViscousForcing(), SolvePressure(), and SolveViscous().
Referenced by v_InitObject().
      
  | 
  inline | 
Definition at line 89 of file VelocityCorrectionScheme.h.
References v_SolveViscous().
Referenced by SolveUnsteadyStokesSystem().
Sets up initial conditions.
Sets the initial conditions.
Reimplemented from Nektar::SolverUtils::UnsteadySystem.
Reimplemented in Nektar::VCSMapping.
Definition at line 228 of file VelocityCorrectionScheme.cpp.
References Nektar::SolverUtils::EquationSystem::m_fieldMetaDataMap, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::IncNavierStokes::m_kinvis, Nektar::IncNavierStokes::m_nConvectiveFields, Nektar::SolverUtils::EquationSystem::m_time, Nektar::SolverUtils::EquationSystem::m_timestep, Nektar::IncNavierStokes::SetBoundaryConditions(), and Nektar::SolverUtils::UnsteadySystem::v_DoInitialise().
      
  | 
  protectedvirtual | 
Explicit part of the method - Advection, Forcing + HOPBCs
Reimplemented in Nektar::VCSMapping.
Definition at line 303 of file VelocityCorrectionScheme.cpp.
References Nektar::IncNavierStokes::EvaluateAdvectionTerms(), Nektar::IncNavierStokes::m_extrapolation, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::IncNavierStokes::m_forcing, Nektar::IncNavierStokes::m_kinvis, Nektar::IncNavierStokes::m_nConvectiveFields, Nektar::IncNavierStokes::m_pressure, and Nektar::IncNavierStokes::m_SmoothAdvection.
Referenced by EvaluateAdvection_SetPressureBCs().
      
  | 
  protectedvirtual | 
Print a summary of time stepping parameters.
Prints a summary with some information regards the time-stepping.
Reimplemented from Nektar::SolverUtils::UnsteadySystem.
Definition at line 188 of file VelocityCorrectionScheme.cpp.
References Nektar::SolverUtils::AddSummaryItem(), Nektar::LibUtilities::eNoTimeIntegrationMethod, Nektar::IncNavierStokes::m_extrapolation, Nektar::SolverUtils::EquationSystem::m_homogen_dealiasing, Nektar::SolverUtils::EquationSystem::m_specHP_dealiasing, m_sVVCutoffRatio, m_sVVDiffCoeff, m_useHomo1DSpecVanVisc, m_useSpecVanVisc, Nektar::LibUtilities::TimeIntegrationMethodMap, and Nektar::SolverUtils::UnsteadySystem::v_GenerateSummary().
      
  | 
  protectedvirtual | 
Implements Nektar::IncNavierStokes.
Definition at line 295 of file VelocityCorrectionScheme.cpp.
References Nektar::SolverUtils::EquationSystem::m_session.
      
  | 
  protectedvirtual | 
Reimplemented from Nektar::SolverUtils::EquationSystem.
Definition at line 284 of file VelocityCorrectionScheme.cpp.
References Nektar::SolverUtils::EquationSystem::m_session.
      
  | 
  virtual | 
Init object for UnsteadySystem class.
Initialization object for UnsteadySystem class.
Reimplemented from Nektar::IncNavierStokes.
Reimplemented in Nektar::VCSMapping.
Definition at line 66 of file VelocityCorrectionScheme.cpp.
References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, >::CreateInstance(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineImplicitSolve(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineOdeRhs(), Nektar::SolverUtils::EquationSystem::eHomogeneous1D, Nektar::eUnsteadyNavierStokes, EvaluateAdvection_SetPressureBCs(), Nektar::GetExtrapolateFactory(), Nektar::NekConstants::kNekUnsetDouble, Nektar::SolverUtils::AdvectionSystem::m_advObject, Nektar::IncNavierStokes::m_equationType, Nektar::SolverUtils::UnsteadySystem::m_explicitDiffusion, Nektar::IncNavierStokes::m_extrapolation, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_HomogeneousType, Nektar::SolverUtils::UnsteadySystem::m_intScheme, Nektar::SolverUtils::UnsteadySystem::m_intVariables, Nektar::IncNavierStokes::m_kinvis, Nektar::IncNavierStokes::m_nConvectiveFields, Nektar::SolverUtils::UnsteadySystem::m_ode, Nektar::IncNavierStokes::m_pressure, m_saved_aii_Dt, Nektar::SolverUtils::EquationSystem::m_session, Nektar::IncNavierStokes::m_SmoothAdvection, Nektar::SolverUtils::EquationSystem::m_specHP_dealiasing, m_sVVCutoffRatio, m_sVVDiffCoeff, m_useHomo1DSpecVanVisc, m_useSpecVanVisc, Nektar::IncNavierStokes::m_velocity, SolveUnsteadyStokesSystem(), and Nektar::IncNavierStokes::v_InitObject().
Referenced by Nektar::VCSMapping::v_InitObject().
      
  | 
  protectedvirtual | 
Forcing term for Poisson solver solver
Reimplemented in Nektar::VCSMapping.
Definition at line 370 of file VelocityCorrectionScheme.cpp.
References Nektar::MultiRegions::DirCartesianMap, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::IncNavierStokes::m_velocity, Vmath::Smul(), Vmath::Vadd(), and Vmath::Zero().
Referenced by SetUpPressureForcing(), and Nektar::VCSMapping::v_SetUpPressureForcing().
      
  | 
  protectedvirtual | 
Forcing term for Helmholtz solver
Reimplemented in Nektar::VCSMapping.
Definition at line 394 of file VelocityCorrectionScheme.cpp.
References Nektar::SolverUtils::EquationSystem::m_fields, Nektar::IncNavierStokes::m_kinvis, Nektar::IncNavierStokes::m_pressure, and Nektar::IncNavierStokes::m_velocity.
Referenced by SetUpViscousForcing().
      
  | 
  protectedvirtual | 
Solve pressure system
Reimplemented in Nektar::VCSMapping.
Definition at line 429 of file VelocityCorrectionScheme.cpp.
References Nektar::StdRegions::eFactorLambda, Nektar::IncNavierStokes::m_pressure, and Nektar::NullFlagList.
Referenced by SolvePressure(), and Nektar::VCSMapping::v_SolvePressure().
      
  | 
  protectedvirtual | 
Solve velocity system
Reimplemented in Nektar::VCSMapping.
Definition at line 444 of file VelocityCorrectionScheme.cpp.
References Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::eFactorSVVCutoffRatio, Nektar::StdRegions::eFactorSVVDiffCoeff, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::IncNavierStokes::m_kinvis, Nektar::IncNavierStokes::m_nConvectiveFields, m_sVVCutoffRatio, m_sVVDiffCoeff, m_useSpecVanVisc, and Nektar::NullFlagList.
Referenced by SolveViscous(), and Nektar::VCSMapping::v_SolveViscous().
Virtual function for transformation to physical space.
Reimplemented from Nektar::IncNavierStokes.
Definition at line 256 of file VelocityCorrectionScheme.cpp.
References Nektar::SolverUtils::EquationSystem::m_fields.
Virtual function for transformation to coefficient space.
Reimplemented from Nektar::IncNavierStokes.
Definition at line 270 of file VelocityCorrectionScheme.cpp.
References Nektar::SolverUtils::EquationSystem::m_fields.
      
  | 
  static | 
Name of class.
Definition at line 58 of file VelocityCorrectionScheme.h.
Save aiiDt value to use as a trip to reset global matrix setup.
Definition at line 122 of file VelocityCorrectionScheme.h.
Referenced by v_InitObject().
      
  | 
  protected | 
cutt off ratio from which to start decayhing modes
Definition at line 117 of file VelocityCorrectionScheme.h.
Referenced by v_GenerateSummary(), v_InitObject(), Nektar::VCSMapping::v_SolveViscous(), and v_SolveViscous().
      
  | 
  protected | 
Diffusion coefficient of SVV modes.
Definition at line 119 of file VelocityCorrectionScheme.h.
Referenced by v_GenerateSummary(), v_InitObject(), Nektar::VCSMapping::v_SolveViscous(), and v_SolveViscous().
      
  | 
  protected | 
bool to identify if spectral vanishing viscosity is active.
Definition at line 113 of file VelocityCorrectionScheme.h.
Referenced by v_GenerateSummary(), and v_InitObject().
      
  | 
  protected | 
bool to identify if spectral vanishing viscosity is active.
Definition at line 115 of file VelocityCorrectionScheme.h.
Referenced by v_GenerateSummary(), v_InitObject(), Nektar::VCSMapping::v_SolveViscous(), and v_SolveViscous().
      
  | 
  protected | 
Variable Coefficient map for the Laplacian which can be activated as part of SVV or otherwise.
Definition at line 125 of file VelocityCorrectionScheme.h.
 1.8.8