Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | Static Protected Attributes | Private Member Functions | Private Attributes | List of all members
Nektar::VCSMapping Class Reference

#include <VCSMapping.h>

Inheritance diagram for Nektar::VCSMapping:
[legend]

Public Member Functions

 VCSMapping (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Constructor. More...
 
void ApplyIncNSMappingForcing (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
virtual ~VCSMapping ()
 
virtual void v_InitObject (bool DeclareField=true) override
 Init object for UnsteadySystem class. More...
 
- Public Member Functions inherited from Nektar::VelocityCorrectionScheme
 VelocityCorrectionScheme (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Constructor. More...
 
virtual ~VelocityCorrectionScheme ()
 
virtual void v_InitObject (bool DeclareField=true) override
 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 ()
 
virtual void v_InitObject (bool DeclareField=true) override
 Init object for UnsteadySystem class. More...
 
int GetNConvectiveFields (void)
 
void AddForcing (const SolverUtils::ForcingSharedPtr &pForce)
 
virtual void v_GetPressure (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &pressure) override
 
virtual void v_GetDensity (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &density) override
 
virtual bool v_HasConstantDensity () override
 
virtual void v_GetVelocity (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &velocity) override
 
virtual void v_SetMovingFrameVelocities (const Array< OneD, NekDouble > &vFrameVels) override
 
virtual void v_GetMovingFrameVelocities (Array< OneD, NekDouble > &vFrameVels) override
 
virtual void v_SetMovingFrameAngles (const Array< OneD, NekDouble > &vFrameTheta) override
 
virtual void v_GetMovingFrameAngles (Array< OneD, NekDouble > &vFrameTheta) override
 
virtual void v_SetMovingFrameProjectionMat (const bnu::matrix< NekDouble > &vProjMat) override
 
virtual void v_GetMovingFrameProjectionMat (bnu::matrix< NekDouble > &vProjMat) override
 
bool DefinedForcing (const std::string &sForce)
 
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 ()
 
virtual SOLVER_UTILS_EXPORT void v_InitObject (bool DeclareField=true) override
 Init object for UnsteadySystem class. More...
 
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)
 
SOLVER_UTILS_EXPORT LibUtilities::TimeIntegrationSchemeSharedPtrGetTimeIntegrationScheme ()
 Returns the time integration scheme. More...
 
SOLVER_UTILS_EXPORT LibUtilities::TimeIntegrationSchemeOperatorsGetTimeIntegrationSchemeOperators ()
 Returns the time integration scheme operators. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT void InitObject (bool DeclareField=true)
 Initialises the members of this object. More...
 
SOLVER_UTILS_EXPORT void DoInitialise (bool dumpInitialConditions=true)
 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 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 NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation. 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, const Array< OneD, const NekDouble > &input)
 
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 int GetInfoSteps ()
 
SOLVER_UTILS_EXPORT void SetInfoSteps (int num)
 
SOLVER_UTILS_EXPORT void SetIterationNumberPIT (int num)
 
SOLVER_UTILS_EXPORT void SetWindowNumberPIT (int num)
 
SOLVER_UTILS_EXPORT Array< OneD, const Array< OneD, NekDouble > > GetTraceNormals ()
 
SOLVER_UTILS_EXPORT void SetTime (const NekDouble time)
 
SOLVER_UTILS_EXPORT void SetTimeStep (const NekDouble timestep)
 
SOLVER_UTILS_EXPORT void SetInitialStep (const int step)
 
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time. More...
 
SOLVER_UTILS_EXPORT bool NegatedOp ()
 Identify if operator is negated in DoSolve. More...
 
SOLVER_UTILS_EXPORT bool ParallelInTime ()
 Check if solver use Parallel-in-Time. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::FluidInterface
virtual ~FluidInterface ()=default
 
SOLVER_UTILS_EXPORT void GetVelocity (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &velocity)
 Extract array with velocity from physfield. More...
 
SOLVER_UTILS_EXPORT bool HasConstantDensity ()
 
SOLVER_UTILS_EXPORT void GetDensity (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &density)
 Extract array with density from physfield. More...
 
SOLVER_UTILS_EXPORT void GetPressure (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &pressure)
 Extract array with pressure from physfield. More...
 
SOLVER_UTILS_EXPORT void SetMovingFrameVelocities (const Array< OneD, NekDouble > &vFrameVels)
 
SOLVER_UTILS_EXPORT void GetMovingFrameVelocities (Array< OneD, NekDouble > &vFrameVels)
 
SOLVER_UTILS_EXPORT void SetMovingFrameProjectionMat (const boost::numeric::ublas::matrix< NekDouble > &vProjMat)
 
SOLVER_UTILS_EXPORT void GetMovingFrameProjectionMat (boost::numeric::ublas::matrix< NekDouble > &vProjMat)
 
SOLVER_UTILS_EXPORT void SetMovingFrameAngles (const Array< OneD, NekDouble > &vFrameTheta)
 
SOLVER_UTILS_EXPORT void GetMovingFrameAngles (Array< OneD, NekDouble > &vFrameTheta)
 

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 Member Functions inherited from Nektar::VelocityCorrectionScheme
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::VelocityCorrectionScheme
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

virtual void v_DoInitialise (bool dumpInitialConditions=true) override
 Sets up initial conditions. More...
 
virtual void v_SetUpPressureForcing (const Array< OneD, const Array< OneD, NekDouble > > &fields, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt) override
 
virtual void v_SetUpViscousForcing (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt) override
 
virtual void v_SolvePressure (const Array< OneD, NekDouble > &Forcing) override
 
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) override
 
virtual void v_EvaluateAdvection_SetPressureBCs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time) override
 
- Protected Member Functions inherited from Nektar::VelocityCorrectionScheme
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) override
 
virtual void v_GenerateSummary (SolverUtils::SummaryList &s) override
 Print a summary of time stepping parameters. More...
 
virtual void v_TransCoeffToPhys (void) override
 Virtual function for transformation to physical space. More...
 
virtual void v_TransPhysToCoeff (void) override
 Virtual function for transformation to coefficient space. More...
 
virtual void v_DoInitialise (bool dumpInitialConditions=true) override
 Sets up initial conditions. More...
 
virtual Array< OneD, bool > v_GetSystemSingularChecks () override
 
virtual int v_GetForceDimension () override
 
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_SolveUnsteadyStokesSystem (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time, const NekDouble a_iixDt)
 
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 () override
 
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)
 
void ComputeGJPNormalVelocity (const Array< OneD, const Array< OneD, NekDouble > > &inarray, StdRegions::VarCoeffMap &varcoeffs)
 
- 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 void v_TransCoeffToPhys (void) override
 Virtual function for transformation to physical space. More...
 
virtual void v_TransPhysToCoeff (void) override
 Virtual function for transformation to coefficient space. More...
 
virtual int v_GetForceDimension ()=0
 
virtual Array< OneD, NekDoublev_GetMaxStdVelocity (const NekDouble SpeedSoundFactor) override
 
virtual bool v_PreIntegrate (int step) override
 
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step) override
 
virtual SOLVER_UTILS_EXPORT Array< OneD, NekDoublev_GetMaxStdVelocity (const NekDouble SpeedSoundFactor=1.0)
 
- 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...
 
virtual SOLVER_UTILS_EXPORT void v_InitObject (bool DeclareField=true) override
 Init object for UnsteadySystem class. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve () override
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_PrintStatusInformation (const int step, const NekDouble cpuTime)
 Print Status Information. More...
 
virtual SOLVER_UTILS_EXPORT void v_PrintSummaryStatistics (const NekDouble intTime)
 Print Summary Statistics. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise (bool dumpInitialConditions=true) override
 Sets up initial conditions. More...
 
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary (SummaryList &s) override
 Print a summary of time stepping parameters. 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 bool v_PreIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT bool v_RequireFwdTrans ()
 
virtual SOLVER_UTILS_EXPORT void v_SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
virtual SOLVER_UTILS_EXPORT bool v_UpdateTimeStepCheck ()
 
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control. More...
 
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...
 
SOLVER_UTILS_EXPORT void DoDummyProjection (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 Perform dummy projection. More...
 
- 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 void v_InitObject (bool DeclareFeld=true)
 Initialisation object for EquationSystem. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise (bool dumpInitialConditions=true)
 Virtual function for initialisation implementation. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve ()
 Virtual function for solve implementation. 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_TransCoeffToPhys ()
 Virtual function for transformation to physical space. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransPhysToCoeff ()
 Virtual function for transformation to coefficient space. More...
 
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary (SummaryList &l)
 Virtual function for generating summary information. 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 MultiRegions::ExpListSharedPtr v_GetPressure (void)
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp (void)
 Virtual function to identify if operator is negated in DoSolve. More...
 
virtual SOLVER_UTILS_EXPORT void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 
- Protected Member Functions inherited from Nektar::SolverUtils::FluidInterface
virtual SOLVER_UTILS_EXPORT void v_GetVelocity (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &velocity)=0
 
virtual SOLVER_UTILS_EXPORT bool v_HasConstantDensity ()=0
 
virtual SOLVER_UTILS_EXPORT void v_GetDensity (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &density)=0
 
virtual SOLVER_UTILS_EXPORT void v_GetPressure (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &pressure)=0
 
virtual SOLVER_UTILS_EXPORT void v_SetMovingFrameVelocities (const Array< OneD, NekDouble > &vFrameVels)
 
virtual SOLVER_UTILS_EXPORT void v_GetMovingFrameVelocities (Array< OneD, NekDouble > &vFrameVels)
 
virtual SOLVER_UTILS_EXPORT void v_SetMovingFrameProjectionMat (const boost::numeric::ublas::matrix< NekDouble > &vProjMat)
 
virtual SOLVER_UTILS_EXPORT void v_GetMovingFrameProjectionMat (boost::numeric::ublas::matrix< NekDouble > &vProjMat)
 
virtual SOLVER_UTILS_EXPORT void v_SetMovingFrameAngles (const Array< OneD, NekDouble > &vFrameTheta)
 
virtual SOLVER_UTILS_EXPORT void v_GetMovingFrameAngles (Array< OneD, NekDouble > &vFrameTheta)
 

Protected Attributes

GlobalMapping::MappingSharedPtr m_mapping
 
bool m_verbose
 
bool m_implicitPressure
 
bool m_implicitViscous
 
bool m_neglectViscous
 
NekDouble m_pressureTolerance
 
NekDouble m_viscousTolerance
 
NekDouble m_pressureRelaxation
 
NekDouble m_viscousRelaxation
 
Array< OneD, Array< OneD, NekDouble > > m_gradP
 
- Protected Attributes inherited from Nektar::VelocityCorrectionScheme
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
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
Array< OneD, Array< OneD, NekDouble > > m_previousSolution
 Storage for previous solution for steady-state check. More...
 
std::vector< int > m_intVariables
 
NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1). More...
 
NekDouble m_CFLGrowth
 CFL growth rate. More...
 
NekDouble m_CFLEnd
 Maximun cfl in cfl growth. More...
 
int m_abortSteps
 Number of steps between checks for abort conditions. More...
 
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...
 
int m_steadyStateSteps
 Check for steady state at step interval. More...
 
NekDouble m_steadyStateTol
 Tolerance to which steady state should be evaluated at. More...
 
int m_filtersInfosteps
 Number of time steps between outputting filters information. More...
 
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
 
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state. More...
 
std::ofstream m_errFile
 
NekDouble m_epsilon
 Diffusion coefficient. More...
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
bool m_verbose
 
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_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_infosteps
 Number of time steps between outputting status information. More...
 
int m_iterPIT = 0
 Number of parallel-in-time time iteration. More...
 
int m_windowPIT = 0
 Index of windows for parallel-in-time time iteration. 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...
 

Static Protected Attributes

static std::string solverTypeLookupId
 
- Static Protected Attributes inherited from Nektar::VelocityCorrectionScheme
static std::string solverTypeLookupId
 
- Static Protected Attributes inherited from Nektar::IncNavierStokes
static std::string eqTypeLookupIds []
 
- Static Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
static std::string equationSystemTypeLookupIds []
 
static std::string projectionTypeLookupIds []
 

Private Member Functions

void MappingAdvectionCorrection (const Array< OneD, const Array< OneD, NekDouble > > &velPhys, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
void MappingAccelerationCorrection (const Array< OneD, const Array< OneD, NekDouble > > &vel, const Array< OneD, const Array< OneD, NekDouble > > &velPhys, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
void MappingPressureCorrection (Array< OneD, Array< OneD, NekDouble > > &outarray)
 
void MappingViscousCorrection (const Array< OneD, const Array< OneD, NekDouble > > &velPhys, Array< OneD, Array< OneD, NekDouble > > &outarray)
 

Private Attributes

Array< OneD, Array< OneD, NekDouble > > m_presForcingCorrection
 

Additional Inherited Members

- Protected Types inherited from Nektar::SolverUtils::EquationSystem
enum  HomogeneousType { eHomogeneous1D , eHomogeneous2D , eHomogeneous3D , eNotHomogeneous }
 Parameter for homogeneous expansions. More...
 

Detailed Description

Definition at line 43 of file VCSMapping.h.

Constructor & Destructor Documentation

◆ VCSMapping()

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

Constructor.

Constructor. Creates ...

Parameters

param

Definition at line 59 of file VCSMapping.cpp.

61 : UnsteadySystem(pSession, pGraph),
62 VelocityCorrectionScheme(pSession, pGraph)
63{
64}
SOLVER_UTILS_EXPORT UnsteadySystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Initialises UnsteadySystem class members.
VelocityCorrectionScheme(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Constructor.

◆ ~VCSMapping()

Nektar::VCSMapping::~VCSMapping ( void  )
virtual

Destructor

Definition at line 131 of file VCSMapping.cpp.

132{
133}

Member Function Documentation

◆ ApplyIncNSMappingForcing()

void Nektar::VCSMapping::ApplyIncNSMappingForcing ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)

Explicit terms of the mapping

Definition at line 723 of file VCSMapping.cpp.

726{
727 size_t physTot = m_fields[0]->GetTotPoints();
728 Array<OneD, Array<OneD, NekDouble>> vel(m_nConvectiveFields);
729 Array<OneD, Array<OneD, NekDouble>> velPhys(m_nConvectiveFields);
730 Array<OneD, Array<OneD, NekDouble>> Forcing(m_nConvectiveFields);
731 Array<OneD, Array<OneD, NekDouble>> tmp(m_nConvectiveFields);
732 for (int i = 0; i < m_nConvectiveFields; ++i)
733 {
734 velPhys[i] = Array<OneD, NekDouble>(physTot, 0.0);
735 Forcing[i] = Array<OneD, NekDouble>(physTot, 0.0);
736 tmp[i] = Array<OneD, NekDouble>(physTot, 0.0);
737 }
738
739 // Get fields and store velocity in wavespace and physical space
740 if (m_fields[0]->GetWaveSpace())
741 {
742 for (int i = 0; i < m_nConvectiveFields; ++i)
743 {
744 vel[i] = inarray[i];
745 m_fields[0]->HomogeneousBwdTrans(physTot, vel[i], velPhys[i]);
746 }
747 }
748 else
749 {
750 for (int i = 0; i < m_nConvectiveFields; ++i)
751 {
752 vel[i] = inarray[i];
753 Vmath::Vcopy(physTot, inarray[i], 1, velPhys[i], 1);
754 }
755 }
756
757 // Advection contribution
758 MappingAdvectionCorrection(velPhys, Forcing);
759
760 // Time-derivative contribution
761 if (m_mapping->IsTimeDependent())
762 {
763 MappingAccelerationCorrection(vel, velPhys, tmp);
764 for (int i = 0; i < m_nConvectiveFields; ++i)
765 {
766 Vmath::Vadd(physTot, tmp[i], 1, Forcing[i], 1, Forcing[i], 1);
767 }
768 }
769
770 // Pressure contribution
772 {
774 for (int i = 0; i < m_nConvectiveFields; ++i)
775 {
776 Vmath::Vadd(physTot, tmp[i], 1, Forcing[i], 1, Forcing[i], 1);
777 }
778 }
779 // Viscous contribution
781 {
782 MappingViscousCorrection(velPhys, tmp);
783 for (int i = 0; i < m_nConvectiveFields; ++i)
784 {
785 Vmath::Smul(physTot, m_kinvis, tmp[i], 1, tmp[i], 1);
786 Vmath::Vadd(physTot, tmp[i], 1, Forcing[i], 1, Forcing[i], 1);
787 }
788 }
789
790 // If necessary, transform to wavespace
791 if (m_fields[0]->GetWaveSpace())
792 {
793 for (int i = 0; i < m_nConvectiveFields; ++i)
794 {
795 m_fields[0]->HomogeneousFwdTrans(physTot, Forcing[i], Forcing[i]);
796 }
797 }
798
799 // Add to outarray
800 for (int i = 0; i < m_nConvectiveFields; ++i)
801 {
802 Vmath::Vadd(physTot, outarray[i], 1, Forcing[i], 1, outarray[i], 1);
803 }
804}
NekDouble m_kinvis
Kinematic viscosity.
int m_nConvectiveFields
Number of fields to be convected;.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
void MappingViscousCorrection(const Array< OneD, const Array< OneD, NekDouble > > &velPhys, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: VCSMapping.cpp:936
GlobalMapping::MappingSharedPtr m_mapping
Definition: VCSMapping.h:75
void MappingPressureCorrection(Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: VCSMapping.cpp:910
void MappingAccelerationCorrection(const Array< OneD, const Array< OneD, NekDouble > > &vel, const Array< OneD, const Array< OneD, NekDouble > > &velPhys, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: VCSMapping.cpp:831
void MappingAdvectionCorrection(const Array< OneD, const Array< OneD, NekDouble > > &velPhys, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: VCSMapping.cpp:806
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:354
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:245
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191

References Nektar::SolverUtils::EquationSystem::m_fields, m_implicitPressure, m_implicitViscous, Nektar::IncNavierStokes::m_kinvis, m_mapping, Nektar::IncNavierStokes::m_nConvectiveFields, m_neglectViscous, MappingAccelerationCorrection(), MappingAdvectionCorrection(), MappingPressureCorrection(), MappingViscousCorrection(), Vmath::Smul(), Vmath::Vadd(), and Vmath::Vcopy().

Referenced by v_EvaluateAdvection_SetPressureBCs().

◆ create()

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

Creates an instance of this class.

Definition at line 47 of file VCSMapping.h.

50 {
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.

◆ MappingAccelerationCorrection()

void Nektar::VCSMapping::MappingAccelerationCorrection ( const Array< OneD, const Array< OneD, NekDouble > > &  vel,
const Array< OneD, const Array< OneD, NekDouble > > &  velPhys,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
private

Definition at line 831 of file VCSMapping.cpp.

835{
836 size_t physTot = m_fields[0]->GetTotPoints();
837 size_t nvel = m_nConvectiveFields;
838
839 Array<OneD, Array<OneD, NekDouble>> wk(nvel * nvel);
840 Array<OneD, Array<OneD, NekDouble>> tmp(nvel);
841 Array<OneD, Array<OneD, NekDouble>> coordVel(nvel);
842 for (size_t i = 0; i < nvel; i++)
843 {
844 tmp[i] = Array<OneD, NekDouble>(physTot, 0.0);
845 coordVel[i] = Array<OneD, NekDouble>(physTot, 0.0);
846 }
847 // Get coordinates velocity in transformed system
848 m_mapping->GetCoordVelocity(tmp);
849 m_mapping->ContravarFromCartesian(tmp, coordVel);
850
851 // Calculate first term: U^j u^i,j = U^j (du^i/dx^j + {i,kj}u^k)
852 m_mapping->ApplyChristoffelContravar(velPhys, wk);
853 for (size_t i = 0; i < nvel; i++)
854 {
855 Vmath::Zero(physTot, outarray[i], 1);
856
857 m_fields[0]->PhysDeriv(velPhys[i], tmp[0], tmp[1]);
858 for (size_t j = 0; j < nvel; j++)
859 {
860 if (j == 2)
861 {
862 m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[j], vel[i],
863 tmp[2]);
864 if (m_fields[0]->GetWaveSpace())
865 {
866 m_fields[0]->HomogeneousBwdTrans(physTot, tmp[2], tmp[2]);
867 }
868 }
869
870 Vmath::Vadd(physTot, wk[i * nvel + j], 1, tmp[j], 1,
871 wk[i * nvel + j], 1);
872
873 Vmath::Vvtvp(physTot, coordVel[j], 1, wk[i * nvel + j], 1,
874 outarray[i], 1, outarray[i], 1);
875 }
876 }
877
878 // Set wavespace to false and store current value
879 bool wavespace = m_fields[0]->GetWaveSpace();
880 m_fields[0]->SetWaveSpace(false);
881
882 // Add -u^j U^i,j
883 m_mapping->ApplyChristoffelContravar(coordVel, wk);
884 for (size_t i = 0; i < nvel; i++)
885 {
886 if (nvel == 2)
887 {
888 m_fields[0]->PhysDeriv(coordVel[i], tmp[0], tmp[1]);
889 }
890 else
891 {
892 m_fields[0]->PhysDeriv(coordVel[i], tmp[0], tmp[1], tmp[2]);
893 }
894
895 for (size_t j = 0; j < nvel; j++)
896 {
897 Vmath::Vadd(physTot, wk[i * nvel + j], 1, tmp[j], 1,
898 wk[i * nvel + j], 1);
899 Vmath::Neg(physTot, wk[i * nvel + j], 1);
900
901 Vmath::Vvtvp(physTot, velPhys[j], 1, wk[i * nvel + j], 1,
902 outarray[i], 1, outarray[i], 1);
903 }
904 }
905
906 // Restore value of wavespace
907 m_fields[0]->SetWaveSpace(wavespace);
908}
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:90
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:513
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:569
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:487

References Nektar::MultiRegions::DirCartesianMap, Nektar::SolverUtils::EquationSystem::m_fields, m_mapping, Nektar::IncNavierStokes::m_nConvectiveFields, Vmath::Neg(), Vmath::Vadd(), Vmath::Vvtvp(), and Vmath::Zero().

Referenced by ApplyIncNSMappingForcing().

◆ MappingAdvectionCorrection()

void Nektar::VCSMapping::MappingAdvectionCorrection ( const Array< OneD, const Array< OneD, NekDouble > > &  velPhys,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
private

Definition at line 806 of file VCSMapping.cpp.

809{
810 size_t physTot = m_fields[0]->GetTotPoints();
811 size_t nvel = m_nConvectiveFields;
812
813 Array<OneD, Array<OneD, NekDouble>> wk(nvel * nvel);
814
815 // Apply Christoffel symbols to obtain {i,kj}vel(k)
816 m_mapping->ApplyChristoffelContravar(velPhys, wk);
817
818 // Calculate correction -U^j*{i,kj}vel(k)
819 for (size_t i = 0; i < nvel; i++)
820 {
821 Vmath::Zero(physTot, outarray[i], 1);
822 for (size_t j = 0; j < nvel; j++)
823 {
824 Vmath::Vvtvp(physTot, wk[i * nvel + j], 1, velPhys[j], 1,
825 outarray[i], 1, outarray[i], 1);
826 }
827 Vmath::Neg(physTot, outarray[i], 1);
828 }
829}

References Nektar::SolverUtils::EquationSystem::m_fields, m_mapping, Nektar::IncNavierStokes::m_nConvectiveFields, Vmath::Neg(), Vmath::Vvtvp(), and Vmath::Zero().

Referenced by ApplyIncNSMappingForcing().

◆ MappingPressureCorrection()

void Nektar::VCSMapping::MappingPressureCorrection ( Array< OneD, Array< OneD, NekDouble > > &  outarray)
private

Definition at line 910 of file VCSMapping.cpp.

912{
913 size_t physTot = m_fields[0]->GetTotPoints();
914 size_t nvel = m_nConvectiveFields;
915
916 // Calculate g^(ij)p_(,j)
917 m_mapping->RaiseIndex(m_gradP, outarray);
918
919 // Calculate correction = (nabla p)/J - g^(ij)p_,j
920 // (Jac is not required if it is constant)
921 if (!m_mapping->HasConstantJacobian())
922 {
923 Array<OneD, NekDouble> Jac(physTot, 0.0);
924 m_mapping->GetJacobian(Jac);
925 for (size_t i = 0; i < nvel; ++i)
926 {
927 Vmath::Vdiv(physTot, m_gradP[i], 1, Jac, 1, m_gradP[i], 1);
928 }
929 }
930 for (size_t i = 0; i < nvel; ++i)
931 {
932 Vmath::Vsub(physTot, m_gradP[i], 1, outarray[i], 1, outarray[i], 1);
933 }
934}
Array< OneD, Array< OneD, NekDouble > > m_gradP
Definition: VCSMapping.h:94
void Vdiv(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:280
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:414

References Nektar::SolverUtils::EquationSystem::m_fields, m_gradP, m_mapping, Nektar::IncNavierStokes::m_nConvectiveFields, Vmath::Vdiv(), and Vmath::Vsub().

Referenced by ApplyIncNSMappingForcing().

◆ MappingViscousCorrection()

void Nektar::VCSMapping::MappingViscousCorrection ( const Array< OneD, const Array< OneD, NekDouble > > &  velPhys,
Array< OneD, Array< OneD, NekDouble > > &  outarray 
)
private

Definition at line 936 of file VCSMapping.cpp.

939{
940 // L(U) - 1.0*d^2(u^i)/dx^jdx^j
941 m_mapping->VelocityLaplacian(velPhys, outarray, 1.0);
942}

References m_mapping.

Referenced by ApplyIncNSMappingForcing().

◆ v_DoInitialise()

void Nektar::VCSMapping::v_DoInitialise ( bool  dumpInitialConditions = true)
overrideprotectedvirtual

Sets up initial conditions.

Sets the initial conditions.

Reimplemented from Nektar::VelocityCorrectionScheme.

Definition at line 135 of file VCSMapping.cpp.

136{
137 UnsteadySystem::v_DoInitialise(dumpInitialConditions);
138
139 // Set up Field Meta Data for output files
140 m_fieldMetaDataMap["Kinvis"] = boost::lexical_cast<std::string>(m_kinvis);
141 m_fieldMetaDataMap["TimeStep"] =
142 boost::lexical_cast<std::string>(m_timestep);
143
144 // Correct Dirichlet boundary conditions to account for mapping
145 m_mapping->UpdateBCs(0.0);
146 //
147 m_F = Array<OneD, Array<OneD, NekDouble>>(m_nConvectiveFields);
148 for (int i = 0; i < m_nConvectiveFields; ++i)
149 {
150 m_fields[i]->ImposeDirichletConditions(m_fields[i]->UpdateCoeffs());
151 m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
152 m_fields[i]->UpdatePhys());
153 m_F[i] = Array<OneD, NekDouble>(m_fields[0]->GetTotPoints(), 0.0);
154 }
155
156 // Initialise m_gradP
157 size_t physTot = m_fields[0]->GetTotPoints();
158 m_gradP = Array<OneD, Array<OneD, NekDouble>>(m_nConvectiveFields);
159 for (int i = 0; i < m_nConvectiveFields; ++i)
160 {
161 m_gradP[i] = Array<OneD, NekDouble>(physTot, 0.0);
163 m_pressure->GetPhys(), m_gradP[i]);
164 if (m_pressure->GetWaveSpace())
165 {
166 m_pressure->HomogeneousBwdTrans(physTot, m_gradP[i], m_gradP[i]);
167 }
168 }
169}
MultiRegions::ExpListSharedPtr m_pressure
Pointer to field holding pressure field.
NekDouble m_timestep
Time step size.
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(bool dumpInitialConditions=true) override
Sets up initial conditions.
Array< OneD, Array< OneD, NekDouble > > m_F

References Nektar::MultiRegions::DirCartesianMap, Nektar::SolverUtils::EquationSystem::GetTotPoints(), Nektar::VelocityCorrectionScheme::m_F, Nektar::SolverUtils::EquationSystem::m_fieldMetaDataMap, Nektar::SolverUtils::EquationSystem::m_fields, m_gradP, Nektar::IncNavierStokes::m_kinvis, m_mapping, Nektar::IncNavierStokes::m_nConvectiveFields, Nektar::IncNavierStokes::m_pressure, Nektar::SolverUtils::EquationSystem::m_timestep, and Nektar::SolverUtils::UnsteadySystem::v_DoInitialise().

◆ v_EvaluateAdvection_SetPressureBCs()

void Nektar::VCSMapping::v_EvaluateAdvection_SetPressureBCs ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
const NekDouble  time 
)
overrideprotectedvirtual

Explicit part of the method - Advection, Forcing + HOPBCs

Reimplemented from Nektar::VelocityCorrectionScheme.

Definition at line 174 of file VCSMapping.cpp.

177{
178 EvaluateAdvectionTerms(inarray, outarray, time);
179
180 // Smooth advection
182 {
183 for (int i = 0; i < m_nConvectiveFields; ++i)
184 {
185 m_pressure->SmoothField(outarray[i]);
186 }
187 }
188
189 // Add forcing terms
190 for (auto &x : m_forcing)
191 {
192 x->Apply(m_fields, inarray, outarray, time);
193 }
194
195 // Add mapping terms
196 ApplyIncNSMappingForcing(inarray, outarray);
197
198 // Calculate High-Order pressure boundary conditions
199 m_extrapolation->EvaluatePressureBCs(inarray, outarray, m_kinvis);
200
201 // Update mapping and deal with Dirichlet boundary conditions
202 if (m_mapping->IsTimeDependent())
203 {
204 if (m_mapping->IsFromFunction())
205 {
206 // If the transformation is explicitly defined, update it here
207 // Otherwise, it will be done somewhere else (ForcingMovingBody)
208 m_mapping->UpdateMapping(time + m_timestep);
209 }
210 m_mapping->UpdateBCs(time + m_timestep);
211 }
212}
bool m_SmoothAdvection
bool to identify if advection term smoothing is requested
ExtrapolateSharedPtr m_extrapolation
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Forcing terms.
void EvaluateAdvectionTerms(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
void ApplyIncNSMappingForcing(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: VCSMapping.cpp:723

References ApplyIncNSMappingForcing(), Nektar::IncNavierStokes::EvaluateAdvectionTerms(), Nektar::IncNavierStokes::m_extrapolation, Nektar::SolverUtils::EquationSystem::m_fields, Nektar::IncNavierStokes::m_forcing, Nektar::IncNavierStokes::m_kinvis, m_mapping, Nektar::IncNavierStokes::m_nConvectiveFields, Nektar::IncNavierStokes::m_pressure, Nektar::IncNavierStokes::m_SmoothAdvection, and Nektar::SolverUtils::EquationSystem::m_timestep.

◆ v_InitObject()

void Nektar::VCSMapping::v_InitObject ( bool  DeclareField = true)
overridevirtual

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::VelocityCorrectionScheme.

Definition at line 66 of file VCSMapping.cpp.

67{
69
71 ASSERTL0(m_mapping, "Could not create mapping in VCSMapping.");
72
73 std::string vExtrapolation = "Mapping";
75 vExtrapolation, m_session, m_fields, m_pressure, m_velocity,
77 m_extrapolation->SubSteppingTimeIntegration(m_intScheme);
78 m_extrapolation->GenerateHOPBCMap(m_session);
79
80 // Storage to extrapolate pressure forcing
81 size_t physTot = m_fields[0]->GetTotPoints();
82 size_t intSteps = 1;
83
84 if (m_intScheme->GetName() == "IMEX" ||
85 m_intScheme->GetName() == "IMEXGear")
86 {
87 m_intSteps = m_intScheme->GetOrder();
88 }
89 else
90 {
92 "Integration method not suitable: "
93 "Options include IMEXGear or IMEXOrder{1,2,3,4}");
94 }
95
96 m_presForcingCorrection = Array<OneD, Array<OneD, NekDouble>>(intSteps);
97 for (size_t i = 0; i < m_presForcingCorrection.size(); i++)
98 {
99 m_presForcingCorrection[i] = Array<OneD, NekDouble>(physTot, 0.0);
100 }
101 m_verbose = (m_session->DefinesCmdLineArgument("verbose")) ? true : false;
102
103 // Load solve parameters related to the mapping
104 // Flags determining if pressure/viscous terms should be treated implicitly
105 m_session->MatchSolverInfo("MappingImplicitPressure", "True",
106 m_implicitPressure, false);
107 m_session->MatchSolverInfo("MappingImplicitViscous", "True",
108 m_implicitViscous, false);
109 m_session->MatchSolverInfo("MappingNeglectViscous", "True",
110 m_neglectViscous, false);
111
113 {
114 m_implicitViscous = false;
115 }
116
117 // Tolerances and relaxation parameters for implicit terms
118 m_session->LoadParameter("MappingPressureTolerance", m_pressureTolerance,
119 1e-12);
120 m_session->LoadParameter("MappingViscousTolerance", m_viscousTolerance,
121 1e-12);
122 m_session->LoadParameter("MappingPressureRelaxation", m_pressureRelaxation,
123 1.0);
124 m_session->LoadParameter("MappingViscousRelaxation", m_viscousRelaxation,
125 1.0);
126}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
static GLOBAL_MAPPING_EXPORT MappingSharedPtr Load(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
Return a pointer to the mapping, creating it on first call.
Definition: Mapping.cpp:272
Array< OneD, int > m_velocity
int which identifies which components of m_fields contains the velocity (u,v,w);
int m_intSteps
Number of time integration steps AND Order of extrapolation for pressure boundary conditions.
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::SessionReaderSharedPtr m_session
The session reader.
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
Wrapper to the time integration scheme.
NekDouble m_pressureTolerance
Definition: VCSMapping.h:88
NekDouble m_viscousRelaxation
Definition: VCSMapping.h:91
NekDouble m_pressureRelaxation
Definition: VCSMapping.h:90
NekDouble m_viscousTolerance
Definition: VCSMapping.h:89
Array< OneD, Array< OneD, NekDouble > > m_presForcingCorrection
Definition: VCSMapping.h:124
virtual void v_InitObject(bool DeclareField=true) override
Init object for UnsteadySystem class.
ExtrapolateFactory & GetExtrapolateFactory()
Definition: Extrapolate.cpp:48

References ASSERTL0, Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::ErrorUtil::efatal, Nektar::GetExtrapolateFactory(), Nektar::GlobalMapping::Mapping::Load(), Nektar::SolverUtils::AdvectionSystem::m_advObject, Nektar::IncNavierStokes::m_extrapolation, Nektar::SolverUtils::EquationSystem::m_fields, m_implicitPressure, m_implicitViscous, Nektar::SolverUtils::UnsteadySystem::m_intScheme, Nektar::IncNavierStokes::m_intSteps, m_mapping, m_neglectViscous, m_presForcingCorrection, Nektar::IncNavierStokes::m_pressure, m_pressureRelaxation, m_pressureTolerance, Nektar::SolverUtils::EquationSystem::m_session, Nektar::IncNavierStokes::m_velocity, m_verbose, m_viscousRelaxation, m_viscousTolerance, NEKERROR, and Nektar::VelocityCorrectionScheme::v_InitObject().

◆ v_SetUpPressureForcing()

void Nektar::VCSMapping::v_SetUpPressureForcing ( const Array< OneD, const Array< OneD, NekDouble > > &  fields,
Array< OneD, Array< OneD, NekDouble > > &  Forcing,
const NekDouble  aii_Dt 
)
overrideprotectedvirtual

Forcing term for Poisson solver solver

Reimplemented from Nektar::VelocityCorrectionScheme.

Definition at line 217 of file VCSMapping.cpp.

220{
221 if (m_mapping->HasConstantJacobian())
222 {
224 aii_Dt);
225 }
226 else
227 {
228 size_t physTot = m_fields[0]->GetTotPoints();
229 size_t nvel = m_nConvectiveFields;
230 Array<OneD, NekDouble> wk(physTot, 0.0);
231
232 Array<OneD, NekDouble> Jac(physTot, 0.0);
233 m_mapping->GetJacobian(Jac);
234
235 // Calculate div(J*u/Dt)
236 Vmath::Zero(physTot, Forcing[0], 1);
237 for (size_t i = 0; i < nvel; ++i)
238 {
239 if (m_fields[i]->GetWaveSpace())
240 {
241 m_fields[i]->HomogeneousBwdTrans(physTot, fields[i], wk);
242 }
243 else
244 {
245 Vmath::Vcopy(physTot, fields[i], 1, wk, 1);
246 }
247 Vmath::Vmul(physTot, wk, 1, Jac, 1, wk, 1);
248 if (m_fields[i]->GetWaveSpace())
249 {
250 m_fields[i]->HomogeneousFwdTrans(physTot, wk, wk);
251 }
252 m_fields[i]->PhysDeriv(MultiRegions::DirCartesianMap[i], wk, wk);
253 Vmath::Vadd(physTot, wk, 1, Forcing[0], 1, Forcing[0], 1);
254 }
255 Vmath::Smul(physTot, 1.0 / aii_Dt, Forcing[0], 1, Forcing[0], 1);
256
257 //
258 // If the mapping viscous terms are being treated explicitly
259 // we need to apply a correction to the forcing
261 {
262 bool wavespace = m_fields[0]->GetWaveSpace();
263 m_fields[0]->SetWaveSpace(false);
264
265 //
266 // Part 1: div(J*grad(U/J . grad(J)))
267 Array<OneD, Array<OneD, NekDouble>> tmp(nvel);
268 Array<OneD, Array<OneD, NekDouble>> velocity(nvel);
269 for (size_t i = 0; i < tmp.size(); i++)
270 {
271 tmp[i] = Array<OneD, NekDouble>(physTot, 0.0);
272 velocity[i] = Array<OneD, NekDouble>(physTot, 0.0);
273 if (wavespace)
274 {
275 m_fields[0]->HomogeneousBwdTrans(
276 physTot, m_fields[i]->GetPhys(), velocity[i]);
277 }
278 else
279 {
280 Vmath::Vcopy(physTot, m_fields[i]->GetPhys(), 1,
281 velocity[i], 1);
282 }
283 }
284 // Calculate wk = U.grad(J)
285 m_mapping->DotGradJacobian(velocity, wk);
286 // Calculate wk = (U.grad(J))/J
287 Vmath::Vdiv(physTot, wk, 1, Jac, 1, wk, 1);
288 // J*grad[(U.grad(J))/J]
289 for (size_t i = 0; i < nvel; ++i)
290 {
291 m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[i], wk,
292 tmp[i]);
293 Vmath::Vmul(physTot, Jac, 1, tmp[i], 1, tmp[i], 1);
294 }
295 // div(J*grad[(U.grad(J))/J])
296 Vmath::Zero(physTot, wk, 1);
297 for (size_t i = 0; i < nvel; ++i)
298 {
299 m_fields[0]->PhysDeriv(MultiRegions::DirCartesianMap[i], tmp[i],
300 tmp[i]);
301 Vmath::Vadd(physTot, wk, 1, tmp[i], 1, wk, 1);
302 }
303
304 // Part 2: grad(J) . curl(curl(U))
305 m_mapping->CurlCurlField(velocity, tmp, m_implicitViscous);
306 // dont need velocity any more, so reuse it
307 m_mapping->DotGradJacobian(tmp, velocity[0]);
308
309 // Add two parts
310 Vmath::Vadd(physTot, velocity[0], 1, wk, 1, wk, 1);
311
312 // Multiply by kinvis and prepare to extrapolate
313 size_t nlevels = m_presForcingCorrection.size();
314 Vmath::Smul(physTot, m_kinvis, wk, 1,
315 m_presForcingCorrection[nlevels - 1], 1);
316
317 // Extrapolate correction
318 m_extrapolation->ExtrapolateArray(m_presForcingCorrection);
319
320 // Put in wavespace
321 if (wavespace)
322 {
323 m_fields[0]->HomogeneousFwdTrans(
324 physTot, m_presForcingCorrection[nlevels - 1], wk);
325 }
326 else
327 {
328 Vmath::Vcopy(physTot, m_presForcingCorrection[nlevels - 1], 1,
329 wk, 1);
330 }
331 // Apply correction: Forcing = Forcing - correction
332 Vmath::Vsub(physTot, Forcing[0], 1, wk, 1, Forcing[0], 1);
333
334 m_fields[0]->SetWaveSpace(wavespace);
335 }
336 }
337}
virtual void v_SetUpPressureForcing(const Array< OneD, const Array< OneD, NekDouble > > &fields, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt)
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:207

References Nektar::MultiRegions::DirCartesianMap, Nektar::IncNavierStokes::m_extrapolation, Nektar::SolverUtils::EquationSystem::m_fields, m_implicitViscous, Nektar::IncNavierStokes::m_kinvis, m_mapping, Nektar::IncNavierStokes::m_nConvectiveFields, m_presForcingCorrection, Vmath::Smul(), Nektar::VelocityCorrectionScheme::v_SetUpPressureForcing(), Vmath::Vadd(), Vmath::Vcopy(), Vmath::Vdiv(), Vmath::Vmul(), Vmath::Vsub(), and Vmath::Zero().

◆ v_SetUpViscousForcing()

void Nektar::VCSMapping::v_SetUpViscousForcing ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  Forcing,
const NekDouble  aii_Dt 
)
overrideprotectedvirtual

Forcing term for Helmholtz solver

Reimplemented from Nektar::VelocityCorrectionScheme.

Definition at line 342 of file VCSMapping.cpp.

345{
346 NekDouble aii_dtinv = 1.0 / aii_Dt;
347 size_t physTot = m_fields[0]->GetTotPoints();
348
349 // Grad p
350 m_pressure->BwdTrans(m_pressure->GetCoeffs(), m_pressure->UpdatePhys());
351
352 size_t nvel = m_velocity.size();
353 if (nvel == 2)
354 {
355 m_pressure->PhysDeriv(m_pressure->GetPhys(), Forcing[0], Forcing[1]);
356 }
357 else
358 {
359 m_pressure->PhysDeriv(m_pressure->GetPhys(), Forcing[0], Forcing[1],
360 Forcing[2]);
361 }
362
363 // Copy grad p in physical space to m_gradP to reuse later
364 if (m_pressure->GetWaveSpace())
365 {
366 for (size_t i = 0; i < nvel; i++)
367 {
368 m_pressure->HomogeneousBwdTrans(physTot, Forcing[i], m_gradP[i]);
369 }
370 }
371 else
372 {
373 for (size_t i = 0; i < nvel; i++)
374 {
375 Vmath::Vcopy(physTot, Forcing[i], 1, m_gradP[i], 1);
376 }
377 }
378
379 if ((!m_mapping->HasConstantJacobian()) || m_implicitPressure)
380 {
381 // If pressure terms are treated explicitly, we need to divide by J
382 // if they are implicit, we need to calculate G(p)
384 {
385 m_mapping->RaiseIndex(m_gradP, Forcing);
386 }
387 else
388 {
389 Array<OneD, NekDouble> Jac(physTot, 0.0);
390 m_mapping->GetJacobian(Jac);
391 for (size_t i = 0; i < nvel; i++)
392 {
393 Vmath::Vdiv(physTot, m_gradP[i], 1, Jac, 1, Forcing[i], 1);
394 }
395 }
396 // Transform back to wavespace
397 if (m_pressure->GetWaveSpace())
398 {
399 for (size_t i = 0; i < nvel; i++)
400 {
401 m_pressure->HomogeneousFwdTrans(physTot, Forcing[i],
402 Forcing[i]);
403 }
404 }
405 }
406
407 // Subtract inarray/(aii_dt) and divide by kinvis. Kinvis will
408 // need to be updated for the convected fields.
409 for (size_t i = 0; i < nvel; ++i)
410 {
411 Blas::Daxpy(physTot, -aii_dtinv, inarray[i], 1, Forcing[i], 1);
412 Blas::Dscal(physTot, 1.0 / m_kinvis, &(Forcing[i])[0], 1);
413 }
414}
static void Dscal(const int &n, const double &alpha, double *x, const int &incx)
BLAS level 1: x = alpha x.
Definition: Blas.hpp:151
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:137
double NekDouble

References Blas::Daxpy(), Blas::Dscal(), Nektar::SolverUtils::EquationSystem::m_fields, m_gradP, m_implicitPressure, Nektar::IncNavierStokes::m_kinvis, m_mapping, Nektar::IncNavierStokes::m_pressure, Nektar::IncNavierStokes::m_velocity, Vmath::Vcopy(), and Vmath::Vdiv().

◆ v_SolvePressure()

void Nektar::VCSMapping::v_SolvePressure ( const Array< OneD, NekDouble > &  Forcing)
overrideprotectedvirtual

Solve pressure system

Reimplemented from Nektar::VelocityCorrectionScheme.

Definition at line 419 of file VCSMapping.cpp.

420{
422 {
424 }
425 else
426 {
427 size_t physTot = m_fields[0]->GetTotPoints();
428 size_t nvel = m_nConvectiveFields;
429 bool converged = false; // flag to mark if system converged
430 size_t s = 0; // iteration counter
431 NekDouble error; // L2 error at current iteration
432 NekDouble forcing_L2 = 0.0; // L2 norm of F
433
434 size_t maxIter;
435 m_session->LoadParameter("MappingMaxIter", maxIter, 5000);
436
437 // rhs of the equation at current iteration
438 Array<OneD, NekDouble> F_corrected(physTot, 0.0);
439 // Pressure field at previous iteration
440 Array<OneD, NekDouble> previous_iter(physTot, 0.0);
441 // Temporary variables
442 Array<OneD, Array<OneD, NekDouble>> wk1(nvel);
443 Array<OneD, Array<OneD, NekDouble>> wk2(nvel);
444 Array<OneD, Array<OneD, NekDouble>> gradP(nvel);
445 for (size_t i = 0; i < nvel; ++i)
446 {
447 wk1[i] = Array<OneD, NekDouble>(physTot, 0.0);
448 wk2[i] = Array<OneD, NekDouble>(physTot, 0.0);
449 gradP[i] = Array<OneD, NekDouble>(physTot, 0.0);
450 }
451
452 // Jacobian
453 Array<OneD, NekDouble> Jac(physTot, 0.0);
454 m_mapping->GetJacobian(Jac);
455
456 // Factors for Laplacian system
459
460 m_pressure->BwdTrans(m_pressure->GetCoeffs(), m_pressure->UpdatePhys());
461 forcing_L2 = m_pressure->L2(Forcing, wk1[0]);
462 while (!converged)
463 {
464 // Update iteration counter and set previous iteration field
465 // (use previous timestep solution for first iteration)
466 s++;
467 ASSERTL0(s < maxIter,
468 "VCSMapping exceeded maximum number of iterations.");
469
470 Vmath::Vcopy(physTot, m_pressure->GetPhys(), 1, previous_iter, 1);
471
472 // Correct pressure bc to account for iteration
473 m_extrapolation->CorrectPressureBCs(previous_iter);
474
475 //
476 // Calculate forcing term for this iteration
477 //
478 for (size_t i = 0; i < nvel; ++i)
479 {
481 previous_iter, gradP[i]);
482 if (m_pressure->GetWaveSpace())
483 {
484 m_pressure->HomogeneousBwdTrans(physTot, gradP[i], wk1[i]);
485 }
486 else
487 {
488 Vmath::Vcopy(physTot, gradP[i], 1, wk1[i], 1);
489 }
490 }
491 m_mapping->RaiseIndex(wk1, wk2); // G(p)
492
493 m_mapping->Divergence(wk2, F_corrected); // div(G(p))
494 if (!m_mapping->HasConstantJacobian())
495 {
496 Vmath::Vmul(physTot, F_corrected, 1, Jac, 1, F_corrected, 1);
497 }
498 // alpha*J*div(G(p))
499 Vmath::Smul(physTot, m_pressureRelaxation, F_corrected, 1,
500 F_corrected, 1);
501 if (m_pressure->GetWaveSpace())
502 {
503 m_pressure->HomogeneousFwdTrans(physTot, F_corrected,
504 F_corrected);
505 }
506 // alpha*J*div(G(p)) - p_ii
507 for (int i = 0; i < m_nConvectiveFields; ++i)
508 {
510 gradP[i], wk1[0]);
511 Vmath::Vsub(physTot, F_corrected, 1, wk1[0], 1, F_corrected, 1);
512 }
513 // p_i,i - J*div(G(p))
514 Vmath::Neg(physTot, F_corrected, 1);
515 // alpha*F - alpha*J*div(G(p)) + p_i,i
516 Vmath::Smul(physTot, m_pressureRelaxation, Forcing, 1, wk1[0], 1);
517 Vmath::Vadd(physTot, wk1[0], 1, F_corrected, 1, F_corrected, 1);
518
519 //
520 // Solve system
521 //
522 m_pressure->HelmSolve(F_corrected, m_pressure->UpdateCoeffs(),
523 factors);
524 m_pressure->BwdTrans(m_pressure->GetCoeffs(),
525 m_pressure->UpdatePhys());
526
527 //
528 // Test convergence
529 //
530 error = m_pressure->L2(m_pressure->GetPhys(), previous_iter);
531 if (forcing_L2 != 0)
532 {
533 if ((error / forcing_L2 < m_pressureTolerance))
534 {
535 converged = true;
536 }
537 }
538 else
539 {
540 if (error < m_pressureTolerance)
541 {
542 converged = true;
543 }
544 }
545 }
546 if (m_verbose && m_session->GetComm()->GetRank() == 0)
547 {
548 std::cout << " Pressure system (mapping) converged in " << s
549 << " iterations with error = " << error << std::endl;
550 }
551 }
552}
virtual void v_SolvePressure(const Array< OneD, NekDouble > &Forcing)
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:408
StdRegions::ConstFactorMap factors

References ASSERTL0, Nektar::MultiRegions::DirCartesianMap, Nektar::StdRegions::eFactorLambda, Nektar::VarcoeffHashingTest::factors, Nektar::IncNavierStokes::m_extrapolation, Nektar::SolverUtils::EquationSystem::m_fields, m_implicitPressure, m_mapping, Nektar::IncNavierStokes::m_nConvectiveFields, Nektar::IncNavierStokes::m_pressure, m_pressureRelaxation, m_pressureTolerance, Nektar::SolverUtils::EquationSystem::m_session, m_verbose, Vmath::Neg(), Vmath::Smul(), Nektar::VelocityCorrectionScheme::v_SolvePressure(), Vmath::Vadd(), Vmath::Vcopy(), Vmath::Vmul(), and Vmath::Vsub().

◆ v_SolveViscous()

void Nektar::VCSMapping::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 
)
overrideprotectedvirtual

Solve velocity system

Reimplemented from Nektar::VelocityCorrectionScheme.

Definition at line 557 of file VCSMapping.cpp.

561{
562 boost::ignore_unused(inarray);
563
565 {
566 VelocityCorrectionScheme::v_SolveViscous(Forcing, inarray, outarray,
567 aii_Dt);
568 }
569 else
570 {
571 size_t physTot = m_fields[0]->GetTotPoints();
572 size_t nvel = m_nConvectiveFields;
573 bool converged = false; // flag to mark if system converged
574 size_t s = 0; // iteration counter
575 NekDouble error, max_error; // L2 error at current iteration
576
577 size_t maxIter;
578 m_session->LoadParameter("MappingMaxIter", maxIter, 5000);
579
580 // L2 norm of F
581 Array<OneD, NekDouble> forcing_L2(m_nConvectiveFields, 0.0);
582
583 // rhs of the equation at current iteration
584 Array<OneD, Array<OneD, NekDouble>> F_corrected(nvel);
585 // Solution at previous iteration
586 Array<OneD, Array<OneD, NekDouble>> previous_iter(nvel);
587 // Working space
588 Array<OneD, Array<OneD, NekDouble>> wk(nvel);
589 for (size_t i = 0; i < nvel; ++i)
590 {
591 F_corrected[i] = Array<OneD, NekDouble>(physTot, 0.0);
592 previous_iter[i] = Array<OneD, NekDouble>(physTot, 0.0);
593 wk[i] = Array<OneD, NekDouble>(physTot, 0.0);
594 }
595
596 // Factors for Helmholtz system
599 1.0 * m_viscousRelaxation / aii_Dt / m_kinvis;
601 {
605 }
606
607 // Calculate L2-norm of F and set initial solution for iteration
608 for (size_t i = 0; i < nvel; ++i)
609 {
610 forcing_L2[i] = m_fields[0]->L2(Forcing[i], wk[0]);
611 m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(), previous_iter[i]);
612 }
613
614 while (!converged)
615 {
616 converged = true;
617 // Iteration counter
618 s++;
619 ASSERTL0(s < maxIter,
620 "VCSMapping exceeded maximum number of iterations.");
621
622 max_error = 0.0;
623
624 //
625 // Calculate forcing term for next iteration
626 //
627
628 // Calculate L(U)- in this parts all components might be coupled
629 if (m_fields[0]->GetWaveSpace())
630 {
631 for (size_t i = 0; i < nvel; ++i)
632 {
633 m_fields[0]->HomogeneousBwdTrans(physTot, previous_iter[i],
634 wk[i]);
635 }
636 }
637 else
638 {
639 for (size_t i = 0; i < nvel; ++i)
640 {
641 Vmath::Vcopy(physTot, previous_iter[i], 1, wk[i], 1);
642 }
643 }
644
645 // (L(U^i) - 1/alpha*U^i_jj)
646 m_mapping->VelocityLaplacian(wk, F_corrected,
647 1.0 / m_viscousRelaxation);
648
649 if (m_fields[0]->GetWaveSpace())
650 {
651 for (size_t i = 0; i < nvel; ++i)
652 {
653 m_fields[0]->HomogeneousFwdTrans(physTot, F_corrected[i],
654 F_corrected[i]);
655 }
656 }
657 else
658 {
659 for (size_t i = 0; i < nvel; ++i)
660 {
661 Vmath::Vcopy(physTot, F_corrected[i], 1, F_corrected[i], 1);
662 }
663 }
664
665 // Loop velocity components
666 for (size_t i = 0; i < nvel; ++i)
667 {
668 // (-alpha*L(U^i) + U^i_jj)
669 Vmath::Smul(physTot, -1.0 * m_viscousRelaxation, F_corrected[i],
670 1, F_corrected[i], 1);
671 // F_corrected = alpha*F + (-alpha*L(U^i) + U^i_jj)
672 Vmath::Smul(physTot, m_viscousRelaxation, Forcing[i], 1, wk[0],
673 1);
674 Vmath::Vadd(physTot, wk[0], 1, F_corrected[i], 1,
675 F_corrected[i], 1);
676
677 //
678 // Solve System
679 //
680 m_fields[i]->HelmSolve(F_corrected[i],
681 m_fields[i]->UpdateCoeffs(), factors);
682 m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(), outarray[i]);
683
684 //
685 // Test convergence
686 //
687 error = m_fields[i]->L2(outarray[i], previous_iter[i]);
688
689 if (forcing_L2[i] != 0)
690 {
691 if ((error / forcing_L2[i] >= m_viscousTolerance))
692 {
693 converged = false;
694 }
695 }
696 else
697 {
698 if (error >= m_viscousTolerance)
699 {
700 converged = false;
701 }
702 }
703 if (error > max_error)
704 {
705 max_error = error;
706 }
707
708 // Copy field to previous_iter
709 Vmath::Vcopy(physTot, outarray[i], 1, previous_iter[i], 1);
710 }
711 }
712 if (m_verbose && m_session->GetComm()->GetRank() == 0)
713 {
714 std::cout << " Velocity system (mapping) converged in " << s
715 << " iterations with error = " << max_error << std::endl;
716 }
717 }
718}
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.
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 ASSERTL0, Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::eFactorSVVCutoffRatio, Nektar::StdRegions::eFactorSVVDiffCoeff, Nektar::VarcoeffHashingTest::factors, Nektar::SolverUtils::EquationSystem::m_fields, m_implicitViscous, Nektar::IncNavierStokes::m_kinvis, m_mapping, Nektar::IncNavierStokes::m_nConvectiveFields, Nektar::SolverUtils::EquationSystem::m_session, Nektar::VelocityCorrectionScheme::m_sVVCutoffRatio, Nektar::VelocityCorrectionScheme::m_sVVDiffCoeff, Nektar::VelocityCorrectionScheme::m_useSpecVanVisc, m_verbose, m_viscousRelaxation, m_viscousTolerance, Vmath::Smul(), Nektar::VelocityCorrectionScheme::v_SolveViscous(), Vmath::Vadd(), and Vmath::Vcopy().

Member Data Documentation

◆ className

string Nektar::VCSMapping::className
static
Initial value:
=
"VCSMapping", VCSMapping::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.
Definition: VCSMapping.h:47
EquationSystemFactory & GetEquationSystemFactory()

Name of class.

Definition at line 58 of file VCSMapping.h.

◆ m_gradP

Array<OneD, Array<OneD, NekDouble> > Nektar::VCSMapping::m_gradP
protected

Definition at line 94 of file VCSMapping.h.

Referenced by MappingPressureCorrection(), v_DoInitialise(), and v_SetUpViscousForcing().

◆ m_implicitPressure

bool Nektar::VCSMapping::m_implicitPressure
protected

◆ m_implicitViscous

bool Nektar::VCSMapping::m_implicitViscous
protected

◆ m_mapping

GlobalMapping::MappingSharedPtr Nektar::VCSMapping::m_mapping
protected

◆ m_neglectViscous

bool Nektar::VCSMapping::m_neglectViscous
protected

Definition at line 85 of file VCSMapping.h.

Referenced by ApplyIncNSMappingForcing(), and v_InitObject().

◆ m_presForcingCorrection

Array<OneD, Array<OneD, NekDouble> > Nektar::VCSMapping::m_presForcingCorrection
private

Definition at line 124 of file VCSMapping.h.

Referenced by v_InitObject(), and v_SetUpPressureForcing().

◆ m_pressureRelaxation

NekDouble Nektar::VCSMapping::m_pressureRelaxation
protected

Definition at line 90 of file VCSMapping.h.

Referenced by v_InitObject(), and v_SolvePressure().

◆ m_pressureTolerance

NekDouble Nektar::VCSMapping::m_pressureTolerance
protected

Definition at line 88 of file VCSMapping.h.

Referenced by v_InitObject(), and v_SolvePressure().

◆ m_verbose

bool Nektar::VCSMapping::m_verbose
protected

Definition at line 79 of file VCSMapping.h.

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

◆ m_viscousRelaxation

NekDouble Nektar::VCSMapping::m_viscousRelaxation
protected

Definition at line 91 of file VCSMapping.h.

Referenced by v_InitObject(), and v_SolveViscous().

◆ m_viscousTolerance

NekDouble Nektar::VCSMapping::m_viscousTolerance
protected

Definition at line 89 of file VCSMapping.h.

Referenced by v_InitObject(), and v_SolveViscous().

◆ solverTypeLookupId

string Nektar::VCSMapping::solverTypeLookupId
staticprotected
Initial value:
=
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.

Definition at line 77 of file VCSMapping.h.