Nektar++
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:
[legend]

Public Member Functions

 VelocityCorrectionScheme (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Constructor. More...
 
virtual ~VelocityCorrectionScheme ()
 
virtual void v_InitObject (bool DeclareField=true)
 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, const Array< OneD, const Array< OneD, NekDouble >> &inarray, 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 ()
 
int GetNConvectiveFields (void)
 
Array< OneD, int > & GetVelocity (void)
 
void AddForcing (const SolverUtils::ForcingSharedPtr &pForce)
 
virtual void GetPressure (const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &pressure) override
 Extract array with pressure from physfield. More...
 
virtual void GetDensity (const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &density) override
 Extract array with density from physfield. More...
 
virtual bool HasConstantDensity () override
 
virtual void GetVelocity (const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, Array< OneD, NekDouble >> &velocity) override
 Extract array with velocity from physfield. More...
 
virtual void SetMovingFrameVelocities (const Array< OneD, NekDouble > &vFrameVels) override
 
virtual void GetMovingFrameVelocities (Array< OneD, NekDouble > &vFrameVels) override
 
virtual void SetMovingFrameAngles (const Array< OneD, NekDouble > &vFrameTheta) override
 
virtual void GetMovingFrameAngles (Array< OneD, NekDouble > &vFrameTheta) override
 
virtual void SetMovingFrameProjectionMat (const bnu::matrix< NekDouble > &vProjMat) override
 
virtual void GetMovingFrameProjectionMat (bnu::matrix< NekDouble > &vProjMat) override
 
virtual bool DefinedForcing (const std::string &sForce)
 
virtual void GetPivotPoint (Array< OneD, NekDouble > &vPivotPoint)
 
- Public Member Functions inherited from Nektar::SolverUtils::AdvectionSystem
SOLVER_UTILS_EXPORT AdvectionSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
virtual SOLVER_UTILS_EXPORT ~AdvectionSystem ()
 
SOLVER_UTILS_EXPORT AdvectionSharedPtr GetAdvObject ()
 Returns the advection object held by this instance. More...
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleGetElmtCFLVals (const bool FlagAcousticCFL=true)
 
SOLVER_UTILS_EXPORT NekDouble GetCFLEstimate (int &elmtid)
 
- 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...
 
SOLVER_UTILS_EXPORT void SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
- 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 (bool DeclareField=true)
 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 >
std::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 ExtraFldOutput (std::vector< Array< OneD, NekDouble >> &fieldcoeffs, std::vector< std::string > &variables)
 
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 SessionFunctionSharedPtr GetFunction (std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
 Get a SessionFunction by name. More...
 
SOLVER_UTILS_EXPORT void SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 Initialise the data in the dependent fields. More...
 
SOLVER_UTILS_EXPORT void EvaluateExactSolution (int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 Evaluates an exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln, bool Normalised=false)
 Compute the L2 error between fields and a given exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, bool Normalised=false)
 Compute the L2 error of the fields. More...
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleErrorExtraPoints (unsigned int field)
 Compute error (L2 and L_inf) over an larger set of quadrature points return [L2 Linf]. More...
 
SOLVER_UTILS_EXPORT void 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 SessionSummary (SummaryList &vSummary)
 Write out a session summary. More...
 
SOLVER_UTILS_EXPORT Array< OneD, MultiRegions::ExpListSharedPtr > & UpdateFields ()
 
SOLVER_UTILS_EXPORT LibUtilities::FieldMetaDataMapUpdateFieldMetaDataMap ()
 Get hold of FieldInfoMap so it can be updated. More...
 
SOLVER_UTILS_EXPORT NekDouble GetFinalTime ()
 Return final time. More...
 
SOLVER_UTILS_EXPORT int GetNcoeffs ()
 
SOLVER_UTILS_EXPORT int GetNcoeffs (const int eid)
 
SOLVER_UTILS_EXPORT int GetNumExpModes ()
 
SOLVER_UTILS_EXPORT const Array< OneD, int > GetNumExpModesPerExp ()
 
SOLVER_UTILS_EXPORT int GetNvariables ()
 
SOLVER_UTILS_EXPORT const std::string GetVariable (unsigned int i)
 
SOLVER_UTILS_EXPORT int GetTraceTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTraceNpoints ()
 
SOLVER_UTILS_EXPORT int GetExpSize ()
 
SOLVER_UTILS_EXPORT int GetPhys_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetCoeff_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTotPoints (int n)
 
SOLVER_UTILS_EXPORT int GetNpoints ()
 
SOLVER_UTILS_EXPORT int 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 SetModifiedBasis (const bool modbasis)
 
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 Array< OneD, const Array< OneD, NekDouble > > GetTraceNormals ()
 
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...
 
- Public Member Functions inherited from Nektar::SolverUtils::FluidInterface
virtual ~FluidInterface ()=default
 
virtual SOLVER_UTILS_EXPORT void SetMovingFrameProjectionMat (const boost::numeric::ublas::matrix< NekDouble > &vProjMat)
 
virtual SOLVER_UTILS_EXPORT void GetMovingFrameProjectionMat (boost::numeric::ublas::matrix< NekDouble > &vProjMat)
 

Static Public Member Functions

static SolverUtils::EquationSystemSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Creates an instance of this class. More...
 

Static Public Attributes

static std::string className
 Name of class. More...
 
- Static Public Attributes inherited from Nektar::SolverUtils::UnsteadySystem
static std::string cmdSetStartTime
 
static std::string cmdSetStartChkNum
 

Protected Member Functions

void SetupFlowrate (NekDouble aii_dt)
 Set up the Stokes solution used to impose constant flowrate through a boundary. More...
 
NekDouble MeasureFlowrate (const Array< OneD, Array< OneD, NekDouble >> &inarray)
 Measure the volumetric flow rate through the volumetric flow rate reference surface. More...
 
virtual bool v_PostIntegrate (int step)
 
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, const Array< OneD, const Array< OneD, NekDouble >> &inarray, 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)
 
void SetUpSVV (void)
 
void SetUpExtrapolation (void)
 
void SVVVarDiffCoeff (const NekDouble velmag, Array< OneD, NekDouble > &diffcoeff, const Array< OneD, Array< OneD, NekDouble >> &vel=NullNekDoubleArrayOfArray)
 
void AppendSVVFactors (StdRegions::ConstFactorMap &factors, MultiRegions::VarFactorsMap &varFactorsMap)
 
- Protected Member Functions inherited from Nektar::IncNavierStokes
 IncNavierStokes (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Constructor. More...
 
EquationType GetEquationType (void)
 
void EvaluateAdvectionTerms (const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
 
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 fldid, const int bndid, std::string womstr)
 Set Up Womersley details. More...
 
void SetMovingReferenceFrameBCs (const NekDouble &time)
 Set the moving reference frame boundary conditions. More...
 
void SetMRFWallBCs (const NekDouble &time)
 
void SetMRFDomainVelBCs (const NekDouble &time)
 
virtual MultiRegions::ExpListSharedPtr v_GetPressure () override
 
virtual Array< OneD, NekDoublev_GetMaxStdVelocity (const NekDouble SpeedSoundFactor) override
 
virtual bool v_PreIntegrate (int step) override
 
- Protected Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
SOLVER_UTILS_EXPORT UnsteadySystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 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 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 void v_SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
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...
 
virtual SOLVER_UTILS_EXPORT bool UpdateTimeStepCheck ()
 
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises EquationSystem class members. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Virtual function for the L_inf error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false)
 Virtual function for the L_2 error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 
virtual SOLVER_UTILS_EXPORT void v_EvaluateExactSolution (unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 
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...
 
bool m_useGJPStabilisation
 bool to identify if GJP semi-implicit is active. More...
 
bool m_useGJPNormalVel
 bool to identify if GJP normal Velocity should be applied in explicit approach More...
 
NekDouble m_GJPJumpScale
 
NekDouble m_sVVCutoffRatio
 cutt off ratio from which to start decayhing modes More...
 
NekDouble m_sVVDiffCoeff
 Diffusion coefficient of SVV modes. More...
 
NekDouble m_sVVCutoffRatioHomo1D
 
NekDouble m_sVVDiffCoeffHomo1D
 Diffusion coefficient of SVV modes in homogeneous 1D Direction. More...
 
Array< OneD, NekDoublem_svvVarDiffCoeff
 Array of coefficient if power kernel is used in SVV. More...
 
bool m_IsSVVPowerKernel
 Identifier for Power Kernel otherwise DG kernel. 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...
 
NekDouble m_flowrate
 Desired volumetric flowrate. More...
 
NekDouble m_flowrateArea
 Area of the boundary through which we are measuring the flowrate. More...
 
bool m_homd1DFlowinPlane
 
NekDouble m_greenFlux
 Flux of the Stokes function solution. More...
 
NekDouble m_alpha
 Current flowrate correction. More...
 
int m_flowrateBndID
 Boundary ID of the flowrate reference surface. More...
 
int m_planeID
 Plane ID for cases with homogeneous expansion. More...
 
MultiRegions::ExpListSharedPtr m_flowrateBnd
 Flowrate reference surface. More...
 
Array< OneD, Array< OneD, NekDouble > > m_flowrateStokes
 Stokes solution used to impose flowrate. More...
 
std::ofstream m_flowrateStream
 Output stream to record flowrate. More...
 
int m_flowrateSteps
 Interval at which to record flowrate data. More...
 
NekDouble m_flowrateAiidt
 Value of aii_dt used to compute Stokes flowrate solution. 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::ForcingSharedPtrm_forcing
 Forcing terms. More...
 
int m_nConvectiveFields
 Number of fields to be convected;. More...
 
Array< OneD, int > m_velocity
 int which identifies which components of m_fields contains the velocity (u,v,w); More...
 
MultiRegions::ExpListSharedPtr m_pressure
 Pointer to field holding pressure field. More...
 
NekDouble m_kinvis
 Kinematic viscosity. More...
 
int m_energysteps
 dump energy to file at steps time More...
 
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...
 
Array< OneD, NekDoublem_pivotPoint
 
std::map< int, 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_abortSteps
 Number of steps between checks for abort conditions. More...
 
int m_filtersInfosteps
 Number of time steps between outputting filters information. More...
 
int m_nanSteps
 
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
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...
 
NekDouble m_steadyStateTol
 Tolerance to which steady state should be evaluated at. More...
 
int m_steadyStateSteps
 Check for steady state at step interval. More...
 
NekDouble m_steadyStateRes = 1.0
 
NekDouble m_steadyStateRes0 = 1.0
 
Array< OneD, Array< OneD, NekDouble > > m_previousSolution
 Storage for previous solution for steady-state check. More...
 
std::ofstream m_errFile
 
std::vector< int > m_intVariables
 
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
 
NekDouble m_filterTimeWarning
 Number of time steps between outputting status information. More...
 
NekDouble m_TimeIntegLambda = 0.0
 coefff of spacial derivatives(rhs or m_F in GLM) in calculating the residual of the whole equation(used in unsteady time integrations) More...
 
bool m_flagImplicitItsStatistics
 
bool m_flagImplicitSolver = false
 
Array< OneD, NekDoublem_magnitdEstimat
 estimate the magnitude of each conserved varibles More...
 
Array< OneD, NekDoublem_locTimeStep
 local time step(notice only for jfnk other see m_cflSafetyFactor) More...
 
NekDouble m_inArrayNorm = -1.0
 
int m_TotLinItePerStep = 0
 
int m_StagesPerStep = 1
 
bool m_flagUpdatePreconMat
 
int m_maxLinItePerNewton
 
int m_TotNewtonIts = 0
 
int m_TotLinIts = 0
 
int m_TotImpStages = 0
 
bool m_CalcPhysicalAV = true
 flag to update artificial viscosity More...
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
bool m_verbose
 
bool m_root
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader. More...
 
std::map< std::string, SolverUtils::SessionFunctionSharedPtrm_sessionFunctions
 Map of known SessionFunctions. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 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_timestepMax = -1.0
 Time step size. More...
 
NekDouble m_lambda
 Lambda constant in real system if one required. More...
 
NekDouble m_checktime
 Time between checkpoints. More...
 
NekDouble m_lastCheckTime
 
NekDouble m_TimeIncrementFactor
 
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, 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...
 
Array< OneD, NekDoublem_movingFrameVelsxyz
 Moving frame of reference velocities. More...
 
Array< OneD, NekDoublem_movingFrameTheta
 Moving frame of reference angles with respect to the. More...
 
boost::numeric::ublas::matrix< NekDoublem_movingFrameProjMat
 Projection matrix for transformation between inertial and moving. 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...
 
NekDouble m_cflNonAcoustic
 
NekDouble m_CFLGrowth
 CFL growth rate. More...
 
NekDouble m_CFLEnd
 maximun cfl in cfl growth More...
 
- Protected Types inherited from Nektar::SolverUtils::EquationSystem
enum  HomogeneousType { eHomogeneous1D , eHomogeneous2D , eHomogeneous3D , eNotHomogeneous }
 Parameter for homogeneous expansions. More...
 
- Static Protected Attributes inherited from Nektar::IncNavierStokes
static std::string eqTypeLookupIds []
 
- Static Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
static std::string equationSystemTypeLookupIds []
 

Detailed Description

Definition at line 42 of file VelocityCorrectionScheme.h.

Constructor & Destructor Documentation

◆ VelocityCorrectionScheme()

Nektar::VelocityCorrectionScheme::VelocityCorrectionScheme ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr pGraph 
)

Constructor.

Constructor. Creates ...

Parameters

Definition at line 61 of file VelocityCorrectionScheme.cpp.

64  : UnsteadySystem(pSession, pGraph), IncNavierStokes(pSession, pGraph),
66 {
67 }
IncNavierStokes(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Constructor.
SOLVER_UTILS_EXPORT UnsteadySystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises UnsteadySystem class members.
StdRegions::VarCoeffMap m_varCoeffLap
Variable Coefficient map for the Laplacian which can be activated as part of SVV or otherwise.
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:241

◆ ~VelocityCorrectionScheme()

Nektar::VelocityCorrectionScheme::~VelocityCorrectionScheme ( void  )
virtual

Destructor

Definition at line 512 of file VelocityCorrectionScheme.cpp.

513 {
514 }

Member Function Documentation

◆ AppendSVVFactors()

void Nektar::VelocityCorrectionScheme::AppendSVVFactors ( StdRegions::ConstFactorMap factors,
MultiRegions::VarFactorsMap varFactorsMap 
)
protected

Definition at line 1154 of file VelocityCorrectionScheme.cpp.

1157 {
1158 
1159  if (m_useSpecVanVisc)
1160  {
1164  {
1165  if (m_IsSVVPowerKernel)
1166  {
1169  }
1170  else
1171  {
1172  varFactorsMap[StdRegions::eFactorSVVDGKerDiffCoeff] =
1174  }
1175  }
1176  }
1177 }
NekDouble m_kinvis
Kinematic viscosity.
Array< OneD, NekDouble > m_svvVarDiffCoeff
Array of coefficient if power kernel is used in SVV.
bool m_IsSVVPowerKernel
Identifier for Power Kernel otherwise DG kernel.
NekDouble m_sVVCutoffRatio
cutt off ratio from which to start decayhing modes
NekDouble m_sVVDiffCoeff
Diffusion coefficient of SVV modes.
bool m_useSpecVanVisc
bool to identify if spectral vanishing viscosity is active.
static Array< OneD, NekDouble > NullNekDouble1DArray

References Nektar::StdRegions::eFactorSVVCutoffRatio, Nektar::StdRegions::eFactorSVVDGKerDiffCoeff, Nektar::StdRegions::eFactorSVVDiffCoeff, Nektar::StdRegions::eFactorSVVPowerKerDiffCoeff, m_IsSVVPowerKernel, Nektar::IncNavierStokes::m_kinvis, m_sVVCutoffRatio, m_sVVDiffCoeff, m_svvVarDiffCoeff, m_useSpecVanVisc, and Nektar::NullNekDouble1DArray.

Referenced by v_SolveViscous().

◆ create()

static SolverUtils::EquationSystemSharedPtr Nektar::VelocityCorrectionScheme::create ( const LibUtilities::SessionReaderSharedPtr pSession,
const SpatialDomains::MeshGraphSharedPtr pGraph 
)
inlinestatic

Creates an instance of this class.

Definition at line 46 of file VelocityCorrectionScheme.h.

49  {
52  pGraph);
53  p->InitObject();
54  return p;
55  }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.

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

◆ EvaluateAdvection_SetPressureBCs()

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 101 of file VelocityCorrectionScheme.h.

104  {
105  v_EvaluateAdvection_SetPressureBCs(inarray, outarray, time);
106  }
virtual void v_EvaluateAdvection_SetPressureBCs(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)

References v_EvaluateAdvection_SetPressureBCs().

Referenced by v_InitObject().

◆ MeasureFlowrate()

NekDouble Nektar::VelocityCorrectionScheme::MeasureFlowrate ( const Array< OneD, Array< OneD, NekDouble >> &  inarray)
protected

Measure the volumetric flow rate through the volumetric flow rate reference surface.

This routine computes the volumetric flow rate

\[ Q(\mathbf{u}) = \frac{1}{\mu(R)} \int_R \mathbf{u} \cdot d\mathbf{s} \]

through the boundary region \( R \).

Definition at line 438 of file VelocityCorrectionScheme.cpp.

440 {
441  NekDouble flowrate = 0.0;
442 
443  if (m_flowrateBnd && m_flowrateBndID >= 0)
444  {
445  // If we're an actual boundary, calculate the vector flux through
446  // the boundary.
447  Array<OneD, Array<OneD, NekDouble>> boundary(m_spacedim);
448 
449  if (!m_homd1DFlowinPlane)
450  {
451  // General case
452  for (int i = 0; i < m_spacedim; ++i)
453  {
454  m_fields[i]->ExtractPhysToBnd(m_flowrateBndID, inarray[i],
455  boundary[i]);
456  }
457  flowrate = m_flowrateBnd->VectorFlux(boundary);
458  }
459  else if (m_planeID == 0)
460  {
461  // Homogeneous with forcing in plane. Calculate flux only on
462  // the meanmode - calculateFlux necessary for hybrid
463  // parallelisation.
464  for (int i = 0; i < m_spacedim; ++i)
465  {
466  m_fields[i]->GetPlane(m_planeID)->ExtractPhysToBnd(
467  m_flowrateBndID, inarray[i], boundary[i]);
468  }
469 
470  // the flowrate is calculated on the mean mode so it needs to be
471  // multiplied by LZ to be consistent with the general case.
472  flowrate = m_flowrateBnd->VectorFlux(boundary) *
473  m_session->GetParameter("LZ");
474  }
475  }
476  else if (m_flowrateBnd && !m_homd1DFlowinPlane)
477  {
478  // 3DH1D case with no Flowrate boundary defined: compute flux
479  // through the zero-th (mean) plane.
480  flowrate = m_flowrateBnd->Integral(inarray[2]);
481  }
482 
483  // Communication to obtain the total flowrate
485  {
486  m_comm->GetColumnComm()->AllReduce(flowrate, LibUtilities::ReduceSum);
487  }
488  else
489  {
490  m_comm->AllReduce(flowrate, LibUtilities::ReduceSum);
491  }
492  return flowrate / m_flowrateArea;
493 }
int m_spacedim
Spatial dimension (>= expansion dim).
LibUtilities::CommSharedPtr m_comm
Communicator.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
enum HomogeneousType m_HomogeneousType
MultiRegions::ExpListSharedPtr m_flowrateBnd
Flowrate reference surface.
NekDouble m_flowrateArea
Area of the boundary through which we are measuring the flowrate.
int m_planeID
Plane ID for cases with homogeneous expansion.
int m_flowrateBndID
Boundary ID of the flowrate reference surface.
double NekDouble

References Nektar::SolverUtils::EquationSystem::eHomogeneous1D, Nektar::SolverUtils::EquationSystem::m_comm, Nektar::SolverUtils::EquationSystem::m_fields, m_flowrateArea, m_flowrateBnd, m_flowrateBndID, m_homd1DFlowinPlane, Nektar::SolverUtils::EquationSystem::m_HomogeneousType, m_planeID, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_spacedim, and Nektar::LibUtilities::ReduceSum.

Referenced by SetupFlowrate(), and SolveUnsteadyStokesSystem().

◆ SetUpExtrapolation()

void Nektar::VelocityCorrectionScheme::SetUpExtrapolation ( void  )
protected

Definition at line 146 of file VelocityCorrectionScheme.cpp.

147 {
148  // creation of the extrapolation object
151  {
152  std::string vExtrapolation = v_GetExtrapolateStr();
153  if (m_session->DefinesSolverInfo("Extrapolation"))
154  {
155  vExtrapolation = v_GetSubSteppingExtrapolateStr(
156  m_session->GetSolverInfo("Extrapolation"));
157  }
159  vExtrapolation, m_session, m_fields, m_pressure, m_velocity,
160  m_advObject);
161 
162  m_extrapolation->SetForcing(m_forcing);
163  m_extrapolation->SubSteppingTimeIntegration(m_intScheme);
164  m_extrapolation->GenerateHOPBCMap(m_session);
165  }
166 }
MultiRegions::ExpListSharedPtr m_pressure
Pointer to field holding pressure field.
ExtrapolateSharedPtr m_extrapolation
Array< OneD, int > m_velocity
int which identifies which components of m_fields contains the velocity (u,v,w);
EquationType m_equationType
equation type;
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Forcing terms.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
Wrapper to the time integration scheme.
virtual std::string v_GetExtrapolateStr(void)
virtual std::string v_GetSubSteppingExtrapolateStr(const std::string &instr)
@ eUnsteadyStokes
@ eUnsteadyNavierStokes
ExtrapolateFactory & GetExtrapolateFactory()
Definition: Extrapolate.cpp:48

References Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::eUnsteadyNavierStokes, Nektar::eUnsteadyStokes, Nektar::GetExtrapolateFactory(), Nektar::SolverUtils::AdvectionSystem::m_advObject, Nektar::IncNavierStokes::m_equationType, Nektar::IncNavierStokes::m_extrapolation, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::IncNavierStokes::m_forcing, Nektar::SolverUtils::UnsteadySystem::m_intScheme, Nektar::IncNavierStokes::m_pressure, Nektar::SolverUtils::EquationSystem::m_session, Nektar::IncNavierStokes::m_velocity, v_GetExtrapolateStr(), and v_GetSubSteppingExtrapolateStr().

Referenced by v_InitObject().

◆ SetupFlowrate()

void Nektar::VelocityCorrectionScheme::SetupFlowrate ( NekDouble  aii_dt)
protected

Set up the Stokes solution used to impose constant flowrate through a boundary.

This routine solves a Stokes equation using a unit forcing direction, specified by the user to be in the desired flow direction. This field can then be used to correct the end of each timestep to impose a constant volumetric flow rate through a user-defined boundary.

There are three modes of operation:

  • Standard two-dimensional or three-dimensional simulations (e.g. pipes or channels)
  • 3DH1D simulations where the forcing is not in the homogeneous direction (e.g. channel flow, where the y-direction of the 2D mesh is perpendicular to the wall);
  • 3DH1D simulations where the forcing is in the homogeneous direction (e.g. pipe flow in the z-direction).

In the first two cases, the user should define:

  • the Flowrate parameter, which dictates the volumetric flux through the reference area
  • tag a boundary region with the Flowrate user-defined type to define the reference area
  • define a FlowrateForce function with components ForceX, ForceY and ForceZ that defines a unit forcing in the appropriate direction.

In the latter case, the user should define only the Flowrate; the reference area is taken to be the homogeneous plane and the force is assumed to be the unit z-vector \( \hat{e}_z \).

This routine solves a single timestep of the Stokes problem (premultiplied by the backwards difference coefficient):

\[ \frac{\partial\mathbf{u}}{\partial t} = -\nabla p + \nu\nabla^2\mathbf{u} + \mathbf{f} \]

with a zero initial condition to obtain a field \( \mathbf{u}_s \). The flowrate is then corrected at each timestep \( n \) by adding the correction \( \alpha\mathbf{u}_s \) where

\[ \alpha = \frac{\overline{Q} - Q(\mathbf{u^n})}{Q(\mathbf{u}_s)} \]

where \( Q(\cdot)\) is the volumetric flux through the appropriate surface or line, which is implemented in VelocityCorrectionScheme::MeasureFlowrate. For more details, see chapter 3.2 of the thesis of D. Moxey (University of Warwick, 2011).

Definition at line 216 of file VelocityCorrectionScheme.cpp.

217 {
218  m_flowrateBndID = -1;
219  m_flowrateArea = 0.0;
220 
221  const Array<OneD, const SpatialDomains::BoundaryConditionShPtr> &bcs =
222  m_fields[0]->GetBndConditions();
223 
224  std::string forces[] = {"X", "Y", "Z"};
225  Array<OneD, NekDouble> flowrateForce(m_spacedim, 0.0);
226 
227  // Set up flowrate forces.
228  bool defined = true;
229  for (int i = 0; i < m_spacedim; ++i)
230  {
231  std::string varName = std::string("Force") + forces[i];
232  defined = m_session->DefinesFunction("FlowrateForce", varName);
233 
234  if (!defined && m_HomogeneousType == eHomogeneous1D)
235  {
236  break;
237  }
238 
239  ASSERTL0(defined,
240  "A 'FlowrateForce' function must defined with components "
241  "[ForceX, ...] to define direction of flowrate forcing");
242 
244  m_session->GetFunction("FlowrateForce", varName);
245  flowrateForce[i] = ffunc->Evaluate();
246  }
247 
248  // Define flag for case with homogeneous expansion and forcing not in the
249  // z-direction
250  m_homd1DFlowinPlane = false;
251  if (defined && m_HomogeneousType == eHomogeneous1D)
252  {
253  m_homd1DFlowinPlane = true;
254  }
255 
256  // For 3DH1D simulations, if force isn't defined then assume in
257  // z-direction.
258  if (!defined)
259  {
260  flowrateForce[2] = 1.0;
261  }
262 
263  // Find the boundary condition that is tagged as the flowrate boundary.
264  for (int i = 0; i < bcs.size(); ++i)
265  {
266  if (boost::iequals(bcs[i]->GetUserDefined(), "Flowrate"))
267  {
268  m_flowrateBndID = i;
269  break;
270  }
271  }
272 
273  int tmpBr = m_flowrateBndID;
274  m_comm->AllReduce(tmpBr, LibUtilities::ReduceMax);
275  ASSERTL0(tmpBr >= 0 || m_HomogeneousType == eHomogeneous1D,
276  "One boundary region must be marked using the 'Flowrate' "
277  "user-defined type to monitor the volumetric flowrate.");
278 
279  // Extract an appropriate expansion list to represents the boundary.
280  if (m_flowrateBndID >= 0)
281  {
282  // For a boundary, extract the boundary itself.
283  m_flowrateBnd = m_fields[0]->GetBndCondExpansions()[m_flowrateBndID];
284  }
286  {
287  // For 3DH1D simulations with no force specified, find the mean
288  // (0th) plane.
289  Array<OneD, unsigned int> zIDs = m_fields[0]->GetZIDs();
290  int tmpId = -1;
291 
292  for (int i = 0; i < zIDs.size(); ++i)
293  {
294  if (zIDs[i] == 0)
295  {
296  tmpId = i;
297  break;
298  }
299  }
300 
301  ASSERTL1(tmpId <= 0, "Should be either at location 0 or -1 if not "
302  "found");
303 
304  if (tmpId != -1)
305  {
306  m_flowrateBnd = m_fields[0]->GetPlane(tmpId);
307  }
308  }
309 
310  // At this point, some processors may not have m_flowrateBnd
311  // set if they don't contain the appropriate boundary. To
312  // calculate the area, we integrate 1.0 over the boundary
313  // (which has been set up with the appropriate subcommunicator
314  // to avoid deadlock), and then communicate this to the other
315  // processors with an AllReduce.
316  if (m_flowrateBnd)
317  {
318  Array<OneD, NekDouble> inArea(m_flowrateBnd->GetNpoints(), 1.0);
319  m_flowrateArea = m_flowrateBnd->Integral(inArea);
320  }
322 
323  // In homogeneous case with forcing not aligned to the z-direction,
324  // redefine m_flowrateBnd so it is a 1D expansion
327  {
328  // For 3DH1D simulations with no force specified, find the mean
329  // (0th) plane.
330  Array<OneD, unsigned int> zIDs = m_fields[0]->GetZIDs();
331  m_planeID = -1;
332 
333  for (int i = 0; i < zIDs.size(); ++i)
334  {
335  if (zIDs[i] == 0)
336  {
337  m_planeID = i;
338  break;
339  }
340  }
341 
342  ASSERTL1(m_planeID <= 0, "Should be either at location 0 or -1 "
343  "if not found");
344 
345  if (m_planeID != -1)
346  {
347  m_flowrateBnd =
348  m_fields[0]->GetBndCondExpansions()[m_flowrateBndID]->GetPlane(
349  m_planeID);
350  }
351  }
352 
353  // Set up some storage for the Stokes solution (to be stored in
354  // m_flowrateStokes) and its initial condition (inTmp), which holds the
355  // unit forcing.
356  int nqTot = m_fields[0]->GetNpoints();
357  Array<OneD, Array<OneD, NekDouble>> inTmp(m_spacedim);
358  m_flowrateStokes = Array<OneD, Array<OneD, NekDouble>>(m_spacedim);
359 
360  for (int i = 0; i < m_spacedim; ++i)
361  {
362  inTmp[i] = Array<OneD, NekDouble>(nqTot, flowrateForce[i] * aii_dt);
363  m_flowrateStokes[i] = Array<OneD, NekDouble>(nqTot, 0.0);
364 
366  {
367  Array<OneD, NekDouble> inTmp2(nqTot);
368  m_fields[i]->HomogeneousFwdTrans(inTmp[i], inTmp2);
369  m_fields[i]->SetWaveSpace(true);
370  inTmp[i] = inTmp2;
371  }
372 
373  Vmath::Zero(m_fields[i]->GetNcoeffs(), m_fields[i]->UpdateCoeffs(), 1);
374  }
375 
376  // Create temporary extrapolation object to avoid issues with
377  // m_extrapolation for HOPBCs using higher order timestepping schemes.
378  // Zero pressure BCs in Neumann boundaries that may have been
379  // set in the advection step.
380  Array<OneD, const SpatialDomains::BoundaryConditionShPtr> PBndConds =
381  m_pressure->GetBndConditions();
382  Array<OneD, MultiRegions::ExpListSharedPtr> PBndExp =
383  m_pressure->GetBndCondExpansions();
384  for (int n = 0; n < PBndConds.size(); ++n)
385  {
386  if (PBndConds[n]->GetBoundaryConditionType() ==
388  {
389  Vmath::Zero(PBndExp[n]->GetNcoeffs(), PBndExp[n]->UpdateCoeffs(),
390  1);
391  }
392  }
393 
394  // Finally, calculate the solution and the flux of the Stokes
395  // solution. We set m_greenFlux to maximum numeric limit, which signals
396  // to SolveUnsteadyStokesSystem that we don't need to apply a flowrate
397  // force.
398  m_greenFlux = numeric_limits<NekDouble>::max();
399  m_flowrateAiidt = aii_dt;
400 
401  // Save the number of convective field in case it is not set
402  // to spacedim. Only need velocity components for stokes forcing
403  int SaveNConvectiveFields = m_nConvectiveFields;
405  SolveUnsteadyStokesSystem(inTmp, m_flowrateStokes, 0.0, aii_dt);
406  m_nConvectiveFields = SaveNConvectiveFields;
408 
409  // If the user specified IO_FlowSteps, open a handle to store output.
410  if (m_comm->GetRank() == 0 && m_flowrateSteps &&
411  !m_flowrateStream.is_open())
412  {
413  std::string filename = m_session->GetSessionName();
414  filename += ".prs";
415  m_flowrateStream.open(filename.c_str());
416  m_flowrateStream.setf(ios::scientific, ios::floatfield);
417  m_flowrateStream << "# step time dP" << endl
418  << "# -------------------------------------------"
419  << endl;
420  }
421 
422  // Replace pressure BCs with those evaluated from advection step
423  m_extrapolation->CopyPressureHBCsToPbndExp();
424 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
int m_nConvectiveFields
Number of fields to be convected;.
SOLVER_UTILS_EXPORT int GetNcoeffs()
NekDouble m_greenFlux
Flux of the Stokes function solution.
NekDouble MeasureFlowrate(const Array< OneD, Array< OneD, NekDouble >> &inarray)
Measure the volumetric flow rate through the volumetric flow rate reference surface.
void SolveUnsteadyStokesSystem(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time, const NekDouble a_iixDt)
Array< OneD, Array< OneD, NekDouble > > m_flowrateStokes
Stokes solution used to impose flowrate.
int m_flowrateSteps
Interval at which to record flowrate data.
std::ofstream m_flowrateStream
Output stream to record flowrate.
NekDouble m_flowrateAiidt
Value of aii_dt used to compute Stokes flowrate solution.
std::shared_ptr< Equation > EquationSharedPtr
Definition: Equation.h:130
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492

References ASSERTL0, ASSERTL1, Nektar::SolverUtils::EquationSystem::eHomogeneous1D, Nektar::SpatialDomains::eNeumann, Nektar::SolverUtils::EquationSystem::GetNcoeffs(), Nektar::SolverUtils::EquationSystem::m_comm, Nektar::IncNavierStokes::m_extrapolation, Nektar::SolverUtils::EquationSystem::m_fields, m_flowrateAiidt, m_flowrateArea, m_flowrateBnd, m_flowrateBndID, m_flowrateSteps, m_flowrateStokes, m_flowrateStream, m_greenFlux, m_homd1DFlowinPlane, Nektar::SolverUtils::EquationSystem::m_HomogeneousType, Nektar::IncNavierStokes::m_nConvectiveFields, m_planeID, Nektar::IncNavierStokes::m_pressure, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_spacedim, MeasureFlowrate(), Nektar::LibUtilities::ReduceMax, SolveUnsteadyStokesSystem(), and Vmath::Zero().

Referenced by SolveUnsteadyStokesSystem().

◆ SetUpPressureForcing()

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 69 of file VelocityCorrectionScheme.h.

72  {
73  v_SetUpPressureForcing(fields, Forcing, aii_Dt);
74  }
virtual void v_SetUpPressureForcing(const Array< OneD, const Array< OneD, NekDouble >> &fields, Array< OneD, Array< OneD, NekDouble >> &Forcing, const NekDouble aii_Dt)

References v_SetUpPressureForcing().

Referenced by SolveUnsteadyStokesSystem().

◆ SetUpSVV()

void Nektar::VelocityCorrectionScheme::SetUpSVV ( void  )
protected

Definition at line 929 of file VelocityCorrectionScheme.cpp.

930 {
931 
932  m_session->MatchSolverInfo("SpectralVanishingViscosity", "PowerKernel",
933  m_useSpecVanVisc, false);
934 
935  if (m_useSpecVanVisc)
936  {
937  m_useHomo1DSpecVanVisc = true;
938  }
939  else
940  {
941  m_session->MatchSolverInfo("SpectralVanishingViscositySpectralHP",
942  "PowerKernel", m_useSpecVanVisc, false);
943  }
944 
945  if (m_useSpecVanVisc)
946  {
947  m_IsSVVPowerKernel = true;
948  }
949  else
950  {
951  m_session->MatchSolverInfo("SpectralVanishingViscosity", "DGKernel",
952  m_useSpecVanVisc, false);
953  if (m_useSpecVanVisc)
954  {
955  m_useHomo1DSpecVanVisc = true;
956  }
957  else
958  {
959  m_session->MatchSolverInfo("SpectralVanishingViscositySpectralHP",
960  "DGKernel", m_useSpecVanVisc, false);
961  }
962 
963  if (m_useSpecVanVisc)
964  {
965  m_IsSVVPowerKernel = false;
966  }
967  }
968 
969  // set up varcoeff kernel if PowerKernel or DG is specified
970  if (m_useSpecVanVisc)
971  {
972  Array<OneD, Array<OneD, NekDouble>> SVVVelFields =
974  if (m_session->DefinesFunction("SVVVelocityMagnitude"))
975  {
976  if (m_comm->GetRank() == 0)
977  {
978  cout << "Seting up SVV velocity from "
979  "SVVVelocityMagnitude section in session file"
980  << endl;
981  }
982  int nvel = m_velocity.size();
983  int phystot = m_fields[0]->GetTotPoints();
984  SVVVelFields = Array<OneD, Array<OneD, NekDouble>>(nvel);
985  vector<string> vars;
986  for (int i = 0; i < nvel; ++i)
987  {
988  SVVVelFields[i] = Array<OneD, NekDouble>(phystot);
989  vars.push_back(m_session->GetVariable(m_velocity[i]));
990  }
991 
992  // Load up files into m_fields;
993  GetFunction("SVVVelocityMagnitude")->Evaluate(vars, SVVVelFields);
994  }
995 
996  m_svvVarDiffCoeff = Array<OneD, NekDouble>(m_fields[0]->GetNumElmts());
997  SVVVarDiffCoeff(1.0, m_svvVarDiffCoeff, SVVVelFields);
998  m_session->LoadParameter("SVVDiffCoeff", m_sVVDiffCoeff, 1.0);
999  }
1000  else
1001  {
1003  m_session->LoadParameter("SVVDiffCoeff", m_sVVDiffCoeff, 0.1);
1004  }
1005 
1006  // Load parameters for Spectral Vanishing Viscosity
1007  if (m_useSpecVanVisc == false)
1008  {
1009  m_session->MatchSolverInfo("SpectralVanishingViscosity", "True",
1010  m_useSpecVanVisc, false);
1011  if (m_useSpecVanVisc == false)
1012  {
1013  m_session->MatchSolverInfo("SpectralVanishingViscosity",
1014  "ExpKernel", m_useSpecVanVisc, false);
1015  }
1017 
1018  if (m_useSpecVanVisc == false)
1019  {
1020  m_session->MatchSolverInfo("SpectralVanishingViscositySpectralHP",
1021  "True", m_useSpecVanVisc, false);
1022  if (m_useSpecVanVisc == false)
1023  {
1024  m_session->MatchSolverInfo(
1025  "SpectralVanishingViscositySpectralHP", "ExpKernel",
1026  m_useSpecVanVisc, false);
1027  }
1028  }
1029  }
1030 
1031  // Case of only Homo1D kernel
1032  if (m_session->DefinesSolverInfo("SpectralVanishingViscosityHomo1D"))
1033  {
1034  m_session->MatchSolverInfo("SpectralVanishingViscosityHomo1D", "True",
1035  m_useHomo1DSpecVanVisc, false);
1036  if (m_useHomo1DSpecVanVisc == false)
1037  {
1038  m_session->MatchSolverInfo("SpectralVanishingViscosityHomo1D",
1039  "ExpKernel", m_useHomo1DSpecVanVisc,
1040  false);
1041  }
1042  }
1043 
1044  m_session->LoadParameter("SVVCutoffRatio", m_sVVCutoffRatio, 0.75);
1045  m_session->LoadParameter("SVVCutoffRatioHomo1D", m_sVVCutoffRatioHomo1D,
1047  m_session->LoadParameter("SVVDiffCoeffHomo1D", m_sVVDiffCoeffHomo1D,
1048  m_sVVDiffCoeff);
1049 
1051  {
1052  ASSERTL0(
1053  m_nConvectiveFields > 2,
1054  "Expect to have three velocity fields with homogenous expansion");
1055 
1057  {
1058  Array<OneD, unsigned int> planes;
1059  planes = m_fields[0]->GetZIDs();
1060 
1061  int num_planes = planes.size();
1062  Array<OneD, NekDouble> SVV(num_planes, 0.0);
1063  NekDouble fac;
1064  int kmodes = m_fields[0]->GetHomogeneousBasis()->GetNumModes();
1065  int pstart;
1066 
1067  pstart = m_sVVCutoffRatioHomo1D * kmodes;
1068 
1069  for (int n = 0; n < num_planes; ++n)
1070  {
1071  if (planes[n] > pstart)
1072  {
1073  fac = (NekDouble)((planes[n] - kmodes) *
1074  (planes[n] - kmodes)) /
1075  ((NekDouble)((planes[n] - pstart) *
1076  (planes[n] - pstart)));
1077  SVV[n] = m_sVVDiffCoeffHomo1D * exp(-fac) / m_kinvis;
1078  }
1079  }
1080 
1081  for (int i = 0; i < m_velocity.size(); ++i)
1082  {
1083  m_fields[m_velocity[i]]->SetHomo1DSpecVanVisc(SVV);
1084  }
1085  }
1086  }
1087 }
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
Get a SessionFunction by name.
NekDouble m_sVVDiffCoeffHomo1D
Diffusion coefficient of SVV modes in homogeneous 1D Direction.
void SVVVarDiffCoeff(const NekDouble velmag, Array< OneD, NekDouble > &diffcoeff, const Array< OneD, Array< OneD, NekDouble >> &vel=NullNekDoubleArrayOfArray)
bool m_useHomo1DSpecVanVisc
bool to identify if spectral vanishing viscosity is active.
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray

References ASSERTL0, Nektar::SolverUtils::EquationSystem::eHomogeneous1D, Nektar::SolverUtils::EquationSystem::GetFunction(), Nektar::SolverUtils::EquationSystem::m_comm, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::SolverUtils::EquationSystem::m_HomogeneousType, m_IsSVVPowerKernel, Nektar::IncNavierStokes::m_kinvis, Nektar::IncNavierStokes::m_nConvectiveFields, Nektar::SolverUtils::EquationSystem::m_session, m_sVVCutoffRatio, m_sVVCutoffRatioHomo1D, m_sVVDiffCoeff, m_sVVDiffCoeffHomo1D, m_svvVarDiffCoeff, m_useHomo1DSpecVanVisc, m_useSpecVanVisc, Nektar::IncNavierStokes::m_velocity, Nektar::NullNekDouble1DArray, Nektar::NullNekDoubleArrayOfArray, and SVVVarDiffCoeff().

Referenced by v_InitObject().

◆ SetUpViscousForcing()

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 76 of file VelocityCorrectionScheme.h.

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)

References v_SetUpViscousForcing().

Referenced by SolveUnsteadyStokesSystem().

◆ SolvePressure()

void Nektar::VelocityCorrectionScheme::SolvePressure ( const Array< OneD, NekDouble > &  Forcing)
inline

Definition at line 83 of file VelocityCorrectionScheme.h.

84  {
85  v_SolvePressure(Forcing);
86  }
virtual void v_SolvePressure(const Array< OneD, NekDouble > &Forcing)

References v_SolvePressure().

Referenced by SolveUnsteadyStokesSystem().

◆ 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 726 of file VelocityCorrectionScheme.cpp.

730 {
731  // Set up flowrate if we're starting for the first time or the value of
732  // aii_Dt has changed.
733  if (m_flowrate > 0.0 && (aii_Dt != m_flowrateAiidt))
734  {
735  SetupFlowrate(aii_Dt);
736  }
737 
738  int physTot = m_fields[0]->GetTotPoints();
739 
740  // Substep the pressure boundary condition if using substepping
741  m_extrapolation->SubStepSetPressureBCs(inarray, aii_Dt, m_kinvis);
742 
743  // Set up forcing term for pressure Poisson equation
744  LibUtilities::Timer timer;
745  timer.Start();
746  SetUpPressureForcing(inarray, m_F, aii_Dt);
747  timer.Stop();
748  timer.AccumulateRegion("Pressure Forcing");
749 
750  // Solve Pressure System
751  timer.Start();
752  SolvePressure(m_F[0]);
753  timer.Stop();
754  timer.AccumulateRegion("Pressure Solve");
755 
756  // Set up forcing term for Helmholtz problems
757  timer.Start();
758  SetUpViscousForcing(inarray, m_F, aii_Dt);
759  timer.Stop();
760  timer.AccumulateRegion("Viscous Forcing");
761 
762  // Solve velocity system
763  timer.Start();
764  SolveViscous(m_F, inarray, outarray, aii_Dt);
765  timer.Stop();
766  timer.AccumulateRegion("Viscous Solve");
767 
768  // Apply flowrate correction
769  if (m_flowrate > 0.0 && m_greenFlux != numeric_limits<NekDouble>::max())
770  {
771  NekDouble currentFlux = MeasureFlowrate(outarray);
772  m_alpha = (m_flowrate - currentFlux) / m_greenFlux;
773 
774  for (int i = 0; i < m_spacedim; ++i)
775  {
776  Vmath::Svtvp(physTot, m_alpha, m_flowrateStokes[i], 1, outarray[i],
777  1, outarray[i], 1);
778  // Enusre coeff space is updated for next time step
779  m_fields[i]->FwdTransLocalElmt(outarray[i],
780  m_fields[i]->UpdateCoeffs());
781  // Impsoe symmetry of flow on coeff space (good to enfore
782  // periodicity).
783  m_fields[i]->LocalToGlobal();
784  m_fields[i]->GlobalToLocal();
785  }
786  }
787 }
void SetUpPressureForcing(const Array< OneD, const Array< OneD, NekDouble >> &fields, Array< OneD, Array< OneD, NekDouble >> &Forcing, const NekDouble aii_Dt)
NekDouble m_flowrate
Desired volumetric flowrate.
NekDouble m_alpha
Current flowrate correction.
void SolveViscous(const Array< OneD, const Array< OneD, NekDouble >> &Forcing, const Array< OneD, const Array< OneD, NekDouble >> &inarray, 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)
Array< OneD, Array< OneD, NekDouble > > m_F
void SetupFlowrate(NekDouble aii_dt)
Set up the Stokes solution used to impose constant flowrate through a boundary.
void SolvePressure(const Array< OneD, NekDouble > &Forcing)
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
svtvp (scalar times vector plus vector): z = alpha*x + y
Definition: Vmath.cpp:622

References Nektar::LibUtilities::Timer::AccumulateRegion(), m_alpha, Nektar::IncNavierStokes::m_extrapolation, m_F, Nektar::SolverUtils::EquationSystem::m_fields, m_flowrate, m_flowrateAiidt, m_flowrateStokes, m_greenFlux, Nektar::IncNavierStokes::m_kinvis, Nektar::SolverUtils::EquationSystem::m_spacedim, MeasureFlowrate(), SetupFlowrate(), SetUpPressureForcing(), SetUpViscousForcing(), SolvePressure(), SolveViscous(), Nektar::LibUtilities::Timer::Start(), Nektar::LibUtilities::Timer::Stop(), and Vmath::Svtvp().

Referenced by SetupFlowrate(), v_InitObject(), and Nektar::SmoothedProfileMethod::v_SolveUnsteadyStokesSystem().

◆ SolveViscous()

void Nektar::VelocityCorrectionScheme::SolveViscous ( const Array< OneD, const Array< OneD, NekDouble >> &  Forcing,
const Array< OneD, const Array< OneD, NekDouble >> &  inarray,
Array< OneD, Array< OneD, NekDouble >> &  outarray,
const NekDouble  aii_Dt 
)
inline

Definition at line 88 of file VelocityCorrectionScheme.h.

92  {
93  v_SolveViscous(Forcing, inarray, outarray, aii_Dt);
94  }
virtual void v_SolveViscous(const Array< OneD, const Array< OneD, NekDouble >> &Forcing, const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble aii_Dt)

References v_SolveViscous().

Referenced by SolveUnsteadyStokesSystem().

◆ SVVVarDiffCoeff()

void Nektar::VelocityCorrectionScheme::SVVVarDiffCoeff ( const NekDouble  velmag,
Array< OneD, NekDouble > &  diffcoeff,
const Array< OneD, Array< OneD, NekDouble >> &  vel = NullNekDoubleArrayOfArray 
)
protected

Definition at line 1089 of file VelocityCorrectionScheme.cpp.

1092 {
1093  int phystot = m_fields[0]->GetTotPoints();
1094  int nel = m_fields[0]->GetNumElmts();
1095  int nvel, cnt;
1096 
1097  Array<OneD, NekDouble> tmp;
1098 
1099  Vmath::Fill(nel, velmag, diffcoeff, 1);
1100 
1101  if (vel != NullNekDoubleArrayOfArray)
1102  {
1103  Array<OneD, NekDouble> Velmag(phystot);
1104  nvel = vel.size();
1105  // calculate magnitude of v
1106  Vmath::Vmul(phystot, vel[0], 1, vel[0], 1, Velmag, 1);
1107  for (int n = 1; n < nvel; ++n)
1108  {
1109  Vmath::Vvtvp(phystot, vel[n], 1, vel[n], 1, Velmag, 1, Velmag, 1);
1110  }
1111  Vmath::Vsqrt(phystot, Velmag, 1, Velmag, 1);
1112 
1113  cnt = 0;
1114  Array<OneD, NekDouble> tmp;
1115  // calculate mean value of vel mag.
1116  for (int i = 0; i < nel; ++i)
1117  {
1118  int nq = m_fields[0]->GetExp(i)->GetTotPoints();
1119  tmp = Velmag + cnt;
1120  diffcoeff[i] = m_fields[0]->GetExp(i)->Integral(tmp);
1121  Vmath::Fill(nq, 1.0, tmp, 1);
1122  NekDouble area = m_fields[0]->GetExp(i)->Integral(tmp);
1123  diffcoeff[i] = diffcoeff[i] / area;
1124  cnt += nq;
1125  }
1126  }
1127  else
1128  {
1129  nvel = m_expdim;
1130  }
1131 
1132  for (int e = 0; e < nel; e++)
1133  {
1134  LocalRegions::ExpansionSharedPtr exp = m_fields[0]->GetExp(e);
1135  NekDouble h = 0;
1136 
1137  // Find maximum length of edge.
1138  for (int i = 0; i < exp->GetGeom()->GetNumEdges(); ++i)
1139  {
1140  h = max(h, exp->GetGeom()->GetEdge(i)->GetVertex(0)->dist(
1141  *(exp->GetGeom()->GetEdge(i)->GetVertex(1))));
1142  }
1143 
1144  int p = 0;
1145  for (int i = 0; i < m_expdim; ++i)
1146  {
1147  p = max(p, exp->GetBasisNumModes(i) - 1);
1148  }
1149 
1150  diffcoeff[e] *= h / p;
1151  }
1152 }
int m_expdim
Expansion dimension.
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
void Vsqrt(int n, const T *x, const int incx, T *y, const int incy)
sqrt y = sqrt(x)
Definition: Vmath.cpp:534
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:209
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:574
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.cpp:45

References Vmath::Fill(), Nektar::SolverUtils::EquationSystem::m_expdim, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::NullNekDoubleArrayOfArray, CellMLToNektar.cellml_metadata::p, Vmath::Vmul(), Vmath::Vsqrt(), and Vmath::Vvtvp().

Referenced by SetUpSVV().

◆ v_DoInitialise()

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 606 of file VelocityCorrectionScheme.cpp.

607 {
608  m_F = Array<OneD, Array<OneD, NekDouble>>(m_nConvectiveFields);
609 
610  for (int i = 0; i < m_nConvectiveFields; ++i)
611  {
612  m_F[i] = Array<OneD, NekDouble>(m_fields[0]->GetTotPoints(), 0.0);
613  }
614 
615  m_flowrateAiidt = 0.0;
616 
618 
619  // Set up Field Meta Data for output files
620  m_fieldMetaDataMap["Kinvis"] = boost::lexical_cast<std::string>(m_kinvis);
621  m_fieldMetaDataMap["TimeStep"] =
622  boost::lexical_cast<std::string>(m_timestep);
623 
624  // set boundary conditions here so that any normal component
625  // correction are imposed before they are imposed on initial
626  // field below
628 
629  // Ensure the initial conditions have correct BCs
630  for (int i = 0; i < m_fields.size(); ++i)
631  {
632  m_fields[i]->ImposeDirichletConditions(m_fields[i]->UpdateCoeffs());
633  m_fields[i]->LocalToGlobal();
634  m_fields[i]->GlobalToLocal();
635  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
636  m_fields[i]->UpdatePhys());
637  }
638 }
void SetBoundaryConditions(NekDouble time)
time dependent boundary conditions updating
NekDouble m_timestep
Time step size.
NekDouble m_time
Current time of simulation.
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
Map to identify relevant solver info to dump in output fields.
SOLVER_UTILS_EXPORT int GetTotPoints()
virtual SOLVER_UTILS_EXPORT void v_DoInitialise()
Sets up initial conditions.

References Nektar::SolverUtils::EquationSystem::GetTotPoints(), m_F, Nektar::SolverUtils::EquationSystem::m_fieldMetaDataMap, Nektar::SolverUtils::EquationSystem::m_fields, m_flowrateAiidt, 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().

◆ v_EvaluateAdvection_SetPressureBCs()

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 691 of file VelocityCorrectionScheme.cpp.

694 {
695  LibUtilities::Timer timer;
696  timer.Start();
697  EvaluateAdvectionTerms(inarray, outarray, time);
698  timer.Stop();
699  timer.AccumulateRegion("Advection Terms");
700 
701  // Smooth advection
702  if (m_SmoothAdvection)
703  {
704  for (int i = 0; i < m_nConvectiveFields; ++i)
705  {
706  m_pressure->SmoothField(outarray[i]);
707  }
708  }
709 
710  // Add forcing terms
711  for (auto &x : m_forcing)
712  {
713  x->Apply(m_fields, inarray, outarray, time);
714  }
715 
716  // Calculate High-Order pressure boundary conditions
717  timer.Start();
718  m_extrapolation->EvaluatePressureBCs(inarray, outarray, m_kinvis);
719  timer.Stop();
720  timer.AccumulateRegion("Pressure BCs");
721 }
bool m_SmoothAdvection
bool to identify if advection term smoothing is requested
void EvaluateAdvectionTerms(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)

References Nektar::LibUtilities::Timer::AccumulateRegion(), 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, Nektar::IncNavierStokes::m_SmoothAdvection, Nektar::LibUtilities::Timer::Start(), and Nektar::LibUtilities::Timer::Stop().

Referenced by EvaluateAdvection_SetPressureBCs().

◆ v_GenerateSummary()

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, and Nektar::SmoothedProfileMethod.

Definition at line 519 of file VelocityCorrectionScheme.cpp.

520 {
522  SolverUtils::AddSummaryItem(s, "Splitting Scheme",
523  "Velocity correction (strong press. form)");
524 
525  if (m_extrapolation->GetSubStepName().size())
526  {
527  SolverUtils::AddSummaryItem(s, "Substepping",
528  m_extrapolation->GetSubStepName());
529  }
530 
531  string dealias = m_homogen_dealiasing ? "Homogeneous1D" : "";
533  {
534  dealias += (dealias == "" ? "" : " + ") + string("spectral/hp");
535  }
536  if (dealias != "")
537  {
538  SolverUtils::AddSummaryItem(s, "Dealiasing", dealias);
539  }
540 
541  string smoothing = m_useSpecVanVisc ? "spectral/hp" : "";
542  if (smoothing != "")
543  {
545  {
547  s, "Smoothing-SpecHP",
548  "SVV (" + smoothing + " Exp Kernel(cut-off = " +
549  boost::lexical_cast<string>(m_sVVCutoffRatio) +
550  ", diff coeff = " +
551  boost::lexical_cast<string>(m_sVVDiffCoeff) + "))");
552  }
553  else
554  {
555  if (m_IsSVVPowerKernel)
556  {
558  s, "Smoothing-SpecHP",
559  "SVV (" + smoothing + " Power Kernel (Power ratio =" +
560  boost::lexical_cast<string>(m_sVVCutoffRatio) +
561  ", diff coeff = " +
562  boost::lexical_cast<string>(m_sVVDiffCoeff) +
563  "*Uh/p))");
564  }
565  else
566  {
568  s, "Smoothing-SpecHP",
569  "SVV (" + smoothing + " DG Kernel (diff coeff = " +
570  boost::lexical_cast<string>(m_sVVDiffCoeff) +
571  "*Uh/p))");
572  }
573  }
574  }
575 
577  {
579  s, "Smoothing-Homo1D",
580  "SVV (Homogeneous1D - Exp Kernel(cut-off = " +
581  boost::lexical_cast<string>(m_sVVCutoffRatioHomo1D) +
582  ", diff coeff = " +
583  boost::lexical_cast<string>(m_sVVDiffCoeffHomo1D) + "))");
584  }
585 
587  {
589  s, "GJP Stab. Impl. ",
590  m_session->GetSolverInfo("GJPStabilisation"));
591  SolverUtils::AddSummaryItem(s, "GJP Stab. JumpScale", m_GJPJumpScale);
592 
593  if (boost::iequals(m_session->GetSolverInfo("GJPStabilisation"),
594  "Explicit"))
595  {
597  s, "GJP Normal Velocity",
598  m_session->GetSolverInfo("GJPNormalVelocity"));
599  }
600  }
601 }
bool m_specHP_dealiasing
Flag to determine if dealisising is usde for the Spectral/hp element discretisation.
bool m_homogen_dealiasing
Flag to determine if dealiasing is used for homogeneous simulations.
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s)
Print a summary of time stepping parameters.
bool m_useGJPStabilisation
bool to identify if GJP semi-implicit is active.
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:49

References Nektar::SolverUtils::AddSummaryItem(), Nektar::SolverUtils::EquationSystem::eHomogeneous1D, Nektar::IncNavierStokes::m_extrapolation, m_GJPJumpScale, Nektar::SolverUtils::EquationSystem::m_homogen_dealiasing, Nektar::SolverUtils::EquationSystem::m_HomogeneousType, m_IsSVVPowerKernel, Nektar::SolverUtils::EquationSystem::m_session, Nektar::SolverUtils::EquationSystem::m_specHP_dealiasing, m_sVVCutoffRatio, m_sVVCutoffRatioHomo1D, m_sVVDiffCoeff, m_sVVDiffCoeffHomo1D, m_svvVarDiffCoeff, m_useGJPStabilisation, m_useHomo1DSpecVanVisc, m_useSpecVanVisc, Nektar::NullNekDouble1DArray, and Nektar::SolverUtils::UnsteadySystem::v_GenerateSummary().

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

◆ v_GetExtrapolateStr()

virtual std::string Nektar::VelocityCorrectionScheme::v_GetExtrapolateStr ( void  )
inlineprotectedvirtual

Reimplemented in Nektar::VCSWeakPressure.

Definition at line 206 of file VelocityCorrectionScheme.h.

207  {
208  return "Standard";
209  }

Referenced by SetUpExtrapolation().

◆ v_GetForceDimension()

int Nektar::VelocityCorrectionScheme::v_GetForceDimension ( void  )
protectedvirtual

Implements Nektar::IncNavierStokes.

Definition at line 683 of file VelocityCorrectionScheme.cpp.

684 {
685  return m_session->GetVariables().size() - 1;
686 }

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

◆ v_GetSubSteppingExtrapolateStr()

virtual std::string Nektar::VelocityCorrectionScheme::v_GetSubSteppingExtrapolateStr ( const std::string &  instr)
inlineprotectedvirtual

Reimplemented in Nektar::VCSWeakPressure.

Definition at line 211 of file VelocityCorrectionScheme.h.

212  {
213  return instr;
214  }

Referenced by SetUpExtrapolation().

◆ v_GetSystemSingularChecks()

Array< OneD, bool > Nektar::VelocityCorrectionScheme::v_GetSystemSingularChecks ( )
protectedvirtual

Reimplemented from Nektar::SolverUtils::EquationSystem.

Definition at line 672 of file VelocityCorrectionScheme.cpp.

673 {
674  int vVar = m_session->GetVariables().size();
675  Array<OneD, bool> vChecks(vVar, false);
676  vChecks[vVar - 1] = true;
677  return vChecks;
678 }

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

◆ v_InitObject()

void Nektar::VelocityCorrectionScheme::v_InitObject ( bool  DeclareField = true)
virtual

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::IncNavierStokes.

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

Definition at line 69 of file VelocityCorrectionScheme.cpp.

70 {
71  int n;
72 
73  IncNavierStokes::v_InitObject(DeclareField);
74  m_explicitDiffusion = false;
75 
76  // Set m_pressure to point to last field of m_fields;
77  if (boost::iequals(m_session->GetVariable(m_fields.size() - 1), "p"))
78  {
79  m_nConvectiveFields = m_fields.size() - 1;
81  }
82  else
83  {
84  ASSERTL0(false, "Need to set up pressure field definition");
85  }
86 
87  // Determine diffusion coefficients for each field
88  m_diffCoeff = Array<OneD, NekDouble>(m_nConvectiveFields, m_kinvis);
89  for (n = 0; n < m_nConvectiveFields; ++n)
90  {
91  std::string varName = m_session->GetVariable(n);
92  if (m_session->DefinesFunction("DiffusionCoefficient", varName))
93  {
95  m_session->GetFunction("DiffusionCoefficient", varName);
96  m_diffCoeff[n] = ffunc->Evaluate();
97  }
98  }
99 
100  // Integrate only the convective fields
101  for (n = 0; n < m_nConvectiveFields; ++n)
102  {
103  m_intVariables.push_back(n);
104  }
105 
107  SetUpSVV();
108 
109  // check to see if it is explicity turned off
110  m_session->MatchSolverInfo("GJPStabilisation", "False",
111  m_useGJPStabilisation, true);
112 
113  // if GJPStabilisation set to False bool will be true and
114  // if not false so negate/revese bool
116 
117  m_session->MatchSolverInfo("GJPNormalVelocity", "True", m_useGJPNormalVel,
118  false);
119 
120  if (m_useGJPNormalVel)
121  {
122  ASSERTL0(boost::iequals(m_session->GetSolverInfo("GJPStabilisation"),
123  "Explicit"),
124  "Can only specify GJPNormalVelocity with"
125  " GJPStabilisation set to Explicit currently");
126  }
127 
128  m_session->LoadParameter("GJPJumpScale", m_GJPJumpScale, 1.0);
129 
130  m_session->MatchSolverInfo("SmoothAdvection", "True", m_SmoothAdvection,
131  false);
132 
133  // set explicit time-intregration class operators
136 
137  // set implicit time-intregration class operators
140 
141  // Set up bits for flowrate.
142  m_session->LoadParameter("Flowrate", m_flowrate, 0.0);
143  m_session->LoadParameter("IO_FlowSteps", m_flowrateSteps, 0);
144 }
virtual void v_InitObject(bool DeclareField=true) override
Init object for UnsteadySystem class.
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
bool m_explicitDiffusion
Indicates if explicit or implicit treatment of diffusion is used.
void EvaluateAdvection_SetPressureBCs(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
Array< OneD, NekDouble > m_diffCoeff
Diffusion coefficients (will be kinvis for velocities)
bool m_useGJPNormalVel
bool to identify if GJP normal Velocity should be applied in explicit approach

References ASSERTL0, Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineImplicitSolve(), Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineOdeRhs(), EvaluateAdvection_SetPressureBCs(), m_diffCoeff, Nektar::SolverUtils::UnsteadySystem::m_explicitDiffusion, Nektar::SolverUtils::EquationSystem::m_fields, m_flowrate, m_flowrateSteps, m_GJPJumpScale, 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, m_useGJPNormalVel, m_useGJPStabilisation, SetUpExtrapolation(), SetUpSVV(), SolveUnsteadyStokesSystem(), and Nektar::IncNavierStokes::v_InitObject().

Referenced by Nektar::SmoothedProfileMethod::v_InitObject(), and Nektar::VCSMapping::v_InitObject().

◆ v_PostIntegrate()

bool Nektar::VelocityCorrectionScheme::v_PostIntegrate ( int  step)
protectedvirtual

Reimplemented from Nektar::SolverUtils::AdvectionSystem.

Definition at line 495 of file VelocityCorrectionScheme.cpp.

496 {
497  if (m_flowrateSteps > 0)
498  {
499  if (m_comm->GetRank() == 0 && (step + 1) % m_flowrateSteps == 0)
500  {
501  m_flowrateStream << setw(8) << step << setw(16) << m_time
502  << setw(16) << m_alpha << endl;
503  }
504  }
505 
507 }
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate(int step)

References m_alpha, Nektar::SolverUtils::EquationSystem::m_comm, m_flowrateSteps, m_flowrateStream, Nektar::SolverUtils::EquationSystem::m_time, and Nektar::SolverUtils::AdvectionSystem::v_PostIntegrate().

◆ v_RequireFwdTrans()

virtual bool Nektar::VelocityCorrectionScheme::v_RequireFwdTrans ( )
inlineprotectedvirtual

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 201 of file VelocityCorrectionScheme.h.

202  {
203  return false;
204  }

◆ v_SetUpPressureForcing()

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::VCSWeakPressure, and Nektar::VCSMapping.

Definition at line 792 of file VelocityCorrectionScheme.cpp.

795 {
796  int i;
797  int physTot = m_fields[0]->GetTotPoints();
798  int nvel = m_velocity.size();
799 
800  m_fields[0]->PhysDeriv(eX, fields[0], Forcing[0]);
801 
802  for (i = 1; i < nvel; ++i)
803  {
804  // Use Forcing[1] as storage since it is not needed for the pressure
805  m_fields[i]->PhysDeriv(DirCartesianMap[i], fields[i], Forcing[1]);
806  Vmath::Vadd(physTot, Forcing[1], 1, Forcing[0], 1, Forcing[0], 1);
807  }
808 
809  Vmath::Smul(physTot, 1.0 / aii_Dt, Forcing[0], 1, Forcing[0], 1);
810 }
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:89
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:359
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:248

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().

◆ v_SetUpViscousForcing()

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 815 of file VelocityCorrectionScheme.cpp.

818 {
819  NekDouble aii_dtinv = 1.0 / aii_Dt;
820  int phystot = m_fields[0]->GetTotPoints();
821 
822  // Grad p
823  m_pressure->BwdTrans(m_pressure->GetCoeffs(), m_pressure->UpdatePhys());
824 
825  int nvel = m_velocity.size();
826  if (nvel == 2)
827  {
828  m_pressure->PhysDeriv(m_pressure->GetPhys(), Forcing[m_velocity[0]],
829  Forcing[m_velocity[1]]);
830  }
831  else
832  {
833  m_pressure->PhysDeriv(m_pressure->GetPhys(), Forcing[m_velocity[0]],
834  Forcing[m_velocity[1]], Forcing[m_velocity[2]]);
835  }
836 
837  // zero convective fields.
838  for (int i = nvel; i < m_nConvectiveFields; ++i)
839  {
840  Vmath::Zero(phystot, Forcing[i], 1);
841  }
842 
843  // Subtract inarray/(aii_dt) and divide by kinvis. Kinvis will
844  // need to be updated for the convected fields.
845  for (int i = 0; i < m_nConvectiveFields; ++i)
846  {
847  Blas::Daxpy(phystot, -aii_dtinv, inarray[i], 1, Forcing[i], 1);
848  Blas::Dscal(phystot, 1.0 / m_diffCoeff[i], &(Forcing[i])[0], 1);
849  }
850 }
static void Dscal(const int &n, const double &alpha, double *x, const int &incx)
BLAS level 1: x = alpha x.
Definition: Blas.hpp:168
static void Daxpy(const int &n, const double &alpha, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: y = alpha x plus y.
Definition: Blas.hpp:154

References Blas::Daxpy(), Blas::Dscal(), 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().

◆ v_SolvePressure()

void Nektar::VelocityCorrectionScheme::v_SolvePressure ( const Array< OneD, NekDouble > &  Forcing)
protectedvirtual

Solve pressure system

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

Definition at line 855 of file VelocityCorrectionScheme.cpp.

857 {
859  // Setup coefficient for equation
860  factors[StdRegions::eFactorLambda] = 0.0;
861 
862  // Solver Pressure Poisson Equation
863  m_pressure->HelmSolve(Forcing, m_pressure->UpdateCoeffs(), factors);
864 
865  // Add presure to outflow bc if using convective like BCs
866  m_extrapolation->AddPressureToOutflowBCs(m_kinvis);
867 }
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:282

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

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

◆ v_SolveViscous()

void Nektar::VelocityCorrectionScheme::v_SolveViscous ( const Array< OneD, const Array< OneD, NekDouble >> &  Forcing,
const Array< OneD, const Array< OneD, NekDouble >> &  inarray,
Array< OneD, Array< OneD, NekDouble >> &  outarray,
const NekDouble  aii_Dt 
)
protectedvirtual

Solve velocity system

Reimplemented in Nektar::VCSMapping.

Definition at line 872 of file VelocityCorrectionScheme.cpp.

876 {
880 
881  AppendSVVFactors(factors, varFactorsMap);
882 
883  // Calculate Normal velocity at Trace for GJP explicit stabiliation
884  if (m_useGJPNormalVel)
885  {
887  std::dynamic_pointer_cast<MultiRegions::ContField>(m_fields[0]);
888 
890  cfield->GetGJPForcing();
891 
892  int nTracePts = GJPData->GetNumTracePts();
893  Array<OneD, NekDouble> unorm(nTracePts, 1.0);
894  Array<OneD, NekDouble> Fwd(nTracePts), Bwd(nTracePts);
895  Array<OneD, Array<OneD, NekDouble>> traceNormals =
896  GJPData->GetTraceNormals();
897 
898  m_fields[0]->GetFwdBwdTracePhys(inarray[0], Fwd, Bwd, true, true);
899  Vmath::Vmul(nTracePts, Fwd, 1, traceNormals[0], 1, unorm, 1);
900 
901  // Evaluate u.n on trace
902  for (int f = 1; f < m_fields[0]->GetCoordim(0); ++f)
903  {
904  m_fields[0]->GetFwdBwdTracePhys(inarray[f], Fwd, Bwd, true, true);
905  Vmath::Vvtvp(nTracePts, Fwd, 1, traceNormals[f], 1, unorm, 1, unorm,
906  1);
907  }
908  Vmath::Vabs(nTracePts, unorm, 1, unorm, 1);
909  varCoeffMap[StdRegions::eVarCoeffGJPNormVel] = unorm;
910  }
911 
912  // Solve Helmholtz system and put in Physical space
913  for (int i = 0; i < m_nConvectiveFields; ++i)
914  {
915  // test by adding GJP implicit
917  {
919  }
920 
921  // Setup coefficients for equation
922  factors[StdRegions::eFactorLambda] = 1.0 / aii_Dt / m_diffCoeff[i];
923  m_fields[i]->HelmSolve(Forcing[i], m_fields[i]->UpdateCoeffs(), factors,
924  varCoeffMap, varFactorsMap);
925  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(), outarray[i]);
926  }
927 }
void AppendSVVFactors(StdRegions::ConstFactorMap &factors, MultiRegions::VarFactorsMap &varFactorsMap)
static VarFactorsMap NullVarFactorsMap
std::map< StdRegions::ConstFactorType, Array< OneD, NekDouble > > VarFactorsMap
std::shared_ptr< GJPStabilisation > GJPStabilisationSharedPtr
std::shared_ptr< ContField > ContFieldSharedPtr
Definition: ContField.h:277
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:240
void Vabs(int n, const T *x, const int incx, T *y, const int incy)
vabs: y = |x|
Definition: Vmath.cpp:553

References AppendSVVFactors(), Nektar::StdRegions::eFactorGJP, Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::eVarCoeffGJPNormVel, m_diffCoeff, Nektar::SolverUtils::EquationSystem::m_fields, m_GJPJumpScale, Nektar::IncNavierStokes::m_nConvectiveFields, m_useGJPNormalVel, m_useGJPStabilisation, Nektar::StdRegions::NullVarCoeffMap, Nektar::MultiRegions::NullVarFactorsMap, Vmath::Vabs(), Vmath::Vmul(), and Vmath::Vvtvp().

Referenced by SolveViscous(), and Nektar::VCSMapping::v_SolveViscous().

◆ v_TransCoeffToPhys()

void Nektar::VelocityCorrectionScheme::v_TransCoeffToPhys ( void  )
protectedvirtual

Virtual function for transformation to physical space.

Reimplemented from Nektar::IncNavierStokes.

Definition at line 643 of file VelocityCorrectionScheme.cpp.

644 {
645  int nfields = m_fields.size() - 1;
646  for (int k = 0; k < nfields; ++k)
647  {
648  // Backward Transformation in physical space for time evolution
649  m_fields[k]->BwdTrans(m_fields[k]->GetCoeffs(),
650  m_fields[k]->UpdatePhys());
651  }
652 }

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

◆ v_TransPhysToCoeff()

void Nektar::VelocityCorrectionScheme::v_TransPhysToCoeff ( void  )
protectedvirtual

Virtual function for transformation to coefficient space.

Reimplemented from Nektar::IncNavierStokes.

Definition at line 657 of file VelocityCorrectionScheme.cpp.

658 {
659 
660  int nfields = m_fields.size() - 1;
661  for (int k = 0; k < nfields; ++k)
662  {
663  // Forward Transformation in physical space for time evolution
664  m_fields[k]->FwdTransLocalElmt(m_fields[k]->GetPhys(),
665  m_fields[k]->UpdateCoeffs());
666  }
667 }

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

Member Data Documentation

◆ className

string Nektar::VelocityCorrectionScheme::className
static
Initial value:
=
"VelocityCorrectionScheme", VelocityCorrectionScheme::create)
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
EquationSystemFactory & GetEquationSystemFactory()

Name of class.

Definition at line 58 of file VelocityCorrectionScheme.h.

◆ m_alpha

NekDouble Nektar::VelocityCorrectionScheme::m_alpha
protected

Current flowrate correction.

Definition at line 147 of file VelocityCorrectionScheme.h.

Referenced by SolveUnsteadyStokesSystem(), and v_PostIntegrate().

◆ m_diffCoeff

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

Diffusion coefficients (will be kinvis for velocities)

Definition at line 132 of file VelocityCorrectionScheme.h.

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

◆ m_F

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

◆ m_flowrate

NekDouble Nektar::VelocityCorrectionScheme::m_flowrate
protected

Desired volumetric flowrate.

Definition at line 139 of file VelocityCorrectionScheme.h.

Referenced by SolveUnsteadyStokesSystem(), and v_InitObject().

◆ m_flowrateAiidt

NekDouble Nektar::VelocityCorrectionScheme::m_flowrateAiidt
protected

Value of aii_dt used to compute Stokes flowrate solution.

Definition at line 161 of file VelocityCorrectionScheme.h.

Referenced by SetupFlowrate(), SolveUnsteadyStokesSystem(), and v_DoInitialise().

◆ m_flowrateArea

NekDouble Nektar::VelocityCorrectionScheme::m_flowrateArea
protected

Area of the boundary through which we are measuring the flowrate.

Definition at line 141 of file VelocityCorrectionScheme.h.

Referenced by MeasureFlowrate(), and SetupFlowrate().

◆ m_flowrateBnd

MultiRegions::ExpListSharedPtr Nektar::VelocityCorrectionScheme::m_flowrateBnd
protected

Flowrate reference surface.

Definition at line 153 of file VelocityCorrectionScheme.h.

Referenced by MeasureFlowrate(), and SetupFlowrate().

◆ m_flowrateBndID

int Nektar::VelocityCorrectionScheme::m_flowrateBndID
protected

Boundary ID of the flowrate reference surface.

Definition at line 149 of file VelocityCorrectionScheme.h.

Referenced by MeasureFlowrate(), and SetupFlowrate().

◆ m_flowrateSteps

int Nektar::VelocityCorrectionScheme::m_flowrateSteps
protected

Interval at which to record flowrate data.

Definition at line 159 of file VelocityCorrectionScheme.h.

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

◆ m_flowrateStokes

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

Stokes solution used to impose flowrate.

Definition at line 155 of file VelocityCorrectionScheme.h.

Referenced by SetupFlowrate(), and SolveUnsteadyStokesSystem().

◆ m_flowrateStream

std::ofstream Nektar::VelocityCorrectionScheme::m_flowrateStream
protected

Output stream to record flowrate.

Definition at line 157 of file VelocityCorrectionScheme.h.

Referenced by SetupFlowrate(), and v_PostIntegrate().

◆ m_GJPJumpScale

NekDouble Nektar::VelocityCorrectionScheme::m_GJPJumpScale
protected

Definition at line 119 of file VelocityCorrectionScheme.h.

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

◆ m_greenFlux

NekDouble Nektar::VelocityCorrectionScheme::m_greenFlux
protected

Flux of the Stokes function solution.

Definition at line 145 of file VelocityCorrectionScheme.h.

Referenced by SetupFlowrate(), and SolveUnsteadyStokesSystem().

◆ m_homd1DFlowinPlane

bool Nektar::VelocityCorrectionScheme::m_homd1DFlowinPlane
protected

Definition at line 143 of file VelocityCorrectionScheme.h.

Referenced by MeasureFlowrate(), and SetupFlowrate().

◆ m_IsSVVPowerKernel

bool Nektar::VelocityCorrectionScheme::m_IsSVVPowerKernel
protected

Identifier for Power Kernel otherwise DG kernel.

Definition at line 130 of file VelocityCorrectionScheme.h.

Referenced by AppendSVVFactors(), SetUpSVV(), and v_GenerateSummary().

◆ m_planeID

int Nektar::VelocityCorrectionScheme::m_planeID
protected

Plane ID for cases with homogeneous expansion.

Definition at line 151 of file VelocityCorrectionScheme.h.

Referenced by MeasureFlowrate(), and SetupFlowrate().

◆ m_sVVCutoffRatio

NekDouble Nektar::VelocityCorrectionScheme::m_sVVCutoffRatio
protected

cutt off ratio from which to start decayhing modes

Definition at line 121 of file VelocityCorrectionScheme.h.

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

◆ m_sVVCutoffRatioHomo1D

NekDouble Nektar::VelocityCorrectionScheme::m_sVVCutoffRatioHomo1D
protected

Definition at line 124 of file VelocityCorrectionScheme.h.

Referenced by SetUpSVV(), and v_GenerateSummary().

◆ m_sVVDiffCoeff

NekDouble Nektar::VelocityCorrectionScheme::m_sVVDiffCoeff
protected

◆ m_sVVDiffCoeffHomo1D

NekDouble Nektar::VelocityCorrectionScheme::m_sVVDiffCoeffHomo1D
protected

Diffusion coefficient of SVV modes in homogeneous 1D Direction.

Definition at line 126 of file VelocityCorrectionScheme.h.

Referenced by SetUpSVV(), and v_GenerateSummary().

◆ m_svvVarDiffCoeff

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

Array of coefficient if power kernel is used in SVV.

Definition at line 128 of file VelocityCorrectionScheme.h.

Referenced by AppendSVVFactors(), SetUpSVV(), and v_GenerateSummary().

◆ m_useGJPNormalVel

bool Nektar::VelocityCorrectionScheme::m_useGJPNormalVel
protected

bool to identify if GJP normal Velocity should be applied in explicit approach

Definition at line 117 of file VelocityCorrectionScheme.h.

Referenced by v_InitObject(), and v_SolveViscous().

◆ m_useGJPStabilisation

bool Nektar::VelocityCorrectionScheme::m_useGJPStabilisation
protected

bool to identify if GJP semi-implicit is active.

Definition at line 114 of file VelocityCorrectionScheme.h.

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

◆ m_useHomo1DSpecVanVisc

bool Nektar::VelocityCorrectionScheme::m_useHomo1DSpecVanVisc
protected

bool to identify if spectral vanishing viscosity is active.

Definition at line 110 of file VelocityCorrectionScheme.h.

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

◆ m_useSpecVanVisc

bool Nektar::VelocityCorrectionScheme::m_useSpecVanVisc
protected

bool to identify if spectral vanishing viscosity is active.

Definition at line 112 of file VelocityCorrectionScheme.h.

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

◆ m_varCoeffLap

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 136 of file VelocityCorrectionScheme.h.