Nektar++
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. 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
SOLVER_UTILS_EXPORT ~UnsteadySystem () override=default
 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 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. More...
 
SOLVER_UTILS_EXPORT LibUtilities::TimeIntegrationSchemeOperatorsGetTimeIntegrationSchemeOperators ()
 Returns the time integration scheme operators. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT void InitObject (bool DeclareField=true)
 Initialises the members of this object. More...
 
SOLVER_UTILS_EXPORT void DoInitialise (bool dumpInitialConditions=true)
 Perform any initialisation necessary before solving the problem. More...
 
SOLVER_UTILS_EXPORT void DoSolve ()
 Solve the problem. More...
 
SOLVER_UTILS_EXPORT void TransCoeffToPhys ()
 Transform from coefficient to physical space. More...
 
SOLVER_UTILS_EXPORT void TransPhysToCoeff ()
 Transform from physical to coefficient space. More...
 
SOLVER_UTILS_EXPORT void Output ()
 Perform output operations after solve. More...
 
SOLVER_UTILS_EXPORT std::string GetSessionName ()
 Get Session name. More...
 
template<class T >
std::shared_ptr< T > as ()
 
SOLVER_UTILS_EXPORT void ResetSessionName (std::string newname)
 Reset Session name. More...
 
SOLVER_UTILS_EXPORT LibUtilities::SessionReaderSharedPtr GetSession ()
 Get Session name. More...
 
SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr GetPressure ()
 Get pressure field if available. More...
 
SOLVER_UTILS_EXPORT void ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 
SOLVER_UTILS_EXPORT void PrintSummary (std::ostream &out)
 Print a summary of parameters and solver characteristics. More...
 
SOLVER_UTILS_EXPORT void SetLambda (NekDouble lambda)
 Set parameter m_lambda. More...
 
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction (std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
 Get a SessionFunction by name. More...
 
SOLVER_UTILS_EXPORT void SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 Initialise the data in the dependent fields. More...
 
SOLVER_UTILS_EXPORT void EvaluateExactSolution (int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 Evaluates an exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln, bool Normalised=false)
 Compute the L2 error between fields and a given exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, bool Normalised=false)
 Compute the L2 error of the fields. More...
 
SOLVER_UTILS_EXPORT NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation. More...
 
SOLVER_UTILS_EXPORT NekDouble H1Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln, bool Normalised=false)
 Compute the H1 error between fields and a given exact solution. 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 GetTime ()
 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 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. More...
 
SOLVER_UTILS_EXPORT bool NegatedOp ()
 Identify if operator is negated in DoSolve. More...
 
- 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. More...
 
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)
 
const Array< OneD, const Array< OneD, NekDouble > > & GetGridVelocity ()
 
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)
 
- 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, 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. More...
 
SOLVER_UTILS_EXPORT void GetAeroForce (Array< OneD, NekDouble > forces)
 Get aerodynamic force and moment. More...
 

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

 VCSMapping (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
 ~VCSMapping () override=default
 
void v_InitObject (bool DeclareField=true) override
 Initialisation object for EquationSystem. More...
 
void v_DoInitialise (bool dumpInitialConditions=true) override
 Sets up initial conditions. More...
 
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. More...
 
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...
 
bool v_PostIntegrate (int step) override
 
void v_GenerateSummary (SolverUtils::SummaryList &s) override
 Print a summary of time stepping parameters. More...
 
void v_TransCoeffToPhys (void) override
 Virtual function for transformation to physical space. More...
 
void v_TransPhysToCoeff (void) override
 Virtual function for transformation to coefficient space. More...
 
void v_DoInitialise (bool dumpInitialConditions=true) override
 Sets up initial conditions. More...
 
Array< OneD, bool > v_GetSystemSingularChecks () override
 
int v_GetForceDimension () override
 
virtual void v_SetUpPressureForcing (const Array< OneD, const Array< OneD, NekDouble > > &fields, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt)
 
virtual void v_SetUpViscousForcing (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt)
 
virtual void v_SolvePressure (const Array< OneD, NekDouble > &Forcing)
 
virtual void v_SolveViscous (const Array< OneD, const Array< OneD, NekDouble > > &Forcing, const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble aii_Dt)
 
virtual void v_SolveUnsteadyStokesSystem (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time, const NekDouble a_iixDt)
 
virtual void v_EvaluateAdvection_SetPressureBCs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 
bool v_RequireFwdTrans () override
 
virtual std::string v_GetExtrapolateStr (void)
 
virtual std::string v_GetSubSteppingExtrapolateStr (const std::string &instr)
 
void SetUpSVV (void)
 
void SetUpExtrapolation (void)
 
void SVVVarDiffCoeff (const NekDouble velmag, Array< OneD, NekDouble > &diffcoeff, const Array< OneD, Array< OneD, NekDouble > > &vel=NullNekDoubleArrayOfArray)
 
void AppendSVVFactors (StdRegions::ConstFactorMap &factors, MultiRegions::VarFactorsMap &varFactorsMap)
 
void ComputeGJPNormalVelocity (const Array< OneD, const Array< OneD, NekDouble > > &inarray, StdRegions::VarCoeffMap &varcoeffs)
 
- Protected Member Functions inherited from Nektar::IncNavierStokes
 IncNavierStokes (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Constructor. More...
 
 ~IncNavierStokes () override=default
 
void v_InitObject (bool DeclareField=true) override
 Initialisation object for EquationSystem. More...
 
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 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...
 
MultiRegions::ExpListSharedPtr v_GetPressure () override
 
void v_TransCoeffToPhys (void) override
 Virtual function for transformation to physical space. More...
 
void v_TransPhysToCoeff (void) override
 Virtual function for transformation to coefficient space. More...
 
virtual int v_GetForceDimension ()=0
 
Array< OneD, NekDoublev_GetMaxStdVelocity (const NekDouble SpeedSoundFactor) override
 
bool v_PreIntegrate (int step) override
 
SOLVER_UTILS_EXPORT void v_InitObject (bool DeclareField=true) override
 Initialisation object for EquationSystem. More...
 
SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step) override
 
virtual SOLVER_UTILS_EXPORT Array< OneD, NekDoublev_GetMaxStdVelocity (const NekDouble SpeedSoundFactor=1.0)
 
- Protected Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
SOLVER_UTILS_EXPORT UnsteadySystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises UnsteadySystem class members. More...
 
SOLVER_UTILS_EXPORT void v_InitObject (bool DeclareField=true) override
 Init object for UnsteadySystem class. More...
 
SOLVER_UTILS_EXPORT void v_DoSolve () override
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_PrintStatusInformation (const int step, const NekDouble cpuTime)
 Print Status Information. More...
 
virtual SOLVER_UTILS_EXPORT void v_PrintSummaryStatistics (const NekDouble intTime)
 Print Summary Statistics. More...
 
SOLVER_UTILS_EXPORT void v_DoInitialise (bool dumpInitialConditions=true) override
 Sets up initial conditions. More...
 
SOLVER_UTILS_EXPORT void v_GenerateSummary (SummaryList &s) override
 Print a summary of time stepping parameters. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_GetTimeStep (const Array< OneD, const Array< OneD, NekDouble > > &inarray)
 Return the timestep to be used for the next step in the time-marching loop. More...
 
virtual SOLVER_UTILS_EXPORT bool v_PreIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT bool v_PostIntegrate (int step)
 
virtual SOLVER_UTILS_EXPORT bool v_RequireFwdTrans ()
 
virtual SOLVER_UTILS_EXPORT void v_SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
virtual SOLVER_UTILS_EXPORT bool v_UpdateTimeStepCheck ()
 
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control. More...
 
SOLVER_UTILS_EXPORT void CheckForRestartTime (NekDouble &time, int &nchk)
 
SOLVER_UTILS_EXPORT void SVVVarDiffCoeff (const Array< OneD, Array< OneD, NekDouble > > vel, StdRegions::VarCoeffMap &varCoeffMap)
 Evaluate the SVV diffusion coefficient according to Moura's paper where it should proportional to h time velocity. More...
 
SOLVER_UTILS_EXPORT void DoDummyProjection (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 Perform dummy projection. More...
 
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises EquationSystem class members. More...
 
virtual SOLVER_UTILS_EXPORT void v_InitObject (bool DeclareFeld=true)
 Initialisation object for EquationSystem. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoInitialise (bool dumpInitialConditions=true)
 Virtual function for initialisation implementation. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve ()
 Virtual function for solve implementation. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Virtual function for the L_inf error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false)
 Virtual function for the L_2 error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT 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. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransCoeffToPhys ()
 Virtual function for transformation to physical space. More...
 
virtual SOLVER_UTILS_EXPORT void v_TransPhysToCoeff ()
 Virtual function for transformation to coefficient space. More...
 
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary (SummaryList &l)
 Virtual function for generating summary information. More...
 
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 
virtual SOLVER_UTILS_EXPORT void v_EvaluateExactSolution (unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 
virtual SOLVER_UTILS_EXPORT void v_Output (void)
 
virtual SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr v_GetPressure (void)
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp (void)
 Virtual function to identify if operator is negated in DoSolve. More...
 
virtual SOLVER_UTILS_EXPORT void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 
virtual SOLVER_UTILS_EXPORT void v_GetVelocity (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &velocity)=0
 
virtual SOLVER_UTILS_EXPORT bool v_HasConstantDensity ()=0
 
virtual SOLVER_UTILS_EXPORT void v_GetDensity (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &density)=0
 
virtual SOLVER_UTILS_EXPORT void v_GetPressure (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &pressure)=0
 
virtual SOLVER_UTILS_EXPORT void v_SetMovingFrameVelocities (const Array< OneD, NekDouble > &vFrameVels, const int step)
 
virtual SOLVER_UTILS_EXPORT bool v_GetMovingFrameVelocities (Array< OneD, NekDouble > &vFrameVels, const int step)
 
virtual SOLVER_UTILS_EXPORT void v_SetMovingFrameDisp (const Array< OneD, NekDouble > &vFrameDisp, const int step)
 
virtual SOLVER_UTILS_EXPORT void v_SetMovingFramePivot (const Array< OneD, NekDouble > &vFramePivot)
 
virtual SOLVER_UTILS_EXPORT bool v_GetMovingFrameDisp (Array< OneD, NekDouble > &vFrameDisp, const int step)
 
virtual SOLVER_UTILS_EXPORT void v_SetAeroForce (Array< OneD, NekDouble > forces)
 
virtual SOLVER_UTILS_EXPORT void v_GetAeroForce (Array< OneD, NekDouble > forces)
 

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
 
IncBoundaryConditionsSharedPtr m_IncNavierStokesBCs
 
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
 pivot point for moving reference frame More...
 
Array< OneD, NekDoublem_aeroForces
 
std::map< int, std::map< int, WomersleyParamsSharedPtr > > m_womersleyParams
 Womersley parameters if required. More...
 
- Protected Attributes inherited from Nektar::SolverUtils::AdvectionSystem
SolverUtils::AdvectionSharedPtr m_advObject
 Advection term. More...
 
- Protected Attributes inherited from Nektar::SolverUtils::UnsteadySystem
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
Array< OneD, Array< OneD, NekDouble > > m_previousSolution
 Storage for previous solution for steady-state check. More...
 
std::vector< int > m_intVariables
 
NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1). More...
 
NekDouble m_CFLGrowth
 CFL growth rate. More...
 
NekDouble m_CFLEnd
 Maximun cfl in cfl growth. More...
 
int m_abortSteps
 Number of steps between checks for abort conditions. More...
 
bool m_explicitDiffusion
 Indicates if explicit or implicit treatment of diffusion is used. More...
 
bool m_explicitAdvection
 Indicates if explicit or implicit treatment of advection is used. More...
 
bool m_explicitReaction
 Indicates if explicit or implicit treatment of reaction is used. More...
 
int m_steadyStateSteps
 Check for steady state at step interval. More...
 
NekDouble m_steadyStateTol
 Tolerance to which steady state should be evaluated at. More...
 
int m_filtersInfosteps
 Number of time steps between outputting filters information. More...
 
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
 
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state. More...
 
std::ofstream m_errFile
 
NekDouble m_epsilon
 Diffusion coefficient. More...
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
bool m_verbose
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader. More...
 
std::map< std::string, SolverUtils::SessionFunctionSharedPtrm_sessionFunctions
 Map of known SessionFunctions. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 Array holding all dependent variables. More...
 
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
 Pointer to boundary conditions object. More...
 
SpatialDomains::MeshGraphSharedPtr m_graph
 Pointer to graph defining mesh. More...
 
std::string m_sessionName
 Name of the session. More...
 
NekDouble m_time
 Current time of simulation. More...
 
int m_initialStep
 Number of the step where the simulation should begin. More...
 
NekDouble m_fintime
 Finish time of the simulation. More...
 
NekDouble m_timestep
 Time step size. More...
 
NekDouble m_lambda
 Lambda constant in real system if one required. More...
 
NekDouble m_checktime
 Time between checkpoints. More...
 
NekDouble m_lastCheckTime
 
NekDouble m_TimeIncrementFactor
 
int m_nchk
 Number of checkpoints written so far. More...
 
int m_steps
 Number of steps to take. More...
 
int m_checksteps
 Number of steps between checkpoints. More...
 
int m_infosteps
 Number of time steps between outputting status information. More...
 
int m_iterPIT = 0
 Number of parallel-in-time time iteration. More...
 
int m_windowPIT = 0
 Index of windows for parallel-in-time time iteration. More...
 
int m_spacedim
 Spatial dimension (>= expansion dim). More...
 
int m_expdim
 Expansion dimension. More...
 
bool m_singleMode
 Flag to determine if single homogeneous mode is used. More...
 
bool m_halfMode
 Flag to determine if half homogeneous mode is used. More...
 
bool m_multipleModes
 Flag to determine if use multiple homogenenous modes are used. More...
 
bool m_useFFT
 Flag to determine if FFT is used for homogeneous transform. More...
 
bool m_homogen_dealiasing
 Flag to determine if dealiasing is used for homogeneous simulations. More...
 
bool m_specHP_dealiasing
 Flag to determine if dealisising is usde for the Spectral/hp element discretisation. More...
 
enum MultiRegions::ProjectionType m_projectionType
 Type of projection; e.g continuous or discontinuous. More...
 
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
 Array holding trace normals for DG simulations in the forwards direction. More...
 
Array< OneD, bool > m_checkIfSystemSingular
 Flag to indicate if the fields should be checked for singularity. More...
 
LibUtilities::FieldMetaDataMap m_fieldMetaDataMap
 Map to identify relevant solver info to dump in output fields. More...
 
Array< OneD, NekDoublem_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. More...
 
std::vector< std::string > m_strFrameData
 variable name in m_movingFrameData 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...
 
- 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_ImplicitALESolver = 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 714 of file VCSMapping.cpp.

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

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

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

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

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

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

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

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

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

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

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 m_intScheme->GetName() == "IMEXGear")
86 {
87 m_intSteps = m_intScheme->GetOrder();
88 }
89 else
90 {
92 "Integration method not suitable: "
93 "Options include IMEXGear or IMEXOrder{1,2,3,4}");
94 }
95
96 m_presForcingCorrection = Array<OneD, Array<OneD, NekDouble>>(intSteps);
97 for (size_t i = 0; i < m_presForcingCorrection.size(); i++)
98 {
99 m_presForcingCorrection[i] = Array<OneD, NekDouble>(physTot, 0.0);
100 }
101 m_verbose = (m_session->DefinesCmdLineArgument("verbose")) ? true : false;
102
103 // Load solve parameters related to the mapping
104 // Flags determining if pressure/viscous terms should be treated implicitly
105 m_session->MatchSolverInfo("MappingImplicitPressure", "True",
106 m_implicitPressure, false);
107 m_session->MatchSolverInfo("MappingImplicitViscous", "True",
108 m_implicitViscous, false);
109 m_session->MatchSolverInfo("MappingNeglectViscous", "True",
110 m_neglectViscous, false);
111
113 {
114 m_implicitViscous = false;
115 }
116
117 // Tolerances and relaxation parameters for implicit terms
118 m_session->LoadParameter("MappingPressureTolerance", m_pressureTolerance,
119 1e-12);
120 m_session->LoadParameter("MappingViscousTolerance", m_viscousTolerance,
121 1e-12);
122 m_session->LoadParameter("MappingPressureRelaxation", m_pressureRelaxation,
123 1.0);
124 m_session->LoadParameter("MappingViscousRelaxation", m_viscousRelaxation,
125 1.0);
126}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
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()
Definition: Extrapolate.cpp:47

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

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

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

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

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

◆ v_SolveViscous()

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

Solve velocity system

Reimplemented from Nektar::VelocityCorrectionScheme.

Definition at line 550 of file VCSMapping.cpp.

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

References ASSERTL0, Nektar::StdRegions::eFactorLambda, Nektar::StdRegions::eFactorSVVCutoffRatio, Nektar::StdRegions::eFactorSVVDiffCoeff, Nektar::VarcoeffHashingTest::factors, Nektar::SolverUtils::EquationSystem::m_fields, m_implicitViscous, Nektar::IncNavierStokes::m_kinvis, m_mapping, Nektar::IncNavierStokes::m_nConvectiveFields, Nektar::SolverUtils::EquationSystem::m_session, Nektar::VelocityCorrectionScheme::m_sVVCutoffRatio, Nektar::VelocityCorrectionScheme::m_sVVDiffCoeff, Nektar::VelocityCorrectionScheme::m_useSpecVanVisc, m_verbose, m_viscousRelaxation, m_viscousTolerance, Vmath::Smul(), Nektar::VelocityCorrectionScheme::v_SolveViscous(), Vmath::Vadd(), and Vmath::Vcopy().

Friends And Related Function 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.