Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | 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 ()
 
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)
 
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 ()
 
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, 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 int GetPararealIterationNumber ()
 
SOLVER_UTILS_EXPORT void SetPararealIterationNumber (int num)
 
SOLVER_UTILS_EXPORT bool GetUseInitialCondition ()
 
SOLVER_UTILS_EXPORT void SetUseInitialCondition (bool 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...
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp ()
 Virtual function to 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 (void) 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 Array< OneD, bool > v_GetSystemSingularChecks () override
 
virtual int v_GetForceDimension () override
 
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)
 
- 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 () override
 Solves an unsteady problem. 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 v_UpdateTimeStepCheck ()
 
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 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 Member Functions inherited from Nektar::SolverUtils::FluidInterface
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)
 

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
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
 
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_infosteps
 Number of time steps between outputting status information. More...
 
int m_pararealIter
 Number of parareal 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_useInitialCondition
 Flag to determine if IC are used. 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...
 

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

- 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 43 of file VCSMapping.h.

Constructor & Destructor Documentation

◆ VCSMapping()

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

Constructor.

Constructor. Creates ...

Parameters

Definition at line 55 of file VCSMapping.cpp.

57  : UnsteadySystem(pSession, pGraph),
58  VelocityCorrectionScheme(pSession, pGraph)
59 {
60 }
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 127 of file VCSMapping.cpp.

128 {
129 }

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 719 of file VCSMapping.cpp.

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

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 827 of file VCSMapping.cpp.

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

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 802 of file VCSMapping.cpp.

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

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 906 of file VCSMapping.cpp.

908 {
909  size_t physTot = m_fields[0]->GetTotPoints();
910  size_t nvel = m_nConvectiveFields;
911 
912  // Calculate g^(ij)p_(,j)
913  m_mapping->RaiseIndex(m_gradP, outarray);
914 
915  // Calculate correction = (nabla p)/J - g^(ij)p_,j
916  // (Jac is not required if it is constant)
917  if (!m_mapping->HasConstantJacobian())
918  {
919  Array<OneD, NekDouble> Jac(physTot, 0.0);
920  m_mapping->GetJacobian(Jac);
921  for (size_t i = 0; i < nvel; ++i)
922  {
923  Vmath::Vdiv(physTot, m_gradP[i], 1, Jac, 1, m_gradP[i], 1);
924  }
925  }
926  for (size_t i = 0; i < nvel; ++i)
927  {
928  Vmath::Vsub(physTot, m_gradP[i], 1, outarray[i], 1, outarray[i], 1);
929  }
930 }
Array< OneD, Array< OneD, NekDouble > > m_gradP
Definition: VCSMapping.h:92
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:284
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:419

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 932 of file VCSMapping.cpp.

935 {
936  // L(U) - 1.0*d^2(u^i)/dx^jdx^j
937  m_mapping->VelocityLaplacian(velPhys, outarray, 1.0);
938 }

References m_mapping.

Referenced by ApplyIncNSMappingForcing().

◆ v_DoInitialise()

void Nektar::VCSMapping::v_DoInitialise ( void  )
overrideprotectedvirtual

Sets up initial conditions.

Sets the initial conditions.

Reimplemented from Nektar::VelocityCorrectionScheme.

Definition at line 131 of file VCSMapping.cpp.

132 {
134 
135  // Set up Field Meta Data for output files
136  m_fieldMetaDataMap["Kinvis"] = boost::lexical_cast<std::string>(m_kinvis);
137  m_fieldMetaDataMap["TimeStep"] =
138  boost::lexical_cast<std::string>(m_timestep);
139 
140  // Correct Dirichlet boundary conditions to account for mapping
141  m_mapping->UpdateBCs(0.0);
142  //
143  m_F = Array<OneD, Array<OneD, NekDouble>>(m_nConvectiveFields);
144  for (int i = 0; i < m_nConvectiveFields; ++i)
145  {
146  m_fields[i]->ImposeDirichletConditions(m_fields[i]->UpdateCoeffs());
147  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
148  m_fields[i]->UpdatePhys());
149  m_F[i] = Array<OneD, NekDouble>(m_fields[0]->GetTotPoints(), 0.0);
150  }
151 
152  // Initialise m_gradP
153  size_t physTot = m_fields[0]->GetTotPoints();
154  m_gradP = Array<OneD, Array<OneD, NekDouble>>(m_nConvectiveFields);
155  for (int i = 0; i < m_nConvectiveFields; ++i)
156  {
157  m_gradP[i] = Array<OneD, NekDouble>(physTot, 0.0);
159  m_pressure->GetPhys(), m_gradP[i]);
160  if (m_pressure->GetWaveSpace())
161  {
162  m_pressure->HomogeneousBwdTrans(physTot, m_gradP[i], m_gradP[i]);
163  }
164  }
165 }
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() 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 170 of file VCSMapping.cpp.

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

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 62 of file VCSMapping.cpp.

63 {
65 
67  ASSERTL0(m_mapping, "Could not create mapping in VCSMapping.");
68 
69  std::string vExtrapolation = "Mapping";
71  vExtrapolation, m_session, m_fields, m_pressure, m_velocity,
72  m_advObject);
73  m_extrapolation->SubSteppingTimeIntegration(m_intScheme);
74  m_extrapolation->GenerateHOPBCMap(m_session);
75 
76  // Storage to extrapolate pressure forcing
77  size_t physTot = m_fields[0]->GetTotPoints();
78  size_t intSteps = 1;
79 
80  if (m_intScheme->GetName() == "IMEX" ||
81  m_intScheme->GetName() == "IMEXGear")
82  {
83  m_intSteps = m_intScheme->GetOrder();
84  }
85  else
86  {
88  "Integration method not suitable: "
89  "Options include IMEXGear or IMEXOrder{1,2,3,4}");
90  }
91 
92  m_presForcingCorrection = Array<OneD, Array<OneD, NekDouble>>(intSteps);
93  for (size_t i = 0; i < m_presForcingCorrection.size(); i++)
94  {
95  m_presForcingCorrection[i] = Array<OneD, NekDouble>(physTot, 0.0);
96  }
97  m_verbose = (m_session->DefinesCmdLineArgument("verbose")) ? true : false;
98 
99  // Load solve parameters related to the mapping
100  // Flags determining if pressure/viscous terms should be treated implicitly
101  m_session->MatchSolverInfo("MappingImplicitPressure", "True",
102  m_implicitPressure, false);
103  m_session->MatchSolverInfo("MappingImplicitViscous", "True",
104  m_implicitViscous, false);
105  m_session->MatchSolverInfo("MappingNeglectViscous", "True",
106  m_neglectViscous, false);
107 
108  if (m_neglectViscous)
109  {
110  m_implicitViscous = false;
111  }
112 
113  // Tolerances and relaxation parameters for implicit terms
114  m_session->LoadParameter("MappingPressureTolerance", m_pressureTolerance,
115  1e-12);
116  m_session->LoadParameter("MappingViscousTolerance", m_viscousTolerance,
117  1e-12);
118  m_session->LoadParameter("MappingPressureRelaxation", m_pressureRelaxation,
119  1.0);
120  m_session->LoadParameter("MappingViscousRelaxation", m_viscousRelaxation,
121  1.0);
122 }
#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:270
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:86
NekDouble m_viscousRelaxation
Definition: VCSMapping.h:89
NekDouble m_pressureRelaxation
Definition: VCSMapping.h:88
NekDouble m_viscousTolerance
Definition: VCSMapping.h:87
Array< OneD, Array< OneD, NekDouble > > m_presForcingCorrection
Definition: VCSMapping.h:122
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 213 of file VCSMapping.cpp.

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

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 338 of file VCSMapping.cpp.

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

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

References ASSERTL0, Nektar::MultiRegions::DirCartesianMap, Nektar::StdRegions::eFactorLambda, 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 553 of file VCSMapping.cpp.

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

Referenced by v_InitObject(), and v_SetUpPressureForcing().

◆ m_pressureRelaxation

NekDouble Nektar::VCSMapping::m_pressureRelaxation
protected

Definition at line 88 of file VCSMapping.h.

Referenced by v_InitObject(), and v_SolvePressure().

◆ m_pressureTolerance

NekDouble Nektar::VCSMapping::m_pressureTolerance
protected

Definition at line 86 of file VCSMapping.h.

Referenced by v_InitObject(), and v_SolvePressure().

◆ m_verbose

bool Nektar::VCSMapping::m_verbose
protected

Definition at line 77 of file VCSMapping.h.

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

◆ m_viscousRelaxation

NekDouble Nektar::VCSMapping::m_viscousRelaxation
protected

Definition at line 89 of file VCSMapping.h.

Referenced by v_InitObject(), and v_SolveViscous().

◆ m_viscousTolerance

NekDouble Nektar::VCSMapping::m_viscousTolerance
protected

Definition at line 87 of file VCSMapping.h.

Referenced by v_InitObject(), and v_SolveViscous().