Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | List of all members
Nektar::VelocityCorrectionScheme Class Reference

#include <VelocityCorrectionScheme.h>

Inheritance diagram for Nektar::VelocityCorrectionScheme:
Inheritance graph
[legend]
Collaboration diagram for Nektar::VelocityCorrectionScheme:
Collaboration graph
[legend]

Public Member Functions

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

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)
 
virtual bool v_RequireFwdTrans ()
 
virtual std::string v_GetExtrapolateStr (void)
 
virtual std::string v_GetSubSteppingExtrapolateStr (const std::string &instr)
 
- 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)
 
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...
 
void SetWomersleyBoundary (const int fldid, const int bndid)
 Set Womersley Profile if specified. More...
 
void SetUpWomersley (const int bndid, std::string womstr)
 Set Up Womersley details. 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, int &nchk)
 
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 EvaluateFunctionExp (std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, const NekDouble &pTime=0.0, const int domain=0)
 
SOLVER_UTILS_EXPORT void EvaluateFunctionFld (std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, const NekDouble &pTime=0.0, const int domain=0)
 
SOLVER_UTILS_EXPORT void EvaluateFunctionPts (std::string pFieldName, Array< OneD, NekDouble > &pArray, const std::string &pFunctionName, const NekDouble &pTime=0.0, const int domain=0)
 
SOLVER_UTILS_EXPORT void LoadPts (std::string funcFilename, std::string filename, Nektar::LibUtilities::PtsFieldSharedPtr &outPts)
 
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, NekDoublem_diffCoeff
 Diffusion coefficients (will be kinvis for velocities) More...
 
StdRegions::VarCoeffMap m_varCoeffLap
 Variable Coefficient map for the Laplacian which can be activated as part of SVV or otherwise. More...
 
Array< OneD, Array< OneD,
NekDouble > > 
m_F
 
- 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...
 
std::map< int,
WomersleyParamsSharedPtr
m_womersleyParams
 Womersley parameters if required. 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...
 
int m_nanSteps
 
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
LibUtilities::TimeIntegrationSolutionSharedPtr m_intSoln
 
NekDouble m_epsilon
 
bool m_explicitDiffusion
 Indicates if explicit or implicit treatment of diffusion is used. More...
 
bool m_explicitAdvection
 Indicates if explicit or implicit treatment of advection is used. More...
 
bool m_explicitReaction
 Indicates if explicit or implicit treatment of reaction is used. More...
 
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state. More...
 
std::vector< int > m_intVariables
 
std::vector< FilterSharedPtrm_filters
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output. More...
 
std::map< std::string,
FieldUtils::Interpolator
m_interpolators
 Map of interpolator objects. More...
 
std::map< std::string,
std::pair< std::string,
LibUtilities::PtsFieldSharedPtr > > 
m_loadedPtsFields
 pts fields we already read from disk: {funcFilename: (filename, ptsfield)} More...
 
std::map< std::string,
std::pair< std::string,
loadedFldField > > 
m_loadedFldFields
 
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...
 
NekDouble m_checktime
 Time between checkpoints. More...
 
int m_nchk
 Number of checkpoints written so far. More...
 
int m_steps
 Number of steps to take. More...
 
int m_checksteps
 Number of steps between checkpoints. More...
 
int m_spacedim
 Spatial dimension (>= expansion dim). More...
 
int m_expdim
 Expansion dimension. More...
 
bool m_singleMode
 Flag to determine if single homogeneous mode is used. More...
 
bool m_halfMode
 Flag to determine if half homogeneous mode is used. More...
 
bool m_multipleModes
 Flag to determine if use multiple homogenenous modes are used. More...
 
bool m_useFFT
 Flag to determine if FFT is used for homogeneous transform. More...
 
bool m_homogen_dealiasing
 Flag to determine if dealiasing is used for homogeneous simulations. More...
 
bool m_specHP_dealiasing
 Flag to determine if dealisising is usde for the Spectral/hp element discretisation. More...
 
enum MultiRegions::ProjectionType m_projectionType
 Type of projection; e.g continuous or discontinuous. More...
 
Array< OneD, Array< OneD,
NekDouble > > 
m_traceNormals
 Array holding trace normals for DG simulations in the forwards direction. More...
 
Array< OneD, Array< OneD,
Array< OneD, NekDouble > > > 
m_gradtan
 1 x nvariable x nq More...
 
Array< OneD, Array< OneD,
Array< OneD, NekDouble > > > 
m_tanbasis
 2 x m_spacedim x nq More...
 
Array< OneD, bool > m_checkIfSystemSingular
 Flag to indicate if the fields should be checked for singularity. More...
 
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
 Map to identify relevant solver info to dump in output fields. More...
 
int m_NumQuadPointsError
 Number of Quadrature points used to work out the error. More...
 
enum HomogeneousType m_HomogeneousType
 
NekDouble m_LhomX
 physical length in X direction (if homogeneous) More...
 
NekDouble m_LhomY
 physical length in Y direction (if homogeneous) More...
 
NekDouble m_LhomZ
 physical length in Z direction (if homogeneous) More...
 
int m_npointsX
 number of points in X direction (if homogeneous) More...
 
int m_npointsY
 number of points in Y direction (if homogeneous) More...
 
int m_npointsZ
 number of points in Z direction (if homogeneous) More...
 
int m_HomoDirec
 number of homogenous directions More...
 

Additional Inherited Members

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

Detailed Description

Definition at line 43 of file VelocityCorrectionScheme.h.

Constructor & Destructor Documentation

Nektar::VelocityCorrectionScheme::VelocityCorrectionScheme ( const LibUtilities::SessionReaderSharedPtr pSession)

Constructor.

Constructor. Creates ...

Parameters
param

Definition at line 59 of file VelocityCorrectionScheme.cpp.

61  : UnsteadySystem(pSession),
62  IncNavierStokes(pSession),
64  {
65 
66  }
SOLVER_UTILS_EXPORT UnsteadySystem(const LibUtilities::SessionReaderSharedPtr &pSession)
Initialises UnsteadySystem class members.
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:228
IncNavierStokes(const LibUtilities::SessionReaderSharedPtr &pSession)
Constructor.
StdRegions::VarCoeffMap m_varCoeffLap
Variable Coefficient map for the Laplacian which can be activated as part of SVV or otherwise...
Nektar::VelocityCorrectionScheme::~VelocityCorrectionScheme ( void  )
virtual

Destructor

Definition at line 192 of file VelocityCorrectionScheme.cpp.

193  {
194  }

Member Function Documentation

static SolverUtils::EquationSystemSharedPtr Nektar::VelocityCorrectionScheme::create ( const LibUtilities::SessionReaderSharedPtr pSession)
inlinestatic

Creates an instance of this class.

Definition at line 48 of file VelocityCorrectionScheme.h.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and CellMLToNektar.cellml_metadata::p.

49  {
52  AllocateSharedPtr(pSession);
53  p->InitObject();
54  return p;
55  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
boost::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.
void Nektar::VelocityCorrectionScheme::EvaluateAdvection_SetPressureBCs ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time 
)
inline

Definition at line 102 of file VelocityCorrectionScheme.h.

References v_EvaluateAdvection_SetPressureBCs().

Referenced by v_InitObject().

106  {
107  v_EvaluateAdvection_SetPressureBCs( inarray, outarray, time);
108  }
virtual void v_EvaluateAdvection_SetPressureBCs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
void Nektar::VelocityCorrectionScheme::SetUpPressureForcing ( const Array< OneD, const Array< OneD, NekDouble > > &  fields,
Array< OneD, Array< OneD, NekDouble > > &  Forcing,
const NekDouble  aii_Dt 
)
inline

Definition at line 67 of file VelocityCorrectionScheme.h.

References v_SetUpPressureForcing().

Referenced by SolveUnsteadyStokesSystem().

71  {
72  v_SetUpPressureForcing( fields, Forcing, aii_Dt);
73  }
virtual void v_SetUpPressureForcing(const Array< OneD, const Array< OneD, NekDouble > > &fields, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt)
void Nektar::VelocityCorrectionScheme::SetUpViscousForcing ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  Forcing,
const NekDouble  aii_Dt 
)
inline

Definition at line 75 of file VelocityCorrectionScheme.h.

References v_SetUpViscousForcing().

Referenced by SolveUnsteadyStokesSystem().

79  {
80  v_SetUpViscousForcing( inarray, Forcing, aii_Dt);
81  }
virtual void v_SetUpViscousForcing(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt)
void Nektar::VelocityCorrectionScheme::SolvePressure ( const Array< OneD, NekDouble > &  Forcing)
inline

Definition at line 83 of file VelocityCorrectionScheme.h.

References v_SolvePressure().

Referenced by SolveUnsteadyStokesSystem().

84  {
85  v_SolvePressure( Forcing);
86  }
virtual void v_SolvePressure(const Array< OneD, NekDouble > &Forcing)
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 347 of file VelocityCorrectionScheme.cpp.

References Nektar::IncNavierStokes::m_extrapolation, m_F, Nektar::IncNavierStokes::m_kinvis, SetUpPressureForcing(), SetUpViscousForcing(), SolvePressure(), and SolveViscous().

Referenced by v_InitObject().

352  {
353  // Substep the pressure boundary condition if using substepping
354  m_extrapolation->SubStepSetPressureBCs(inarray,aii_Dt,m_kinvis);
355 
356  // Set up forcing term for pressure Poisson equation
357  SetUpPressureForcing(inarray, m_F, aii_Dt);
358 
359  // Solve Pressure System
360  SolvePressure (m_F[0]);
361 
362  // Set up forcing term for Helmholtz problems
363  SetUpViscousForcing(inarray, m_F, aii_Dt);
364 
365  // Solve velocity system
366  SolveViscous( m_F, outarray, aii_Dt);
367  }
NekDouble m_kinvis
Kinematic viscosity.
void SolvePressure(const Array< OneD, NekDouble > &Forcing)
ExtrapolateSharedPtr m_extrapolation
void SolveViscous(const Array< OneD, const Array< OneD, NekDouble > > &Forcing, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble aii_Dt)
void SetUpViscousForcing(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt)
void SetUpPressureForcing(const Array< OneD, const Array< OneD, NekDouble > > &fields, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt)
Array< OneD, Array< OneD, NekDouble > > m_F
void Nektar::VelocityCorrectionScheme::SolveViscous ( const Array< OneD, const Array< OneD, NekDouble > > &  Forcing,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  aii_Dt 
)
inline

Definition at line 88 of file VelocityCorrectionScheme.h.

References v_SolveViscous().

Referenced by SolveUnsteadyStokesSystem().

92  {
93  v_SolveViscous( Forcing, outarray, aii_Dt);
94  }
virtual void v_SolveViscous(const Array< OneD, const Array< OneD, NekDouble > > &Forcing, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble aii_Dt)
void Nektar::VelocityCorrectionScheme::v_DoInitialise ( void  )
protectedvirtual

Sets up initial conditions.

Sets the initial conditions.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Reimplemented in Nektar::VCSMapping.

Definition at line 240 of file VelocityCorrectionScheme.cpp.

References m_F, 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().

241  {
242 
244 
245  // Set up Field Meta Data for output files
246  m_fieldMetaDataMap["Kinvis"] = boost::lexical_cast<std::string>(m_kinvis);
247  m_fieldMetaDataMap["TimeStep"] = boost::lexical_cast<std::string>(m_timestep);
248 
249  // set boundary conditions here so that any normal component
250  // correction are imposed before they are imposed on intiial
251  // field below
253 
254  m_F = Array<OneD, Array< OneD, NekDouble> > (m_nConvectiveFields);
255  for(int i = 0; i < m_nConvectiveFields; ++i)
256  {
257  m_fields[i]->LocalToGlobal();
258  m_fields[i]->ImposeDirichletConditions(m_fields[i]->UpdateCoeffs());
259  m_fields[i]->GlobalToLocal();
260  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
261  m_fields[i]->UpdatePhys());
262  m_F[i] = Array< OneD, NekDouble> (m_fields[0]->GetTotPoints(), 0.0);
263  }
264  }
void SetBoundaryConditions(NekDouble time)
time dependent boundary conditions updating
NekDouble m_time
Current time of simulation.
NekDouble m_kinvis
Kinematic viscosity.
NekDouble m_timestep
Time step size.
int m_nConvectiveFields
Number of fields to be convected;.
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
Map to identify relevant solver info to dump in output fields.
Array< OneD, Array< OneD, NekDouble > > m_F
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
virtual SOLVER_UTILS_EXPORT void v_DoInitialise()
Sets up initial conditions.
void Nektar::VelocityCorrectionScheme::v_EvaluateAdvection_SetPressureBCs ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time 
)
protectedvirtual

Explicit part of the method - Advection, Forcing + HOPBCs

Reimplemented in Nektar::VCSMapping.

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

321  {
322  EvaluateAdvectionTerms(inarray, outarray);
323 
324  // Smooth advection
326  {
327  for(int i = 0; i < m_nConvectiveFields; ++i)
328  {
329  m_pressure->SmoothField(outarray[i]);
330  }
331  }
332 
333  // Add forcing terms
334  std::vector<SolverUtils::ForcingSharedPtr>::const_iterator x;
335  for (x = m_forcing.begin(); x != m_forcing.end(); ++x)
336  {
337  (*x)->Apply(m_fields, inarray, outarray, time);
338  }
339 
340  // Calculate High-Order pressure boundary conditions
341  m_extrapolation->EvaluatePressureBCs(inarray,outarray,m_kinvis);
342  }
NekDouble m_kinvis
Kinematic viscosity.
ExtrapolateSharedPtr m_extrapolation
void EvaluateAdvectionTerms(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
int m_nConvectiveFields
Number of fields to be convected;.
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Forcing terms.
bool m_SmoothAdvection
bool to identify if advection term smoothing is requested
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
MultiRegions::ExpListSharedPtr m_pressure
Pointer to field holding pressure field.
void Nektar::VelocityCorrectionScheme::v_GenerateSummary ( SolverUtils::SummaryList s)
protectedvirtual

Print a summary of time stepping parameters.

Prints a summary with some information regards the time-stepping.

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Reimplemented in Nektar::VCSWeakPressure.

Definition at line 199 of file VelocityCorrectionScheme.cpp.

References Nektar::SolverUtils::AddSummaryItem(), Nektar::SolverUtils::EquationSystem::eHomogeneous1D, Nektar::LibUtilities::eNoTimeIntegrationMethod, Nektar::IncNavierStokes::m_extrapolation, Nektar::SolverUtils::EquationSystem::m_homogen_dealiasing, Nektar::SolverUtils::EquationSystem::m_HomogeneousType, Nektar::SolverUtils::EquationSystem::m_specHP_dealiasing, m_sVVCutoffRatio, m_sVVDiffCoeff, m_useHomo1DSpecVanVisc, m_useSpecVanVisc, Nektar::LibUtilities::TimeIntegrationMethodMap, and Nektar::SolverUtils::UnsteadySystem::v_GenerateSummary().

200  {
202  SolverUtils::AddSummaryItem(s, "Splitting Scheme", "Velocity correction (strong press. form)");
203 
204  if (m_extrapolation->GetSubStepIntegrationMethod() !=
206  {
207  SolverUtils::AddSummaryItem(s, "Substepping",
209  m_extrapolation->GetSubStepIntegrationMethod()]);
210  }
211 
212  string dealias = m_homogen_dealiasing ? "Homogeneous1D" : "";
214  {
215  dealias += (dealias == "" ? "" : " + ") + string("spectral/hp");
216  }
217  if (dealias != "")
218  {
219  SolverUtils::AddSummaryItem(s, "Dealiasing", dealias);
220  }
221 
222  string smoothing = m_useSpecVanVisc ? "spectral/hp" : "";
224  {
225  smoothing += (smoothing == "" ? "" : " + ") + string("Homogeneous1D");
226  }
227  if (smoothing != "")
228  {
230  s, "Smoothing", "SVV (" + smoothing + " SVV (cut-off = "
231  + boost::lexical_cast<string>(m_sVVCutoffRatio)
232  + ", diff coeff = "
233  + boost::lexical_cast<string>(m_sVVDiffCoeff)+")");
234  }
235  }
bool m_useHomo1DSpecVanVisc
bool to identify if spectral vanishing viscosity is active.
NekDouble m_sVVDiffCoeff
Diffusion coefficient of SVV modes.
ExtrapolateSharedPtr m_extrapolation
bool m_specHP_dealiasing
Flag to determine if dealisising is usde for the Spectral/hp element discretisation.
const char *const TimeIntegrationMethodMap[]
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s)
Print a summary of time stepping parameters.
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:50
bool m_homogen_dealiasing
Flag to determine if dealiasing is used for homogeneous simulations.
NekDouble m_sVVCutoffRatio
cutt off ratio from which to start decayhing modes
bool m_useSpecVanVisc
bool to identify if spectral vanishing viscosity is active.
enum HomogeneousType m_HomogeneousType
virtual std::string Nektar::VelocityCorrectionScheme::v_GetExtrapolateStr ( void  )
inlineprotectedvirtual

Reimplemented in Nektar::VCSWeakPressure.

Definition at line 165 of file VelocityCorrectionScheme.h.

Referenced by v_InitObject().

166  {
167  return "Standard";
168  }
int Nektar::VelocityCorrectionScheme::v_GetForceDimension ( void  )
protectedvirtual

Implements Nektar::IncNavierStokes.

Definition at line 309 of file VelocityCorrectionScheme.cpp.

References Nektar::SolverUtils::EquationSystem::m_session.

310  {
311  return m_session->GetVariables().size() - 1;
312  }
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
virtual std::string Nektar::VelocityCorrectionScheme::v_GetSubSteppingExtrapolateStr ( const std::string &  instr)
inlineprotectedvirtual

Reimplemented in Nektar::VCSWeakPressure.

Definition at line 170 of file VelocityCorrectionScheme.h.

Referenced by v_InitObject().

172  {
173  return instr;
174  }
Array< OneD, bool > Nektar::VelocityCorrectionScheme::v_GetSystemSingularChecks ( )
protectedvirtual

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 298 of file VelocityCorrectionScheme.cpp.

References Nektar::SolverUtils::EquationSystem::m_session.

299  {
300  int vVar = m_session->GetVariables().size();
301  Array<OneD, bool> vChecks(vVar, false);
302  vChecks[vVar-1] = true;
303  return vChecks;
304  }
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
void Nektar::VelocityCorrectionScheme::v_InitObject ( )
virtual

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::IncNavierStokes.

Reimplemented in Nektar::VCSMapping.

Definition at line 68 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::SolverUtils::AdvectionSystem::m_advObject, m_diffCoeff, 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, 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(), v_GetExtrapolateStr(), v_GetSubSteppingExtrapolateStr(), and Nektar::IncNavierStokes::v_InitObject().

Referenced by Nektar::VCSMapping::v_InitObject().

69  {
70  int n;
71 
73  m_explicitDiffusion = false;
74 
75  // Set m_pressure to point to last field of m_fields;
76  if (boost::iequals(m_session->GetVariable(m_fields.num_elements()-1), "p"))
77  {
78  m_nConvectiveFields = m_fields.num_elements()-1;
80  }
81  else
82  {
83  ASSERTL0(false,"Need to set up pressure field definition");
84  }
85 
86  // Determine diffusion coefficients for each field
87  m_diffCoeff = Array<OneD, NekDouble> (m_nConvectiveFields, m_kinvis);
88  for (n = 0; n < m_nConvectiveFields; ++n)
89  {
90  std::string varName = m_session->GetVariable(n);
91  if ( m_session->DefinesFunction("DiffusionCoefficient", varName))
92  {
94  = m_session->GetFunction("DiffusionCoefficient", varName);
95  m_diffCoeff[n] = ffunc->Evaluate();
96  }
97  }
98 
99  // creation of the extrapolation object
101  {
102  std::string vExtrapolation = v_GetExtrapolateStr();
103 
104  if (m_session->DefinesSolverInfo("Extrapolation"))
105  {
106  vExtrapolation = v_GetSubSteppingExtrapolateStr(
107  m_session->GetSolverInfo("Extrapolation"));
108  }
109 
111  vExtrapolation,
112  m_session,
113  m_fields,
114  m_pressure,
115  m_velocity,
116  m_advObject);
117  }
118 
119  // Integrate only the convective fields
120  for (n = 0; n < m_nConvectiveFields; ++n)
121  {
122  m_intVariables.push_back(n);
123  }
124 
125  // Load parameters for Spectral Vanishing Viscosity
126  m_session->MatchSolverInfo("SpectralVanishingViscosity","True",
127  m_useSpecVanVisc, false);
129  if(m_useSpecVanVisc == false)
130  {
131  m_session->MatchSolverInfo("SpectralVanishingViscositySpectralHP",
132  "True", m_useSpecVanVisc, false);
133  m_session->MatchSolverInfo("SpectralVanishingViscosityHomo1D",
134  "True", m_useHomo1DSpecVanVisc, false);
135  }
136  m_session->LoadParameter("SVVCutoffRatio",m_sVVCutoffRatio,0.75);
137  m_session->LoadParameter("SVVDiffCoeff", m_sVVDiffCoeff, 0.1);
138 
139  m_session->MatchSolverInfo("SPECTRALHPDEALIASING","True",
140  m_specHP_dealiasing,false);
141 
143  {
144  ASSERTL0(m_nConvectiveFields > 2,"Expect to have three velocity fields with homogenous expansion");
145 
147  {
148  Array<OneD, unsigned int> planes;
149  planes = m_fields[0]->GetZIDs();
150 
151  int num_planes = planes.num_elements();
152  Array<OneD, NekDouble> SVV(num_planes,0.0);
153  NekDouble fac;
154  int kmodes = m_fields[0]->GetHomogeneousBasis()->GetNumModes();
155  int pstart;
156 
157  pstart = m_sVVCutoffRatio*kmodes;
158 
159  for(n = 0; n < num_planes; ++n)
160  {
161  if(planes[n] > pstart)
162  {
163  fac = (NekDouble)((planes[n] - kmodes)*(planes[n] - kmodes))/
164  ((NekDouble)((planes[n] - pstart)*(planes[n] - pstart)));
165  SVV[n] = m_sVVDiffCoeff*exp(-fac)/m_kinvis;
166  }
167  }
168 
169  for(int i = 0; i < m_velocity.num_elements(); ++i)
170  {
171  m_fields[m_velocity[i]]->SetHomo1DSpecVanVisc(SVV);
172  }
173  }
174 
175  }
176 
177  m_session->MatchSolverInfo("SmoothAdvection", "True", m_SmoothAdvection, false);
178 
179  // set explicit time-intregration class operators
181 
182  m_extrapolation->SubSteppingTimeIntegration(m_intScheme->GetIntegrationMethod(), m_intScheme);
183  m_extrapolation->GenerateHOPBCMap(m_session);
184 
185  // set implicit time-intregration class operators
187  }
EquationType m_equationType
equation type;
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:198
void EvaluateAdvection_SetPressureBCs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
tBaseSharedPtr CreateInstance(tKey idKey BOOST_PP_COMMA_IF(MAX_PARAM) BOOST_PP_ENUM_BINARY_PARAMS(MAX_PARAM, tParam, x))
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:162
void SolveUnsteadyStokesSystem(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time, const NekDouble a_iixDt)
bool m_explicitDiffusion
Indicates if explicit or implicit treatment of diffusion is used.
NekDouble m_kinvis
Kinematic viscosity.
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
ExtrapolateFactory & GetExtrapolateFactory()
Definition: Extrapolate.cpp:50
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
Array< OneD, int > m_velocity
int which identifies which components of m_fields contains the velocity (u,v,w);
bool m_useHomo1DSpecVanVisc
bool to identify if spectral vanishing viscosity is active.
NekDouble m_sVVDiffCoeff
Diffusion coefficient of SVV modes.
ExtrapolateSharedPtr m_extrapolation
virtual std::string v_GetSubSteppingExtrapolateStr(const std::string &instr)
virtual std::string v_GetExtrapolateStr(void)
bool m_specHP_dealiasing
Flag to determine if dealisising is usde for the Spectral/hp element discretisation.
int m_nConvectiveFields
Number of fields to be convected;.
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
bool m_SmoothAdvection
bool to identify if advection term smoothing is requested
double NekDouble
boost::shared_ptr< Equation > EquationSharedPtr
NekDouble m_sVVCutoffRatio
cutt off ratio from which to start decayhing modes
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
bool m_useSpecVanVisc
bool to identify if spectral vanishing viscosity is active.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
Array< OneD, NekDouble > m_diffCoeff
Diffusion coefficients (will be kinvis for velocities)
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
Wrapper to the time integration scheme.
MultiRegions::ExpListSharedPtr m_pressure
Pointer to field holding pressure field.
virtual void v_InitObject()
Init object for UnsteadySystem class.
enum HomogeneousType m_HomogeneousType
virtual bool Nektar::VelocityCorrectionScheme::v_RequireFwdTrans ( )
inlineprotectedvirtual

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 160 of file VelocityCorrectionScheme.h.

161  {
162  return false;
163  }
void Nektar::VelocityCorrectionScheme::v_SetUpPressureForcing ( const Array< OneD, const Array< OneD, NekDouble > > &  fields,
Array< OneD, Array< OneD, NekDouble > > &  Forcing,
const NekDouble  aii_Dt 
)
protectedvirtual

Forcing term for Poisson solver solver

Reimplemented in Nektar::VCSMapping, and Nektar::VCSWeakPressure.

Definition at line 372 of file VelocityCorrectionScheme.cpp.

References Nektar::MultiRegions::DirCartesianMap, Nektar::MultiRegions::eX, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::IncNavierStokes::m_velocity, Vmath::Smul(), and Vmath::Vadd().

Referenced by SetUpPressureForcing(), and Nektar::VCSMapping::v_SetUpPressureForcing().

376  {
377  int i;
378  int physTot = m_fields[0]->GetTotPoints();
379  int nvel = m_velocity.num_elements();
380 
381  m_fields[0]->PhysDeriv(eX,fields[0], Forcing[0]);
382 
383  for(i = 1; i < nvel; ++i)
384  {
385  // Use Forcing[1] as storage since it is not needed for the pressure
386  m_fields[i]->PhysDeriv(DirCartesianMap[i],fields[i],Forcing[1]);
387  Vmath::Vadd(physTot,Forcing[1],1,Forcing[0],1,Forcing[0],1);
388  }
389 
390  Vmath::Smul(physTot,1.0/aii_Dt,Forcing[0],1,Forcing[0],1);
391  }
Array< OneD, int > m_velocity
int which identifies which components of m_fields contains the velocity (u,v,w);
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:213
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:86
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:299
void Nektar::VelocityCorrectionScheme::v_SetUpViscousForcing ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  Forcing,
const NekDouble  aii_Dt 
)
protectedvirtual

Forcing term for Helmholtz solver

Reimplemented in Nektar::VCSMapping.

Definition at line 396 of file VelocityCorrectionScheme.cpp.

References m_diffCoeff, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::IncNavierStokes::m_nConvectiveFields, Nektar::IncNavierStokes::m_pressure, Nektar::IncNavierStokes::m_velocity, and Vmath::Zero().

Referenced by SetUpViscousForcing().

400  {
401  NekDouble aii_dtinv = 1.0/aii_Dt;
402  int phystot = m_fields[0]->GetTotPoints();
403 
404  // Grad p
405  m_pressure->BwdTrans(m_pressure->GetCoeffs(),m_pressure->UpdatePhys());
406 
407  int nvel = m_velocity.num_elements();
408  if(nvel == 2)
409  {
410  m_pressure->PhysDeriv(m_pressure->GetPhys(),
411  Forcing[m_velocity[0]],
412  Forcing[m_velocity[1]]);
413  }
414  else
415  {
416  m_pressure->PhysDeriv(m_pressure->GetPhys(),
417  Forcing[m_velocity[0]],
418  Forcing[m_velocity[1]],
419  Forcing[m_velocity[2]]);
420  }
421 
422  // zero convective fields.
423  for(int i = nvel; i < m_nConvectiveFields; ++i)
424  {
425  Vmath::Zero(phystot,Forcing[i],1);
426  }
427 
428  // Subtract inarray/(aii_dt) and divide by kinvis. Kinvis will
429  // need to be updated for the convected fields.
430  for(int i = 0; i < m_nConvectiveFields; ++i)
431  {
432  Blas::Daxpy(phystot,-aii_dtinv,inarray[i],1,Forcing[i],1);
433  Blas::Dscal(phystot,1.0/m_diffCoeff[i],&(Forcing[i])[0],1);
434  }
435  }
Array< OneD, int > m_velocity
int which identifies which components of m_fields contains the velocity (u,v,w);
int m_nConvectiveFields
Number of fields to be convected;.
double NekDouble
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
Array< OneD, NekDouble > m_diffCoeff
Diffusion coefficients (will be kinvis for velocities)
MultiRegions::ExpListSharedPtr m_pressure
Pointer to field holding pressure field.
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:373
void Nektar::VelocityCorrectionScheme::v_SolvePressure ( const Array< OneD, NekDouble > &  Forcing)
protectedvirtual

Solve pressure system

Reimplemented in Nektar::VCSMapping, and Nektar::VCSWeakPressure.

Definition at line 441 of file VelocityCorrectionScheme.cpp.

References Nektar::StdRegions::eFactorLambda, Nektar::IncNavierStokes::m_extrapolation, Nektar::IncNavierStokes::m_kinvis, Nektar::IncNavierStokes::m_pressure, and Nektar::NullFlagList.

Referenced by SolvePressure(), and Nektar::VCSMapping::v_SolvePressure().

443  {
445  // Setup coefficient for equation
446  factors[StdRegions::eFactorLambda] = 0.0;
447 
448  // Solver Pressure Poisson Equation
449  m_pressure->HelmSolve(Forcing, m_pressure->UpdateCoeffs(),
450  NullFlagList, factors);
451 
452  // Add presure to outflow bc if using convective like BCs
453  m_extrapolation->AddPressureToOutflowBCs(m_kinvis);
454  }
NekDouble m_kinvis
Kinematic viscosity.
ExtrapolateSharedPtr m_extrapolation
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:252
MultiRegions::ExpListSharedPtr m_pressure
Pointer to field holding pressure field.
static FlagList NullFlagList
An empty flag list.
void Nektar::VelocityCorrectionScheme::v_SolveViscous ( const Array< OneD, const Array< OneD, NekDouble > > &  Forcing,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  aii_Dt 
)
protectedvirtual

Solve velocity system

Reimplemented in Nektar::VCSMapping.

Definition at line 459 of file VelocityCorrectionScheme.cpp.

References Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::eFactorSVVCutoffRatio, Nektar::StdRegions::eFactorSVVDiffCoeff, m_diffCoeff, 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().

463  {
465 
466  if(m_useSpecVanVisc)
467  {
470  }
471 
472  // Solve Helmholtz system and put in Physical space
473  for(int i = 0; i < m_nConvectiveFields; ++i)
474  {
475  // Setup coefficients for equation
476  factors[StdRegions::eFactorLambda] = 1.0/aii_Dt/m_diffCoeff[i];
477  m_fields[i]->HelmSolve(Forcing[i], m_fields[i]->UpdateCoeffs(),
478  NullFlagList, factors);
479  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),outarray[i]);
480  }
481  }
NekDouble m_kinvis
Kinematic viscosity.
NekDouble m_sVVDiffCoeff
Diffusion coefficient of SVV modes.
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:252
int m_nConvectiveFields
Number of fields to be convected;.
NekDouble m_sVVCutoffRatio
cutt off ratio from which to start decayhing modes
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
bool m_useSpecVanVisc
bool to identify if spectral vanishing viscosity is active.
Array< OneD, NekDouble > m_diffCoeff
Diffusion coefficients (will be kinvis for velocities)
static FlagList NullFlagList
An empty flag list.
void Nektar::VelocityCorrectionScheme::v_TransCoeffToPhys ( void  )
protectedvirtual

Virtual function for transformation to physical space.

Reimplemented from Nektar::IncNavierStokes.

Definition at line 270 of file VelocityCorrectionScheme.cpp.

References Nektar::SolverUtils::EquationSystem::m_fields.

271  {
272  int nfields = m_fields.num_elements() - 1;
273  for (int k=0 ; k < nfields; ++k)
274  {
275  //Backward Transformation in physical space for time evolution
276  m_fields[k]->BwdTrans_IterPerExp(m_fields[k]->GetCoeffs(),
277  m_fields[k]->UpdatePhys());
278  }
279  }
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
void Nektar::VelocityCorrectionScheme::v_TransPhysToCoeff ( void  )
protectedvirtual

Virtual function for transformation to coefficient space.

Reimplemented from Nektar::IncNavierStokes.

Definition at line 284 of file VelocityCorrectionScheme.cpp.

References Nektar::SolverUtils::EquationSystem::m_fields.

285  {
286 
287  int nfields = m_fields.num_elements() - 1;
288  for (int k=0 ; k < nfields; ++k)
289  {
290  //Forward Transformation in physical space for time evolution
291  m_fields[k]->FwdTrans_IterPerExp(m_fields[k]->GetPhys(),m_fields[k]->UpdateCoeffs());
292  }
293  }
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.

Member Data Documentation

string Nektar::VelocityCorrectionScheme::className
static
Initial value:

Name of class.

Definition at line 58 of file VelocityCorrectionScheme.h.

Array<OneD, NekDouble> Nektar::VelocityCorrectionScheme::m_diffCoeff
protected

Diffusion coefficients (will be kinvis for velocities)

Definition at line 120 of file VelocityCorrectionScheme.h.

Referenced by v_InitObject(), v_SetUpViscousForcing(), and v_SolveViscous().

Array<OneD, Array< OneD, NekDouble> > Nektar::VelocityCorrectionScheme::m_F
protected
NekDouble Nektar::VelocityCorrectionScheme::m_sVVCutoffRatio
protected

cutt off ratio from which to start decayhing modes

Definition at line 116 of file VelocityCorrectionScheme.h.

Referenced by Nektar::VCSWeakPressure::v_GenerateSummary(), v_GenerateSummary(), v_InitObject(), Nektar::VCSMapping::v_SolveViscous(), and v_SolveViscous().

NekDouble Nektar::VelocityCorrectionScheme::m_sVVDiffCoeff
protected
bool Nektar::VelocityCorrectionScheme::m_useHomo1DSpecVanVisc
protected

bool to identify if spectral vanishing viscosity is active.

Definition at line 112 of file VelocityCorrectionScheme.h.

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

bool Nektar::VelocityCorrectionScheme::m_useSpecVanVisc
protected

bool to identify if spectral vanishing viscosity is active.

Definition at line 114 of file VelocityCorrectionScheme.h.

Referenced by Nektar::VCSWeakPressure::v_GenerateSummary(), v_GenerateSummary(), v_InitObject(), Nektar::VCSMapping::v_SolveViscous(), and v_SolveViscous().

StdRegions::VarCoeffMap Nektar::VelocityCorrectionScheme::m_varCoeffLap
protected

Variable Coefficient map for the Laplacian which can be activated as part of SVV or otherwise.

Definition at line 123 of file VelocityCorrectionScheme.h.