Nektar++
Loading...
Searching...
No Matches
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | Static Protected Attributes | Private Member Functions | Private Attributes | Friends | List of all members
Nektar::VCSMapping Class Reference

#include <VCSMapping.h>

Inheritance diagram for Nektar::VCSMapping:
[legend]

Public Member Functions

void ApplyIncNSMappingForcing (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
- Public Member Functions inherited from Nektar::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
int GetNConvectiveFields (void)
 
void AddForcing (const SolverUtils::ForcingSharedPtr &pForce)
 
bool DefinedForcing (const std::string &sForce)
 
- Public Member Functions inherited from Nektar::SolverUtils::AdvectionSystem
SOLVER_UTILS_EXPORT AdvectionSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
SOLVER_UTILS_EXPORT ~AdvectionSystem () override=default
 
SOLVER_UTILS_EXPORT AdvectionSharedPtr GetAdvObject ()
 Returns the advection object held by this instance.
 
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
SOLVER_UTILS_EXPORT ~UnsteadySystem () override=default
 Destructor.
 
SOLVER_UTILS_EXPORT NekDouble GetTimeStep (const Array< OneD, const Array< OneD, NekDouble > > &inarray)
 Calculate the larger time-step mantaining the problem stable.
 
SOLVER_UTILS_EXPORT NekDouble GetTimeStep ()
 
SOLVER_UTILS_EXPORT void SetTimeStep (const NekDouble timestep)
 
SOLVER_UTILS_EXPORT void SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
SOLVER_UTILS_EXPORT LibUtilities::TimeIntegrationSchemeSharedPtrGetTimeIntegrationScheme ()
 Returns the time integration scheme.
 
SOLVER_UTILS_EXPORT LibUtilities::TimeIntegrationSchemeOperatorsGetTimeIntegrationSchemeOperators ()
 Returns the time integration scheme operators.
 
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor.
 
SOLVER_UTILS_EXPORT void InitObject (bool DeclareField=true)
 Initialises the members of this object.
 
SOLVER_UTILS_EXPORT void DoInitialise (bool dumpInitialConditions=true)
 Perform any initialisation necessary before solving the problem.
 
SOLVER_UTILS_EXPORT void DoSolve ()
 Solve the problem.
 
SOLVER_UTILS_EXPORT void TransCoeffToPhys ()
 Transform from coefficient to physical space.
 
SOLVER_UTILS_EXPORT void TransPhysToCoeff ()
 Transform from physical to coefficient space.
 
SOLVER_UTILS_EXPORT void Output ()
 Perform output operations after solve.
 
SOLVER_UTILS_EXPORT std::string GetSessionName ()
 Get Session name.
 
template<class T >
std::shared_ptr< T > as ()
 
SOLVER_UTILS_EXPORT void ResetSessionName (std::string newname)
 Reset Session name.
 
SOLVER_UTILS_EXPORT LibUtilities::SessionReaderSharedPtr GetSession ()
 Get Session name.
 
SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr GetPressure ()
 Get pressure field if available.
 
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.
 
SOLVER_UTILS_EXPORT void SetLambda (NekDouble lambda)
 Set parameter m_lambda.
 
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction (std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
 Get a SessionFunction by name.
 
SOLVER_UTILS_EXPORT void SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 Initialise the data in the dependent fields.
 
SOLVER_UTILS_EXPORT void EvaluateExactSolution (int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 Evaluates an exact solution.
 
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.
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, bool Normalised=false)
 Compute the L2 error of the fields.
 
SOLVER_UTILS_EXPORT NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation.
 
SOLVER_UTILS_EXPORT NekDouble H1Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln, bool Normalised=false)
 Compute the H1 error between fields and a given exact solution.
 
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].
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n)
 Write checkpoint file of m_fields.
 
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.
 
SOLVER_UTILS_EXPORT void Checkpoint_BaseFlow (const int n)
 Write base flow file of m_fields.
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname)
 Write field data to the given filename.
 
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.
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 Input field data from the given file.
 
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.
 
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.
 
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.
 
SOLVER_UTILS_EXPORT void SessionSummary (SummaryList &vSummary)
 Write out a session summary.
 
SOLVER_UTILS_EXPORT Array< OneD, MultiRegions::ExpListSharedPtr > & UpdateFields ()
 
SOLVER_UTILS_EXPORT LibUtilities::FieldMetaDataMapUpdateFieldMetaDataMap ()
 Get hold of FieldInfoMap so it can be updated.
 
SOLVER_UTILS_EXPORT NekDouble GetTime ()
 Return final time.
 
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 void SetSteps (const int steps)
 
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 Array< OneD, NekDouble > & UpdatePhysField (const int i)
 
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.
 
SOLVER_UTILS_EXPORT bool NegatedOp ()
 Identify if operator is negated in DoSolve.
 
- Public Member Functions inherited from Nektar::SolverUtils::ALEHelper
virtual ~ALEHelper ()=default
 
virtual SOLVER_UTILS_EXPORT void v_ALEInitObject (int spaceDim, Array< OneD, MultiRegions::ExpListSharedPtr > &fields)
 
SOLVER_UTILS_EXPORT void InitObject (int spaceDim, Array< OneD, MultiRegions::ExpListSharedPtr > &fields)
 
virtual SOLVER_UTILS_EXPORT void v_UpdateGridVelocity (const NekDouble &time)
 
virtual SOLVER_UTILS_EXPORT void v_ALEPreMultiplyMass (Array< OneD, Array< OneD, NekDouble > > &fields)
 
SOLVER_UTILS_EXPORT void ALEDoElmtInvMass (Array< OneD, Array< OneD, NekDouble > > &traceNormals, Array< OneD, Array< OneD, NekDouble > > &fields, NekDouble time)
 Update m_fields with u^n by multiplying by inverse mass matrix. That's then used in e.g. checkpoint output and L^2 error calculation.
 
SOLVER_UTILS_EXPORT void ALEDoElmtInvMassBwdTrans (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
SOLVER_UTILS_EXPORT void MoveMesh (const NekDouble &time, Array< OneD, Array< OneD, NekDouble > > &traceNormals)
 
SOLVER_UTILS_EXPORT void ResetMatricesNormal (Array< OneD, Array< OneD, NekDouble > > &traceNormals)
 
SOLVER_UTILS_EXPORT void UpdateNormalsFlag ()
 
const Array< OneD, const Array< OneD, NekDouble > > & GetGridVelocity ()
 
bool & GetUpdateNormalsFlag ()
 
SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & GetGridVelocityTrace ()
 
SOLVER_UTILS_EXPORT void ExtraFldOutputGridVelocity (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 
SOLVER_UTILS_EXPORT void ExtraFldOutputGrid (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 
- 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.
 
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.
 
SOLVER_UTILS_EXPORT void GetPressure (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &pressure)
 Extract array with pressure from physfield.
 
SOLVER_UTILS_EXPORT void SetMovingFrameVelocities (const Array< OneD, NekDouble > &vFrameVels, const int step)
 
SOLVER_UTILS_EXPORT bool GetMovingFrameVelocities (Array< OneD, NekDouble > &vFrameVels, const int step)
 
SOLVER_UTILS_EXPORT void SetMovingFrameDisp (const Array< OneD, NekDouble > &vFrameDisp, const int step)
 
SOLVER_UTILS_EXPORT void SetMovingFramePivot (const Array< OneD, NekDouble > &vFramePivot)
 
SOLVER_UTILS_EXPORT bool GetMovingFrameDisp (Array< OneD, NekDouble > &vFrameDisp, const int step)
 
SOLVER_UTILS_EXPORT void SetAeroForce (Array< OneD, NekDouble > forces)
 Set aerodynamic force and moment.
 
SOLVER_UTILS_EXPORT void GetAeroForce (Array< OneD, NekDouble > forces)
 Get aerodynamic force and moment.
 

Static Public Member Functions

static SolverUtils::EquationSystemSharedPtr create (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Creates an instance of this class.
 
- 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.
 

Static Public Attributes

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

Protected Member Functions

 VCSMapping (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
 ~VCSMapping () override=default
 
void v_InitObject (bool DeclareField=true) override
 Initialisation object for EquationSystem.
 
void v_DoInitialise (bool dumpInitialConditions=true) override
 Sets up initial conditions.
 
void v_SetUpPressureForcing (const Array< OneD, const Array< OneD, NekDouble > > &fields, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt) override
 
void v_SetUpViscousForcing (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt) override
 
void v_SolvePressure (const Array< OneD, NekDouble > &Forcing) override
 
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
 
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
 VelocityCorrectionScheme (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
 ~VelocityCorrectionScheme () override=default
 
void v_InitObject (bool DeclareField=true) override
 Initialisation object for EquationSystem.
 
void SetupFlowrate (NekDouble aii_dt)
 Set up the Stokes solution used to impose constant flowrate through a boundary.
 
NekDouble MeasureFlowrate (const Array< OneD, Array< OneD, NekDouble > > &inarray)
 Measure the volumetric flow rate through the volumetric flow rate reference surface.
 
bool v_PostIntegrate (int step) override
 
void v_GenerateSummary (SolverUtils::SummaryList &s) override
 Print a summary of time stepping parameters.
 
void v_TransCoeffToPhys (void) override
 Virtual function for transformation to physical space.
 
void v_TransPhysToCoeff (void) override
 Virtual function for transformation to coefficient space.
 
void v_DoInitialise (bool dumpInitialConditions=true) override
 Sets up initial conditions.
 
Array< OneD, bool > v_GetSystemSingularChecks () override
 
int v_GetForceDimension () override
 
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)
 
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.
 
 ~IncNavierStokes () override=default
 
void v_GetPressure (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &pressure) override
 
void v_GetDensity (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &density) override
 
bool v_HasConstantDensity () override
 
void v_GetVelocity (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &velocity) override
 
void v_SetMovingFrameVelocities (const Array< OneD, NekDouble > &vFrameVels, const int step) override
 
bool v_GetMovingFrameVelocities (Array< OneD, NekDouble > &vFrameVels, const int step) override
 
void v_SetMovingFrameDisp (const Array< OneD, NekDouble > &vFrameDisp, const int step) override
 
void v_SetMovingFramePivot (const Array< OneD, NekDouble > &vFramePivot) override
 
bool v_GetMovingFrameDisp (Array< OneD, NekDouble > &vFrameDisp, const int step) override
 
void v_SetAeroForce (Array< OneD, NekDouble > forces) override
 
void v_GetAeroForce (Array< OneD, NekDouble > forces) override
 
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
 
void SetRadiationBoundaryForcing (int fieldid)
 Set Radiation forcing term.
 
void SetZeroNormalVelocity ()
 Set Normal Velocity Component to Zero.
 
void SetWomersleyBoundary (const int fldid, const int bndid)
 Set Womersley Profile if specified.
 
void SetUpWomersley (const int fldid, const int bndid, std::string womstr)
 Set Up Womersley details.
 
MultiRegions::ExpListSharedPtr v_GetPressure () override
 
Array< OneD, NekDoublev_GetMaxStdVelocity (const NekDouble SpeedSoundFactor) override
 
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.
 
SOLVER_UTILS_EXPORT void v_DoSolve () override
 Solves an unsteady problem.
 
virtual SOLVER_UTILS_EXPORT void v_PrintStatusInformation (const int step, const NekDouble cpuTime)
 Print Status Information.
 
virtual SOLVER_UTILS_EXPORT void v_PrintSummaryStatistics (const NekDouble intTime)
 Print Summary Statistics.
 
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.
 
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.
 
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.
 
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.
 
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises EquationSystem class members.
 
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.
 
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.
 
virtual SOLVER_UTILS_EXPORT NekDouble v_H1Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false)
 Virtual function for the H_1 error computation between fields and a given exact solution.
 
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 bool v_NegatedOp (void)
 Virtual function to identify if operator is negated in DoSolve.
 
virtual SOLVER_UTILS_EXPORT void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 

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.
 
bool m_useSpecVanVisc
 bool to identify if spectral vanishing viscosity is active.
 
bool m_useGJPStabilisation
 bool to identify if GJP semi-implicit is active.
 
bool m_useGJPNormalVel
 bool to identify if GJP normal Velocity should be applied in explicit approach
 
NekDouble m_GJPJumpScale
 
NekDouble m_sVVCutoffRatio
 cutt off ratio from which to start decayhing modes
 
NekDouble m_sVVDiffCoeff
 Diffusion coefficient of SVV modes.
 
NekDouble m_sVVCutoffRatioHomo1D
 
NekDouble m_sVVDiffCoeffHomo1D
 Diffusion coefficient of SVV modes in homogeneous 1D Direction.
 
Array< OneD, NekDoublem_svvVarDiffCoeff
 Array of coefficient if power kernel is used in SVV.
 
bool m_IsSVVPowerKernel
 Identifier for Power Kernel otherwise DG kernel.
 
Array< OneD, NekDoublem_diffCoeff
 Diffusion coefficients (will be kinvis for velocities)
 
StdRegions::VarCoeffMap m_varCoeffLap
 Variable Coefficient map for the Laplacian which can be activated as part of SVV or otherwise.
 
NekDouble m_flowrate
 Desired volumetric flowrate.
 
NekDouble m_flowrateArea
 Area of the boundary through which we are measuring the flowrate.
 
bool m_homd1DFlowinPlane
 
NekDouble m_greenFlux
 Flux of the Stokes function solution.
 
NekDouble m_alpha
 Current flowrate correction.
 
int m_flowrateBndID
 Boundary ID of the flowrate reference surface.
 
int m_planeID
 Plane ID for cases with homogeneous expansion.
 
MultiRegions::ExpListSharedPtr m_flowrateBnd
 Flowrate reference surface.
 
Array< OneD, Array< OneD, NekDouble > > m_flowrateStokes
 Stokes solution used to impose flowrate.
 
std::ofstream m_flowrateStream
 Output stream to record flowrate.
 
int m_flowrateSteps
 Interval at which to record flowrate data.
 
NekDouble m_flowrateAiidt
 Value of aii_dt used to compute Stokes flowrate solution.
 
Array< OneD, Array< OneD, NekDouble > > m_F
 
- Protected Attributes inherited from Nektar::IncNavierStokes
ExtrapolateSharedPtr m_extrapolation
 
IncBoundaryConditionsSharedPtr m_IncNavierStokesBCs
 
std::ofstream m_mdlFile
 modal energy file
 
bool m_SmoothAdvection
 bool to identify if advection term smoothing is requested
 
std::vector< SolverUtils::ForcingSharedPtrm_forcing
 Forcing terms.
 
int m_nConvectiveFields
 Number of fields to be convected;.
 
Array< OneD, int > m_velocity
 int which identifies which components of m_fields contains the velocity (u,v,w);
 
MultiRegions::ExpListSharedPtr m_pressure
 Pointer to field holding pressure field.
 
NekDouble m_kinvis
 Kinematic viscosity.
 
int m_energysteps
 dump energy to file at steps time
 
EquationType m_equationType
 equation type;
 
Array< OneD, Array< OneD, int > > m_fieldsBCToElmtID
 Mapping from BCs to Elmt IDs.
 
Array< OneD, Array< OneD, int > > m_fieldsBCToTraceID
 Mapping from BCs to Elmt Edge IDs.
 
Array< OneD, Array< OneD, NekDouble > > m_fieldsRadiationFactor
 RHS Factor for Radiation Condition.
 
int m_intSteps
 Number of time integration steps AND Order of extrapolation for pressure boundary conditions.
 
Array< OneD, NekDoublem_pivotPoint
 pivot point for moving reference frame
 
Array< OneD, NekDoublem_aeroForces
 
std::map< int, std::map< int, WomersleyParamsSharedPtr > > m_womersleyParams
 Womersley parameters if required.
 
- Protected Attributes inherited from Nektar::SolverUtils::AdvectionSystem
SolverUtils::AdvectionSharedPtr m_advObject
 Advection term.
 
- Protected Attributes inherited from Nektar::SolverUtils::UnsteadySystem
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
 Wrapper to the time integration scheme.
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use.
 
Array< OneD, Array< OneD, NekDouble > > m_previousSolution
 Storage for previous solution for steady-state check.
 
std::vector< int > m_intVariables
 
NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1).
 
NekDouble m_CFLGrowth
 CFL growth rate.
 
NekDouble m_CFLEnd
 Maximun cfl in cfl growth.
 
int m_abortSteps
 Number of steps between checks for abort conditions.
 
bool m_explicitDiffusion
 Indicates if explicit or implicit treatment of diffusion is used.
 
bool m_explicitAdvection
 Indicates if explicit or implicit treatment of advection is used.
 
bool m_explicitReaction
 Indicates if explicit or implicit treatment of reaction is used.
 
int m_steadyStateSteps
 Check for steady state at step interval.
 
NekDouble m_steadyStateTol
 Tolerance to which steady state should be evaluated at.
 
int m_filtersInfosteps
 Number of time steps between outputting filters information.
 
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
 
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state.
 
std::ofstream m_errFile
 
NekDouble m_epsilon
 Diffusion coefficient.
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator.
 
bool m_verbose
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader.
 
std::map< std::string, SolverUtils::SessionFunctionSharedPtrm_sessionFunctions
 Map of known SessionFunctions.
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output.
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 Array holding all dependent variables.
 
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
 Pointer to boundary conditions object.
 
SpatialDomains::MeshGraphSharedPtr m_graph
 Pointer to graph defining mesh.
 
std::string m_sessionName
 Name of the session.
 
NekDouble m_time
 Current time of simulation.
 
int m_initialStep
 Number of the step where the simulation should begin.
 
NekDouble m_fintime
 Finish time of the simulation.
 
NekDouble m_timestep
 Time step size.
 
NekDouble m_lambda
 Lambda constant in real system if one required.
 
NekDouble m_checktime
 Time between checkpoints.
 
NekDouble m_lastCheckTime
 
NekDouble m_TimeIncrementFactor
 
int m_nchk
 Number of checkpoints written so far.
 
int m_steps
 Number of steps to take.
 
int m_checksteps
 Number of steps between checkpoints.
 
int m_infosteps
 Number of time steps between outputting status information.
 
int m_iterPIT = 0
 Number of parallel-in-time time iteration.
 
int m_windowPIT = 0
 Index of windows for parallel-in-time time iteration.
 
int m_spacedim
 Spatial dimension (>= expansion dim).
 
int m_expdim
 Expansion dimension.
 
bool m_singleMode
 Flag to determine if single homogeneous mode is used.
 
bool m_halfMode
 Flag to determine if half homogeneous mode is used.
 
bool m_multipleModes
 Flag to determine if use multiple homogenenous modes are used.
 
bool m_useFFT
 Flag to determine if FFT is used for homogeneous transform.
 
bool m_homogen_dealiasing
 Flag to determine if dealiasing is used for homogeneous simulations.
 
bool m_specHP_dealiasing
 Flag to determine if dealisising is usde for the Spectral/hp element discretisation.
 
enum MultiRegions::ProjectionType m_projectionType
 Type of projection; e.g continuous or discontinuous.
 
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
 Array holding trace normals for DG simulations in the forwards direction.
 
Array< OneD, bool > m_checkIfSystemSingular
 Flag to indicate if the fields should be checked for singularity.
 
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
 Map to identify relevant solver info to dump in output fields.
 
Array< OneD, NekDoublem_movingFrameData
 Moving reference frame status in the inertial frame X, Y, Z, Theta_x, Theta_y, Theta_z, U, V, W, Omega_x, Omega_y, Omega_z, A_x, A_y, A_z, DOmega_x, DOmega_y, DOmega_z, pivot_x, pivot_y, pivot_z.
 
std::vector< std::string > m_strFrameData
 variable name in m_movingFrameData
 
int m_NumQuadPointsError
 Number of Quadrature points used to work out the error.
 
enum HomogeneousType m_HomogeneousType
 
NekDouble m_LhomX
 physical length in X direction (if homogeneous)
 
NekDouble m_LhomY
 physical length in Y direction (if homogeneous)
 
NekDouble m_LhomZ
 physical length in Z direction (if homogeneous)
 
int m_npointsX
 number of points in X direction (if homogeneous)
 
int m_npointsY
 number of points in Y direction (if homogeneous)
 
int m_npointsZ
 number of points in Z direction (if homogeneous)
 
int m_HomoDirec
 number of homogenous directions
 
- Protected Attributes inherited from Nektar::SolverUtils::ALEHelper
Array< OneD, MultiRegions::ExpListSharedPtrm_fieldsALE
 
Array< OneD, Array< OneD, NekDouble > > m_gridVelocity
 
Array< OneD, Array< OneD, NekDouble > > m_gridVelocityTrace
 
std::vector< ALEBaseShPtrm_ALEs
 
bool m_ALESolver = false
 
bool m_meshDistorted = false
 
bool m_implicitALESolver = false
 
bool m_updateNormals = false
 
NekDouble m_prevStageTime = 0.0
 
int m_spaceDim
 

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
 

Friends

class MemoryManager< VCSMapping >
 

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

Constructor & Destructor Documentation

◆ VCSMapping()

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

Constructor. Creates ...

Parameters

param

Definition at line 58 of file VCSMapping.cpp.

60 : UnsteadySystem(pSession, pGraph),
61 VelocityCorrectionScheme(pSession, pGraph)
62{
63}
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)

◆ ~VCSMapping()

Nektar::VCSMapping::~VCSMapping ( )
overrideprotecteddefault

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

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

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

53 {
56 p->InitObject();
57 return p;
58 }
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.
std::vector< double > p(NPUPPER)

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

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

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

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

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

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

901{
902 size_t physTot = m_fields[0]->GetTotPoints();
903 size_t nvel = m_nConvectiveFields;
904
905 // Calculate g^(ij)p_(,j)
906 m_mapping->RaiseIndex(m_gradP, outarray);
907
908 // Calculate correction = (nabla p)/J - g^(ij)p_,j
909 // (Jac is not required if it is constant)
910 if (!m_mapping->HasConstantJacobian())
911 {
912 Array<OneD, NekDouble> Jac(physTot, 0.0);
913 m_mapping->GetJacobian(Jac);
914 for (size_t i = 0; i < nvel; ++i)
915 {
916 Vmath::Vdiv(physTot, m_gradP[i], 1, Jac, 1, m_gradP[i], 1);
917 }
918 }
919 for (size_t i = 0; i < nvel; ++i)
920 {
921 Vmath::Vsub(physTot, m_gradP[i], 1, outarray[i], 1, outarray[i], 1);
922 }
923}
Array< OneD, Array< OneD, NekDouble > > m_gradP
Definition VCSMapping.h:88
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.hpp:126
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.hpp:220

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

928{
929 // L(U) - 1.0*d^2(u^i)/dx^jdx^j
930 m_mapping->VelocityLaplacian(velPhys, outarray, 1.0);
931}

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::SolverUtils::UnsteadySystem.

Definition at line 126 of file VCSMapping.cpp.

127{
128 UnsteadySystem::v_DoInitialise(dumpInitialConditions);
129
130 // Set up Field Meta Data for output files
131 m_fieldMetaDataMap["Kinvis"] = boost::lexical_cast<std::string>(m_kinvis);
132 m_fieldMetaDataMap["TimeStep"] =
133 boost::lexical_cast<std::string>(m_timestep);
134
135 // Correct Dirichlet boundary conditions to account for mapping
136 m_mapping->UpdateBCs(0.0);
137 //
138 m_F = Array<OneD, Array<OneD, NekDouble>>(m_nConvectiveFields);
139 for (int i = 0; i < m_nConvectiveFields; ++i)
140 {
141 m_fields[i]->ImposeDirichletConditions(m_fields[i]->UpdateCoeffs());
142 m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(),
143 m_fields[i]->UpdatePhys());
144 m_F[i] = Array<OneD, NekDouble>(m_fields[0]->GetTotPoints(), 0.0);
145 }
146
147 // Initialise m_gradP
148 size_t physTot = m_fields[0]->GetTotPoints();
149 m_gradP = Array<OneD, Array<OneD, NekDouble>>(m_nConvectiveFields);
150 for (int i = 0; i < m_nConvectiveFields; ++i)
151 {
152 m_gradP[i] = Array<OneD, NekDouble>(physTot, 0.0);
154 m_pressure->GetPhys(), m_gradP[i]);
155 if (m_pressure->GetWaveSpace())
156 {
157 m_pressure->HomogeneousBwdTrans(physTot, m_gradP[i], m_gradP[i]);
158 }
159 }
160}
MultiRegions::ExpListSharedPtr m_pressure
Pointer to field holding pressure field.
NekDouble m_timestep
Time step size.
SOLVER_UTILS_EXPORT int GetTotPoints()
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
Map to identify relevant solver info to dump in output fields.
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 165 of file VCSMapping.cpp.

168{
169 EvaluateAdvectionTerms(inarray, outarray, time);
170
171 // Smooth advection
173 {
174 for (int i = 0; i < m_nConvectiveFields; ++i)
175 {
176 m_pressure->SmoothField(outarray[i]);
177 }
178 }
179
180 // Add forcing terms
181 for (auto &x : m_forcing)
182 {
183 x->Apply(m_fields, inarray, outarray, time);
184 }
185
186 // Add mapping terms
187 ApplyIncNSMappingForcing(inarray, outarray);
188
189 // Calculate High-Order pressure boundary conditions
190 m_extrapolation->EvaluatePressureBCs(inarray, outarray, m_kinvis);
191
192 // Update mapping and deal with Dirichlet boundary conditions
193 if (m_mapping->IsTimeDependent())
194 {
195 if (m_mapping->IsFromFunction())
196 {
197 // If the transformation is explicitly defined, update it here
198 // Otherwise, it will be done somewhere else (ForcingMovingBody)
199 m_mapping->UpdateMapping(time + m_timestep);
200 }
201 m_mapping->UpdateBCs(time + m_timestep);
202 }
203}
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)

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  DeclareFeld = true)
overrideprotectedvirtual

Initialisation object for EquationSystem.

Continuous field

Setting up the normals

Setting up the normals

Reimplemented from Nektar::IncNavierStokes.

Definition at line 65 of file VCSMapping.cpp.

66{
68
70 ASSERTL0(m_mapping, "Could not create mapping in VCSMapping.");
71
72 std::string vExtrapolation = "Mapping";
74 vExtrapolation, m_session, m_fields, m_pressure, m_velocity,
76 m_extrapolation->SubSteppingTimeIntegration(m_intScheme);
77 m_extrapolation->GenerateBndElmtExpansion();
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 {
86 m_intSteps = m_intScheme->GetOrder();
87 }
88 else
89 {
90 NEKERROR(ErrorUtil::efatal, "Integration method not suitable: "
91 "Options include IMEXOrder{1,2,3,4}");
92 }
93
94 m_presForcingCorrection = Array<OneD, Array<OneD, NekDouble>>(intSteps);
95 for (size_t i = 0; i < m_presForcingCorrection.size(); i++)
96 {
97 m_presForcingCorrection[i] = Array<OneD, NekDouble>(physTot, 0.0);
98 }
99 m_verbose = (m_session->DefinesCmdLineArgument("verbose")) ? true : false;
100
101 // Load solve parameters related to the mapping
102 // Flags determining if pressure/viscous terms should be treated implicitly
103 m_session->MatchSolverInfo("MappingImplicitPressure", "True",
104 m_implicitPressure, false);
105 m_session->MatchSolverInfo("MappingImplicitViscous", "True",
106 m_implicitViscous, false);
107 m_session->MatchSolverInfo("MappingNeglectViscous", "True",
108 m_neglectViscous, false);
109
111 {
112 m_implicitViscous = false;
113 }
114
115 // Tolerances and relaxation parameters for implicit terms
116 m_session->LoadParameter("MappingPressureTolerance", m_pressureTolerance,
117 1e-12);
118 m_session->LoadParameter("MappingViscousTolerance", m_viscousTolerance,
119 1e-12);
120 m_session->LoadParameter("MappingPressureRelaxation", m_pressureRelaxation,
121 1.0);
122 m_session->LoadParameter("MappingViscousRelaxation", m_viscousRelaxation,
123 1.0);
124}
#define ASSERTL0(condition, msg)
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
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:264
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.
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:82
NekDouble m_viscousRelaxation
Definition VCSMapping.h:85
NekDouble m_pressureRelaxation
Definition VCSMapping.h:84
NekDouble m_viscousTolerance
Definition VCSMapping.h:83
Array< OneD, Array< OneD, NekDouble > > m_presForcingCorrection
Definition VCSMapping.h:124
void v_InitObject(bool DeclareField=true) override
Initialisation object for EquationSystem.
ExtrapolateFactory & GetExtrapolateFactory()

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

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

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

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

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

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

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

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

Friends And Related Symbol Documentation

◆ MemoryManager< VCSMapping >

friend class MemoryManager< VCSMapping >
friend

Definition at line 156 of file VCSMapping.h.

Member Data Documentation

◆ className

std::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.
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
Definition VCSMapping.h:50
EquationSystemFactory & GetEquationSystemFactory()

Name of class.

Definition at line 61 of file VCSMapping.h.

◆ m_gradP

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

Definition at line 88 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 79 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 84 of file VCSMapping.h.

Referenced by v_InitObject(), and v_SolvePressure().

◆ m_pressureTolerance

NekDouble Nektar::VCSMapping::m_pressureTolerance
protected

Definition at line 82 of file VCSMapping.h.

Referenced by v_InitObject(), and v_SolvePressure().

◆ m_verbose

bool Nektar::VCSMapping::m_verbose
protected

Definition at line 73 of file VCSMapping.h.

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

◆ m_viscousRelaxation

NekDouble Nektar::VCSMapping::m_viscousRelaxation
protected

Definition at line 85 of file VCSMapping.h.

Referenced by v_InitObject(), and v_SolveViscous().

◆ m_viscousTolerance

NekDouble Nektar::VCSMapping::m_viscousTolerance
protected

Definition at line 83 of file VCSMapping.h.

Referenced by v_InitObject(), and v_SolveViscous().

◆ solverTypeLookupId

std::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 71 of file VCSMapping.h.