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

#include <VCSMapping.h>

Inheritance diagram for Nektar::VCSMapping:
[legend]

Public Member Functions

 VCSMapping (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Constructor. More...
 
void ApplyIncNSMappingForcing (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
 ~VCSMapping () override
 
void v_InitObject (bool DeclareField=true) override
 Initialisation object for EquationSystem. More...
 
- Public Member Functions inherited from Nektar::VelocityCorrectionScheme
 VelocityCorrectionScheme (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Constructor. More...
 
 ~VelocityCorrectionScheme () override
 
void v_InitObject (bool DeclareField=true) override
 Initialisation object for EquationSystem. More...
 
void SetUpPressureForcing (const Array< OneD, const Array< OneD, NekDouble > > &fields, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt)
 
void SetUpViscousForcing (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt)
 
void SolvePressure (const Array< OneD, NekDouble > &Forcing)
 
void SolveViscous (const Array< OneD, const Array< OneD, NekDouble > > &Forcing, const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble aii_Dt)
 
void SolveUnsteadyStokesSystem (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time, const NekDouble a_iixDt)
 
void EvaluateAdvection_SetPressureBCs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 
- Public Member Functions inherited from Nektar::IncNavierStokes
 ~IncNavierStokes () override
 
void v_InitObject (bool DeclareField=true) override
 Initialisation object for EquationSystem. More...
 
int GetNConvectiveFields (void)
 
void AddForcing (const SolverUtils::ForcingSharedPtr &pForce)
 
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
 
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
 
SOLVER_UTILS_EXPORT void v_InitObject (bool DeclareField=true) override
 Initialisation object for EquationSystem. More...
 
SOLVER_UTILS_EXPORT AdvectionSharedPtr GetAdvObject ()
 Returns the advection object held by this instance. More...
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleGetElmtCFLVals (const bool FlagAcousticCFL=true)
 
SOLVER_UTILS_EXPORT NekDouble GetCFLEstimate (int &elmtid)
 
- Public Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
SOLVER_UTILS_EXPORT ~UnsteadySystem () override
 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 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 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 SetSteps (const int steps)
 
SOLVER_UTILS_EXPORT void ZeroPhysFields ()
 
SOLVER_UTILS_EXPORT void FwdTransFields ()
 
SOLVER_UTILS_EXPORT void SetModifiedBasis (const bool modbasis)
 
SOLVER_UTILS_EXPORT int GetCheckpointNumber ()
 
SOLVER_UTILS_EXPORT void SetCheckpointNumber (int num)
 
SOLVER_UTILS_EXPORT int GetCheckpointSteps ()
 
SOLVER_UTILS_EXPORT void SetCheckpointSteps (int num)
 
SOLVER_UTILS_EXPORT int GetInfoSteps ()
 
SOLVER_UTILS_EXPORT void SetInfoSteps (int num)
 
SOLVER_UTILS_EXPORT void SetIterationNumberPIT (int num)
 
SOLVER_UTILS_EXPORT void SetWindowNumberPIT (int num)
 
SOLVER_UTILS_EXPORT Array< OneD, const Array< OneD, NekDouble > > GetTraceNormals ()
 
SOLVER_UTILS_EXPORT void SetTime (const NekDouble time)
 
SOLVER_UTILS_EXPORT void SetTimeStep (const NekDouble timestep)
 
SOLVER_UTILS_EXPORT void SetInitialStep (const int step)
 
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time. More...
 
SOLVER_UTILS_EXPORT bool NegatedOp ()
 Identify if operator is negated in DoSolve. More...
 
- 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

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
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...
 
EquationType GetEquationType (void)
 
void EvaluateAdvectionTerms (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 
void WriteModalEnergy (void)
 
void SetBoundaryConditions (NekDouble time)
 time dependent boundary conditions updating More...
 
void SetRadiationBoundaryForcing (int fieldid)
 Set Radiation forcing term. More...
 
void SetZeroNormalVelocity ()
 Set Normal Velocity Component to Zero. More...
 
void SetWomersleyBoundary (const int fldid, const int bndid)
 Set Womersley Profile if specified. More...
 
void SetUpWomersley (const int fldid, const int bndid, std::string womstr)
 Set Up Womersley details. More...
 
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 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 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
 

Additional Inherited Members

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

Detailed Description

Definition at line 43 of file VCSMapping.h.

Constructor & Destructor Documentation

◆ VCSMapping()

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

Constructor.

Constructor. Creates ...

Parameters

param

Definition at line 59 of file VCSMapping.cpp.

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

◆ ~VCSMapping()

Nektar::VCSMapping::~VCSMapping ( void  )
override

Destructor

Definition at line 132 of file VCSMapping.cpp.

133{
134}

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

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

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

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

◆ MappingAccelerationCorrection()

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

Definition at line 830 of file VCSMapping.cpp.

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

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

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

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

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

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

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

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

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)
overridevirtual

Initialisation object for EquationSystem.

Continuous field

Setting up the normals

Setting up the normals

Reimplemented from Nektar::IncNavierStokes.

Definition at line 66 of file VCSMapping.cpp.

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

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

◆ v_SetUpPressureForcing()

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

Forcing term for Poisson solver solver

Reimplemented from Nektar::VelocityCorrectionScheme.

Definition at line 218 of file VCSMapping.cpp.

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

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

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

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

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

Member Data Documentation

◆ className

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

Name of class.

Definition at line 58 of file VCSMapping.h.

◆ m_gradP

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

Definition at line 94 of file VCSMapping.h.

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

◆ m_implicitPressure

bool Nektar::VCSMapping::m_implicitPressure
protected

◆ m_implicitViscous

bool Nektar::VCSMapping::m_implicitViscous
protected

◆ m_mapping

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

◆ m_neglectViscous

bool Nektar::VCSMapping::m_neglectViscous
protected

Definition at line 85 of file VCSMapping.h.

Referenced by ApplyIncNSMappingForcing(), and v_InitObject().

◆ m_presForcingCorrection

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

Definition at line 123 of file VCSMapping.h.

Referenced by v_InitObject(), and v_SetUpPressureForcing().

◆ m_pressureRelaxation

NekDouble Nektar::VCSMapping::m_pressureRelaxation
protected

Definition at line 90 of file VCSMapping.h.

Referenced by v_InitObject(), and v_SolvePressure().

◆ m_pressureTolerance

NekDouble Nektar::VCSMapping::m_pressureTolerance
protected

Definition at line 88 of file VCSMapping.h.

Referenced by v_InitObject(), and v_SolvePressure().

◆ m_verbose

bool Nektar::VCSMapping::m_verbose
protected

Definition at line 79 of file VCSMapping.h.

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

◆ m_viscousRelaxation

NekDouble Nektar::VCSMapping::m_viscousRelaxation
protected

Definition at line 91 of file VCSMapping.h.

Referenced by v_InitObject(), and v_SolveViscous().

◆ m_viscousTolerance

NekDouble Nektar::VCSMapping::m_viscousTolerance
protected

Definition at line 89 of file VCSMapping.h.

Referenced by v_InitObject(), and v_SolveViscous().

◆ solverTypeLookupId

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

Definition at line 77 of file VCSMapping.h.