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

#include <SmoothedProfileMethod.h>

Inheritance diagram for Nektar::SmoothedProfileMethod:
[legend]

Public Member Functions

 SmoothedProfileMethod (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Construct a new Smoothed Profile Method object. More...
 
 ~SmoothedProfileMethod () override
 Destroy the Smoothed Profile Method object. More...
 
void v_InitObject (bool DeclareField=true) override
 Initialisation object for EquationSystem. More...
 
void v_GenerateSummary (SolverUtils::SummaryList &s) override
 Generates the summary of the current simulation. 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_SolveUnsteadyStokesSystem (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, NekDouble time, NekDouble a_iixDt) override
 Linear terms due to pressure and visosity are calculated here. After solving the velocity filed without taking into account the immersed boundaries, a new correction is applied through the force \(f_s\): More...
 
void SetUpCorrectionPressure (const Array< OneD, const Array< OneD, NekDouble > > &fields, Array< OneD, Array< OneD, NekDouble > > &Forcing)
 Sets the forcing term of the equation for the correction pressure \(p_p\): More...
 
void SolveCorrectionPressure (const Array< OneD, NekDouble > &Forcing)
 Solves the Poisson equation for the correction pressure \(p_p\): More...
 
void SolveCorrectedVelocity (Array< OneD, Array< OneD, NekDouble > > &Forcing, Array< OneD, Array< OneD, NekDouble > > &fields, NekDouble dt)
 Corrects the velocity field so that the IBs are taken into account. Solves the explicit equation: More...
 
void SetCorrectionPressureBCs ()
 Updates the BCs for boundaries with Neumann or periodic BCs in the pressure: More...
 
void UpdatePhiUp (NekDouble time)
 Calculates the values of the shape function. More...
 
void UpdateForcing (const Array< OneD, const Array< OneD, NekDouble > > &fields, NekDouble dt)
 For a body with a velocity \(\mathbf{u_p}\), the force \(\mathbf{f_s}\) applied to the fluid ensures that the IBC are met: More...
 
bool GetVarTimeDependence (std::string funcName, std::string attrName)
 True if the function is timedependent, false otherwise. More...
 
TiXmlElement * GetFunctionHdl (std::string functionName)
 Returns a handle to the requested function. Returns NULL if it does not exist. More...
 
void ReadPhi ()
 
template<typename T >
void SetUpExpansions (int nvel)
 Initialises the expansions for the intermediate pressure, the 'Phi' function and the IB forcing. More...
 
- 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

MultiRegions::ExpListSharedPtr m_pressureP
 Correction pressure field for SPM. More...
 
Array< OneD, Array< OneD, NekDouble > > m_up
 Velocity of the immersed body(ies) More...
 
Array< OneD, Array< OneD, NekDouble > > m_upPrev
 
std::vector< std::string > m_velName
 Vector storing the names of the components of \u_p. More...
 
bool m_timeDependentUp
 Flag signaling if \(\u_p\) depends on time. More...
 
NekDouble m_gamma0
 Stiffly-stable scheme \(\gamma_0\) coefficient. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fs
 Forcing function 'f_s'. More...
 
MultiRegions::ExpListSharedPtr m_phi
 Shape function 'phi' as expansion list. More...
 
SolverUtils::SessionFunctionSharedPtr m_phiEvaluator
 Function that evaluates the values of \Phi. More...
 
bool m_timeDependentPhi
 Flag that is true when phi depends on time. More...
 
bool m_filePhi
 Flag indicating that phi was defined in a file. More...
 
int m_forcesFilter
 Position of "AeroForcesSPM" filter in 'm_session->GetFilters()'. More...
 
- 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 []
 

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 SmoothedProfileMethod.h.

Constructor & Destructor Documentation

◆ SmoothedProfileMethod()

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

Construct a new Smoothed Profile Method object.

Parameters
pSession
pGraph

Definition at line 62 of file SmoothedProfileMethod.cpp.

65 : UnsteadySystem(pSession, pGraph),
66 VelocityCorrectionScheme(pSession, pGraph)
67{
68}
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.

◆ ~SmoothedProfileMethod()

Nektar::SmoothedProfileMethod::~SmoothedProfileMethod ( void  )
override

Destroy the Smoothed Profile Method object.

Definition at line 74 of file SmoothedProfileMethod.cpp.

75{
76}

Member Function Documentation

◆ create()

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

Creates an instance of this class.

Definition at line 48 of file SmoothedProfileMethod.h.

51 {
54 pGraph);
55 p->InitObject();
56 return p;
57 }
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.

◆ GetFunctionHdl()

TiXmlElement * Nektar::SmoothedProfileMethod::GetFunctionHdl ( std::string  functionName)
protected

Returns a handle to the requested function. Returns NULL if it does not exist.

Parameters
functionName
Returns
TiXmlElement*

Definition at line 572 of file SmoothedProfileMethod.cpp.

573{
574 // Get the handler of first function block
575 TiXmlElement *conds = m_session->GetElement("Nektar/Conditions");
576 TiXmlElement *function = conds->FirstChildElement("FUNCTION");
577
578 // Loop over functions until the block 'name' is found
579 string functionType = function->Attribute("NAME");
580 while (function && !boost::iequals(functionType, functionName))
581 {
582 function = function->NextSiblingElement("FUNCTION");
583 functionType = function->Attribute("NAME");
584 }
585
586 return function;
587}
LibUtilities::SessionReaderSharedPtr m_session
The session reader.

References Nektar::SolverUtils::EquationSystem::m_session.

Referenced by GetVarTimeDependence(), and ReadPhi().

◆ GetVarTimeDependence()

bool Nektar::SmoothedProfileMethod::GetVarTimeDependence ( std::string  funcName,
std::string  attrName 
)
protected

True if the function is timedependent, false otherwise.

Parameters
name
type
attribute
Returns
string

Definition at line 532 of file SmoothedProfileMethod.cpp.

534{
535 // Get the handler of the function
536 TiXmlElement *function = GetFunctionHdl(funcName);
537
538 // Go to the first element
539 TiXmlElement *functionDef = function->FirstChildElement();
540 ASSERTL0(functionDef, "At least one element must be defined in " + funcName)
541
542 // And search the element with name 'elemName'
543 string varName = functionDef->Attribute("VAR");
544 while (functionDef && !boost::iequals(varName, elemName))
545 {
546 functionDef = functionDef->NextSiblingElement();
547 varName = functionDef->Attribute("VAR");
548 }
549
550 ASSERTL0(functionDef,
551 "Variable " + elemName + " must be defined in " + funcName + ".");
552
553 // And return the value of USERDEFINEDTYPE
554 string attr;
555 int err = functionDef->QueryStringAttribute("USERDEFINEDTYPE", &attr);
556 bool output = boost::iequals(attr, "TimeDependent");
557
558 ASSERTL0((err == TIXML_NO_ATTRIBUTE) || (err == TIXML_SUCCESS && output),
559 "USERDEFINEDTYPE in " + elemName +
560 " must be TimeDependent if defined");
561
562 return output;
563}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
TiXmlElement * GetFunctionHdl(std::string functionName)
Returns a handle to the requested function. Returns NULL if it does not exist.

References ASSERTL0, and GetFunctionHdl().

Referenced by ReadPhi().

◆ ReadPhi()

void Nektar::SmoothedProfileMethod::ReadPhi ( )
protected

Definition at line 589 of file SmoothedProfileMethod.cpp.

590{
591 // Function evaluator for Phi and Up
592 m_phiEvaluator = GetFunction("ShapeFunction");
593
594 TiXmlElement *function = GetFunctionHdl("ShapeFunction");
595 TiXmlElement *child = function->FirstChildElement();
596 m_filePhi = false;
597
598 // If defined by using a file
599 if (boost::iequals(child->ValueStr(), "F"))
600 {
601 // Get name of STL file
602 string fileName;
603 int status = child->QueryStringAttribute("FILE", &fileName);
604 ASSERTL0(status == TIXML_SUCCESS,
605 "An FLD file with the values "
606 "of the phi function has to be supplied.")
607 ASSERTL0(boost::iequals(fileName.substr(fileName.length() - 4), ".fld"),
608 "A valid FLD file must be supplied in the "
609 "'ShapeFunction' field.")
610
611 // Get phi values from XML file (after "FieldConvert" the STL file)
612 // First, load the data
613 std::vector<LibUtilities::FieldDefinitionsSharedPtr> fieldDef;
614 std::vector<std::vector<NekDouble>> fieldData;
615 LibUtilities::FieldMetaDataMap fieldMetaData;
616 LibUtilities::FieldIOSharedPtr phiFile =
617 LibUtilities::FieldIO::CreateForFile(m_session, fileName);
618 phiFile->Import(fileName, fieldDef, fieldData, fieldMetaData);
619
620 // Only Phi field should be defined in the file
621 ASSERTL0(fieldData.size() == 1, "Only one field (phi) must be "
622 "defined in the FLD file.")
623
624 // Extract Phi field to output
625 string tmp("phi");
626 m_phi->ExtractDataToCoeffs(fieldDef[0], fieldData[0], tmp,
627 m_phi->UpdateCoeffs());
628 m_phi->BwdTrans(m_phi->GetCoeffs(), m_phi->UpdatePhys());
629 m_filePhi = true;
630 m_timeDependentPhi = false;
631 m_timeDependentUp = false;
632 }
633 else
634 {
635 // Check if Phi is timedependent
636 m_timeDependentPhi = GetVarTimeDependence("ShapeFunction", "Phi");
637
638 // If so, check if its velocity changes as well
639 m_timeDependentUp = GetVarTimeDependence("ShapeFunction", "Up");
640 switch (m_velocity.size())
641 {
642 case 2:
644 GetVarTimeDependence("ShapeFunction", "Vp");
645 break;
646 case 3:
648 GetVarTimeDependence("ShapeFunction", "Vp");
650 GetVarTimeDependence("ShapeFunction", "Wp");
651 break;
652 }
653 }
654}
Array< OneD, int > m_velocity
int which identifies which components of m_fields contains the velocity (u,v,w);
bool m_filePhi
Flag indicating that phi was defined in a file.
bool m_timeDependentPhi
Flag that is true when phi depends on time.
SolverUtils::SessionFunctionSharedPtr m_phiEvaluator
Function that evaluates the values of \Phi.
bool m_timeDependentUp
Flag signaling if depends on time.
bool GetVarTimeDependence(std::string funcName, std::string attrName)
True if the function is timedependent, false otherwise.
MultiRegions::ExpListSharedPtr m_phi
Shape function 'phi' as expansion list.
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
Get a SessionFunction by name.
void Import(const std::string &infilename, std::vector< FieldDefinitionsSharedPtr > &fielddefs, std::vector< std::vector< NekDouble > > &fielddata, FieldMetaDataMap &fieldinfomap, const Array< OneD, int > &ElementIDs)
This function allows for data to be imported from an FLD file when a session and/or communicator is n...
Definition: FieldIO.cpp:288
std::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:322
std::map< std::string, std::string > FieldMetaDataMap
Definition: FieldIO.h:50
std::shared_ptr< FieldDefinitions > FieldDefinitionsSharedPtr
Definition: FieldIO.h:184
double NekDouble
STL namespace.

References ASSERTL0, Nektar::LibUtilities::FieldIO::CreateForFile(), Nektar::SolverUtils::EquationSystem::GetFunction(), GetFunctionHdl(), GetVarTimeDependence(), m_filePhi, m_phi, m_phiEvaluator, Nektar::SolverUtils::EquationSystem::m_session, m_timeDependentPhi, m_timeDependentUp, and Nektar::IncNavierStokes::m_velocity.

Referenced by v_InitObject().

◆ SetCorrectionPressureBCs()

void Nektar::SmoothedProfileMethod::SetCorrectionPressureBCs ( )
protected

Updates the BCs for boundaries with Neumann or periodic BCs in the pressure:

\[ \frac{\partial p_p}{\partial\mathbf{n}} = \mathbf{f_s}\cdot\mathbf{n} \]

Definition at line 404 of file SmoothedProfileMethod.cpp.

405{
406 size_t nvel = m_velocity.size();
407 Array<OneD, ExpListSharedPtr> BndExp;
408 Array<OneD, SpatialDomains::BoundaryConditionShPtr> BndCond;
409
410 // Get the BC expansions
411 BndExp = m_pressureP->GetBndCondExpansions();
412 BndCond = m_pressureP->GetBndConditions();
413
414 // For each boundary...
415 for (size_t b = 0; b < BndExp.size(); ++b)
416 {
417 // Only for BCs based on the derivative
418 if (BndCond[b]->GetBoundaryConditionType() ==
420 BndCond[b]->GetBoundaryConditionType() == SpatialDomains::ePeriodic)
421 {
422 // Calculate f_s values
423 Array<OneD, Array<OneD, NekDouble>> f_s(nvel);
424 for (size_t i = 0; i < nvel; ++i)
425 {
426 f_s[i] = m_fs[0]->GetBndCondExpansions()[b]->GetPhys();
427 }
428
429 // BC is f_s * n
430 BndExp[b]->NormVectorIProductWRTBase(f_s, BndExp[b]->UpdatePhys());
431 }
432 }
433}
MultiRegions::ExpListSharedPtr m_pressureP
Correction pressure field for SPM.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fs
Forcing function 'f_s'.

References Nektar::SpatialDomains::eNeumann, Nektar::SpatialDomains::ePeriodic, m_fs, m_pressureP, and Nektar::IncNavierStokes::m_velocity.

Referenced by SetUpCorrectionPressure().

◆ SetUpCorrectionPressure()

void Nektar::SmoothedProfileMethod::SetUpCorrectionPressure ( const Array< OneD, const Array< OneD, NekDouble > > &  fields,
Array< OneD, Array< OneD, NekDouble > > &  Forcing 
)
protected

Sets the forcing term of the equation for the correction pressure \(p_p\):

\[ \nabla\cdot\mathbf{f_s} \]

Parameters
fields
Forcing

Definition at line 301 of file SmoothedProfileMethod.cpp.

304{
305 size_t physTot = m_fs[0]->GetNpoints();
306 size_t nvel = m_velocity.size();
307
308 // Set boundary conditions
310
311 // Divergence of 'fs'
312 m_fields[m_velocity[0]]->PhysDeriv(eX, m_fs[0]->GetPhys(), Forcing[0]);
313
314 // Using 'Forcing[1]' as storage
315 for (size_t i = 1; i < nvel; ++i)
316 {
317 size_t ind = m_velocity[i];
318 m_fields[ind]->PhysDeriv(DirCartesianMap[i], m_fs[i]->GetPhys(),
319 Forcing[1]);
320 Vmath::Vadd(physTot, Forcing[1], 1, Forcing[0], 1, Forcing[0], 1);
321 }
322}
void SetCorrectionPressureBCs()
Updates the BCs for boundaries with Neumann or periodic BCs in the pressure:
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:87
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

References Nektar::MultiRegions::DirCartesianMap, Nektar::MultiRegions::eX, Nektar::SolverUtils::EquationSystem::m_fields, m_fs, Nektar::IncNavierStokes::m_velocity, SetCorrectionPressureBCs(), and Vmath::Vadd().

Referenced by v_SolveUnsteadyStokesSystem().

◆ SetUpExpansions()

template<typename T >
void Nektar::SmoothedProfileMethod::SetUpExpansions ( int  nvel)
inlineprotected

Initialises the expansions for the intermediate pressure, the 'Phi' function and the IB forcing.

Definition at line 132 of file SmoothedProfileMethod.h.

133 {
134 int iVel = m_velocity[0];
136 *std::dynamic_pointer_cast<T>(m_pressure));
138 *std::dynamic_pointer_cast<T>(m_fields[iVel]));
139
140 m_fs = Array<OneD, MultiRegions::ExpListSharedPtr>(nvel);
141 for (int i = 0; i < nvel; ++i)
142 {
144 *std::dynamic_pointer_cast<T>(m_fields[iVel]));
145 }
146
147 // Set to wave space if homogeneous case
149 {
150 m_pressureP->SetWaveSpace(true);
151 m_phi->SetWaveSpace(true);
152 for (int i = 0; i < nvel; ++i)
153 {
154 m_fs[i]->SetWaveSpace(true);
155 }
156 }
157 }
MultiRegions::ExpListSharedPtr m_pressure
Pointer to field holding pressure field.
enum HomogeneousType m_HomogeneousType

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::SolverUtils::EquationSystem::eNotHomogeneous, Nektar::SolverUtils::EquationSystem::m_fields, m_fs, Nektar::SolverUtils::EquationSystem::m_HomogeneousType, m_phi, Nektar::IncNavierStokes::m_pressure, m_pressureP, and Nektar::IncNavierStokes::m_velocity.

◆ SolveCorrectedVelocity()

void Nektar::SmoothedProfileMethod::SolveCorrectedVelocity ( Array< OneD, Array< OneD, NekDouble > > &  Forcing,
Array< OneD, Array< OneD, NekDouble > > &  fields,
NekDouble  dt 
)
protected

Corrects the velocity field so that the IBs are taken into account. Solves the explicit equation:

\[ \frac{\gamma_0(\mathbf{u_p}^{n+1} - \mathbf{u}^*)}{\Delta t} = \mathbf{f_s} - \nabla p_p \]

Parameters
Forcing
fields
dt

Definition at line 357 of file SmoothedProfileMethod.cpp.

360{
361 size_t physTot = m_phi->GetNpoints();
362
363 // Gradient of p_p
364 size_t nvel = m_velocity.size();
365 if (nvel == 2)
366 {
367 m_pressureP->PhysDeriv(m_pressureP->GetPhys(), Forcing[0], Forcing[1]);
368 }
369 else
370 {
371 m_pressureP->PhysDeriv(m_pressureP->GetPhys(), Forcing[0], Forcing[1],
372 Forcing[2]);
373 }
374
375 // Velocity correction
376 for (size_t i = 0; i < nvel; ++i)
377 {
378 // Adding -(1-m_phi)*grad(p_p) instead of -grad(p_p) reduces the
379 // flux through the walls, but the flow is not incompressible
380 if (m_session->DefinesSolverInfo("ForceBoundary") &&
381 boost::iequals(m_session->GetSolverInfo("ForceBoundary"), "True"))
382 {
383 Vmath::Vvtvm(physTot, m_phi->GetPhys(), 1, Forcing[i], 1,
384 Forcing[i], 1, Forcing[i], 1);
385 Vmath::Vadd(physTot, m_fs[i]->GetPhys(), 1, Forcing[i], 1,
386 Forcing[i], 1);
387 }
388 else
389 {
390 Vmath::Vsub(physTot, m_fs[i]->GetPhys(), 1, Forcing[i], 1,
391 Forcing[i], 1);
392 }
393 Blas::Daxpy(physTot, dt / m_gamma0, Forcing[i], 1, fields[i], 1);
394 }
395}
NekDouble m_gamma0
Stiffly-stable scheme coefficient.
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
void Vvtvm(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)
vvtvm (vector times vector minus vector): z = w*x - y
Definition: Vmath.hpp:381
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 Blas::Daxpy(), m_fs, m_gamma0, m_phi, m_pressureP, Nektar::SolverUtils::EquationSystem::m_session, Nektar::IncNavierStokes::m_velocity, Vmath::Vadd(), Vmath::Vsub(), and Vmath::Vvtvm().

Referenced by v_SolveUnsteadyStokesSystem().

◆ SolveCorrectionPressure()

void Nektar::SmoothedProfileMethod::SolveCorrectionPressure ( const Array< OneD, NekDouble > &  Forcing)
protected

Solves the Poisson equation for the correction pressure \(p_p\):

\[ \nabla^2 p_p = \nabla\cdot\mathbf{f_s} \]

Parameters
Forcing

Definition at line 332 of file SmoothedProfileMethod.cpp.

334{
336 // Factor 'lambda=0' in Helmholtz equation to get the Poisson form
338
339 // Solve the Poisson equation
340 m_pressureP->HelmSolve(Forcing, m_pressureP->UpdateCoeffs(), factors);
341
342 // Update node values from coefficients
343 m_pressureP->BwdTrans(m_pressureP->GetCoeffs(), m_pressureP->UpdatePhys());
344}
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:430
StdRegions::ConstFactorMap factors

References Nektar::StdRegions::eFactorLambda, Nektar::VarcoeffHashingTest::factors, and m_pressureP.

Referenced by v_SolveUnsteadyStokesSystem().

◆ UpdateForcing()

void Nektar::SmoothedProfileMethod::UpdateForcing ( const Array< OneD, const Array< OneD, NekDouble > > &  fields,
NekDouble  dt 
)
protected

For a body with a velocity \(\mathbf{u_p}\), the force \(\mathbf{f_s}\) applied to the fluid ensures that the IBC are met:

\[ \mathbf{f_s} = \frac{\Phi^{n+1}\left(\mathbf{u_p}^{n+1} - \mathbf{u^*}\right)}{\Delta t} \]

Parameters
fields
dt
f_s

Definition at line 484 of file SmoothedProfileMethod.cpp.

486{
487 size_t nvel = m_velocity.size();
488 size_t nq = m_phi->GetNpoints();
489
490 for (size_t i = 0; i < nvel; ++i)
491 {
492 // In homogeneous cases, switch out of wave space
493 Array<OneD, NekDouble> tmpField(nq);
494 size_t ind = m_velocity[i];
495
497 m_fields[ind]->GetWaveSpace())
498 {
499 m_fields[ind]->HomogeneousBwdTrans(nq, fields[i], tmpField);
500 m_fs[i]->HomogeneousBwdTrans(nq, m_fs[i]->GetPhys(),
501 m_fs[i]->UpdatePhys());
502 }
503 else
504 {
505 tmpField = fields[i];
506 }
507
508 Vmath::Vsub(nq, m_up[i], 1, tmpField, 1, m_fs[i]->UpdatePhys(), 1);
509 Vmath::Vmul(nq, m_phi->GetPhys(), 1, m_fs[i]->GetPhys(), 1,
510 m_fs[i]->UpdatePhys(), 1);
511 Vmath::Smul(nq, m_gamma0 / dt, m_fs[i]->GetPhys(), 1,
512 m_fs[i]->UpdatePhys(), 1);
513
514 // And go back to wave space if the 'if' was executed
516 m_fields[ind]->GetWaveSpace())
517 {
518 m_fs[i]->HomogeneousFwdTrans(nq, m_fs[i]->GetPhys(),
519 m_fs[i]->UpdatePhys());
520 }
521 }
522}
Array< OneD, Array< OneD, NekDouble > > m_up
Velocity of the immersed body(ies)
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
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

References Nektar::SolverUtils::EquationSystem::eNotHomogeneous, Nektar::SolverUtils::EquationSystem::m_fields, m_fs, m_gamma0, Nektar::SolverUtils::EquationSystem::m_HomogeneousType, m_phi, m_up, Nektar::IncNavierStokes::m_velocity, Vmath::Smul(), Vmath::Vmul(), and Vmath::Vsub().

Referenced by v_SolveUnsteadyStokesSystem().

◆ UpdatePhiUp()

void Nektar::SmoothedProfileMethod::UpdatePhiUp ( NekDouble  t)
protected

Calculates the values of the shape function.

Parameters
t

Definition at line 440 of file SmoothedProfileMethod.cpp.

441{
442 // Initialise 'm_up' and 'm_phi' during first step
443 if (t <= 0.0)
444 {
445 if (!m_filePhi)
446 {
447 // Update 'm_phi' only if it was provided as a function
448 m_phiEvaluator->Evaluate("Phi", m_phi->UpdatePhys(), t);
449 }
450
451 // Initialize both variables for the first step
452 m_phiEvaluator->Evaluate(m_velName, m_up, t);
453
454 // Initialise 'm_upPrev' in all cases
455 m_upPrev = m_up;
456 }
457 // If timedependent 'm_phi'
458 // Phi functions from files are not timedependent
459 else if (m_timeDependentPhi)
460 {
461 m_phiEvaluator->Evaluate("Phi", m_phi->UpdatePhys(), t);
462
463 // And if velocities are timedependent as well
465 {
466 // Store previous value of u_p during simulation
467 m_upPrev = m_up;
468 m_phiEvaluator->Evaluate(m_velName, m_up, t);
469 }
470 }
471}
std::vector< std::string > m_velName
Vector storing the names of the components of \u_p.
Array< OneD, Array< OneD, NekDouble > > m_upPrev

References m_filePhi, m_phi, m_phiEvaluator, m_timeDependentPhi, m_timeDependentUp, m_up, m_upPrev, and m_velName.

Referenced by v_InitObject(), and v_SolveUnsteadyStokesSystem().

◆ v_GenerateSummary()

void Nektar::SmoothedProfileMethod::v_GenerateSummary ( SolverUtils::SummaryList s)
overridevirtual

Generates the summary of the current simulation.

Parameters
s

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 233 of file SmoothedProfileMethod.cpp.

234{
236 SolverUtils::AddSummaryItem(s, "IB formulation",
237 "Smoothed Profile Method (SPM)");
238}
void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print a summary of time stepping parameters.
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:47

References Nektar::SolverUtils::AddSummaryItem(), and Nektar::VelocityCorrectionScheme::v_GenerateSummary().

◆ v_InitObject()

void Nektar::SmoothedProfileMethod::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 78 of file SmoothedProfileMethod.cpp.

79{
81
82 // Update implicit time-intregration class operators
85
86 // Number of dims as number of velocity vectors
87 size_t nvel = m_velocity.size();
88
89 // Initialization of correction pressure and shape function
90 switch (nvel)
91 {
92 case 1:
94 {
95 SetUpExpansions<ContField>(nvel);
96 }
98 {
99 SetUpExpansions<DisContField>(nvel);
100 }
101 break;
102
103 case 2:
105 {
106 SetUpExpansions<ContField>(nvel);
107 }
109 {
110 SetUpExpansions<DisContField>(nvel);
111 }
112 break;
113
114 case 3:
116 {
118 {
119 SetUpExpansions<ContField>(nvel);
120 }
122 {
123 SetUpExpansions<ContField3DHomogeneous1D>(nvel);
124 }
127 {
128 SetUpExpansions<ContField3DHomogeneous2D>(nvel);
129 }
130 }
132 {
134 {
135 SetUpExpansions<DisContField>(nvel);
136 }
138 {
139 SetUpExpansions<DisContField3DHomogeneous1D>(nvel);
140 }
143 {
144 SetUpExpansions<DisContField3DHomogeneous2D>(nvel);
145 }
146 }
147 break;
148 }
149
150 // Read 'm_phi' and its velocity
151 ASSERTL0(m_session->DefinesFunction("ShapeFunction"),
152 "ShapeFunction must be defined in the session file.")
153 ReadPhi();
154
155 // Allocate the vector 'm_up'
156 size_t physTot = m_phi->GetTotPoints();
157 m_velName.push_back("Up");
158 if (nvel > 1)
159 {
160 m_velName.push_back("Vp");
161 }
162 if (nvel == 3)
163 {
164 m_velName.push_back("Wp");
165 }
166
167 m_up = Array<OneD, Array<OneD, NekDouble>>(nvel);
168 for (size_t i = 0; i < nvel; ++i)
169 {
170 m_up[i] = Array<OneD, NekDouble>(physTot, 0.0);
171 }
172
173 // Make sure that m_phi and m_up are defined
174 UpdatePhiUp(0.0);
175
176 // Get the time integration scheme.
177 LibUtilities::TimeIntScheme timeInt;
178 if (m_session->DefinesTimeIntScheme())
179 {
180 timeInt = m_session->GetTimeIntScheme();
181 }
182 else
183 {
184 timeInt.method = m_session->GetSolverInfo("TimeIntegrationMethod");
185 timeInt.order = timeInt.method.back() - '0';
186
187 // Remove everything past the IMEX.
188 timeInt.method = timeInt.method.substr(0, 4);
189 }
190
191 // Select 'm_gamma0' depending on IMEX order
192 ASSERTL0(
193 boost::iequals(timeInt.method, "IMEX") && 1 <= timeInt.order &&
194 timeInt.order <= 4,
195 "The TimeIntegrationMethod scheme must be IMEX with order '1' to '4'.")
196
197 switch (timeInt.order)
198 {
199 case 1:
200 m_gamma0 = 1.0;
201 break;
202
203 case 2:
204 m_gamma0 = 3.0 / 2.0;
205 break;
206
207 case 3:
208 m_gamma0 = 11.0 / 6.0;
209 break;
210
211 case 4:
212 m_gamma0 = 25.0 / 12.0;
213 break;
214 }
215
216 // Check if the aeroforces filter is active, negative if inactive
217 m_forcesFilter = -1;
218 for (size_t i = 0; i < m_session->GetFilters().size(); ++i)
219 {
220 if (boost::iequals(m_session->GetFilters()[i].first, "AeroForcesSPM"))
221 {
222 m_forcesFilter = i;
223 break;
224 }
225 }
226}
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
void v_SolveUnsteadyStokesSystem(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, NekDouble time, NekDouble a_iixDt) override
Linear terms due to pressure and visosity are calculated here. After solving the velocity filed witho...
int m_forcesFilter
Position of "AeroForcesSPM" filter in 'm_session->GetFilters()'.
void UpdatePhiUp(NekDouble time)
Calculates the values of the shape function.
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SOLVER_UTILS_EXPORT int GetTotPoints()
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
void v_InitObject(bool DeclareField=true) override
Initialisation object for EquationSystem.

References ASSERTL0, Nektar::LibUtilities::TimeIntegrationSchemeOperators::DefineImplicitSolve(), Nektar::MultiRegions::eDiscontinuous, Nektar::MultiRegions::eGalerkin, Nektar::SolverUtils::EquationSystem::eHomogeneous1D, Nektar::SolverUtils::EquationSystem::eHomogeneous2D, Nektar::SolverUtils::EquationSystem::eHomogeneous3D, Nektar::SolverUtils::EquationSystem::eNotHomogeneous, m_forcesFilter, m_gamma0, Nektar::SolverUtils::EquationSystem::m_HomogeneousType, Nektar::SolverUtils::UnsteadySystem::m_ode, m_phi, Nektar::SolverUtils::EquationSystem::m_projectionType, Nektar::SolverUtils::EquationSystem::m_session, m_up, m_velName, Nektar::IncNavierStokes::m_velocity, Nektar::LibUtilities::TimeIntScheme::method, Nektar::LibUtilities::TimeIntScheme::order, ReadPhi(), UpdatePhiUp(), Nektar::VelocityCorrectionScheme::v_InitObject(), and v_SolveUnsteadyStokesSystem().

◆ v_SolveUnsteadyStokesSystem()

void Nektar::SmoothedProfileMethod::v_SolveUnsteadyStokesSystem ( const Array< OneD, const Array< OneD, NekDouble > > &  inarray,
Array< OneD, Array< OneD, NekDouble > > &  outarray,
NekDouble  time,
NekDouble  a_iixDt 
)
overrideprotectedvirtual

Linear terms due to pressure and visosity are calculated here. After solving the velocity filed without taking into account the immersed boundaries, a new correction is applied through the force \(f_s\):

\[ \mathbf{f_s} = \frac{\Phi^{n+1}(\mathbf{u_p}-\mathbf{u^*})} {\Delta t} \]

Parameters
inarray
outarray
time
a_iixDt

Reimplemented from Nektar::VelocityCorrectionScheme.

Definition at line 254 of file SmoothedProfileMethod.cpp.

258{
260 time, a_iixDt);
261
262 size_t physTot = m_pressureP->GetNpoints();
263
264 /* SPM correction of velocity */
265 // Update 'm_phi' and 'm_up' if needed (evaluated at next time step)
266 UpdatePhiUp(time + a_iixDt);
267 // Update calculation of IB forcing 'm_fs'
268 UpdateForcing(outarray, a_iixDt);
269 // Estimate forces only if requested
270 if (m_forcesFilter >= 0)
271 {
272 static_pointer_cast<FilterAeroForcesSPM>(
274 ->CalculateForces(outarray, m_upPrev, m_phi, time, a_iixDt);
275 }
276 // Set BC conditions for pressure p_p
277 SetUpCorrectionPressure(outarray, m_F);
278 // Solve Poisson equation for pressure p_p
280 // Solve velocity in the next step with IB
281 SolveCorrectedVelocity(m_F, outarray, a_iixDt);
282
283 // Add pressures to get final value
284 Vmath::Vadd(physTot, m_pressure->GetPhys(), 1, m_pressureP->GetPhys(), 1,
285 m_pressure->UpdatePhys(), 1);
286 m_pressure->FwdTrans(m_pressure->GetPhys(), m_pressure->UpdateCoeffs());
287
288 // Add presure to outflow bc if using convective like BCs
289 m_extrapolation->AddPressureToOutflowBCs(m_kinvis);
290}
NekDouble m_kinvis
Kinematic viscosity.
ExtrapolateSharedPtr m_extrapolation
void UpdateForcing(const Array< OneD, const Array< OneD, NekDouble > > &fields, NekDouble dt)
For a body with a velocity , the force applied to the fluid ensures that the IBC are met:
void SolveCorrectedVelocity(Array< OneD, Array< OneD, NekDouble > > &Forcing, Array< OneD, Array< OneD, NekDouble > > &fields, NekDouble dt)
Corrects the velocity field so that the IBs are taken into account. Solves the explicit equation:
void SolveCorrectionPressure(const Array< OneD, NekDouble > &Forcing)
Solves the Poisson equation for the correction pressure :
void SetUpCorrectionPressure(const Array< OneD, const Array< OneD, NekDouble > > &fields, Array< OneD, Array< OneD, NekDouble > > &Forcing)
Sets the forcing term of the equation for the correction pressure :
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
Array< OneD, Array< OneD, NekDouble > > m_F
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)

References Nektar::IncNavierStokes::m_extrapolation, Nektar::VelocityCorrectionScheme::m_F, Nektar::SolverUtils::UnsteadySystem::m_filters, m_forcesFilter, Nektar::IncNavierStokes::m_kinvis, m_phi, Nektar::IncNavierStokes::m_pressure, m_pressureP, m_upPrev, SetUpCorrectionPressure(), SolveCorrectedVelocity(), SolveCorrectionPressure(), UpdateForcing(), UpdatePhiUp(), Nektar::VelocityCorrectionScheme::v_SolveUnsteadyStokesSystem(), and Vmath::Vadd().

Referenced by v_InitObject().

Member Data Documentation

◆ className

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

Name of class.

Definition at line 60 of file SmoothedProfileMethod.h.

◆ m_filePhi

bool Nektar::SmoothedProfileMethod::m_filePhi
protected

Flag indicating that phi was defined in a file.

Definition at line 94 of file SmoothedProfileMethod.h.

Referenced by ReadPhi(), and UpdatePhiUp().

◆ m_forcesFilter

int Nektar::SmoothedProfileMethod::m_forcesFilter
protected

Position of "AeroForcesSPM" filter in 'm_session->GetFilters()'.

Definition at line 96 of file SmoothedProfileMethod.h.

Referenced by v_InitObject(), and v_SolveUnsteadyStokesSystem().

◆ m_fs

Array<OneD, MultiRegions::ExpListSharedPtr> Nektar::SmoothedProfileMethod::m_fs
protected

◆ m_gamma0

NekDouble Nektar::SmoothedProfileMethod::m_gamma0
protected

Stiffly-stable scheme \(\gamma_0\) coefficient.

Definition at line 84 of file SmoothedProfileMethod.h.

Referenced by SolveCorrectedVelocity(), UpdateForcing(), and v_InitObject().

◆ m_phi

MultiRegions::ExpListSharedPtr Nektar::SmoothedProfileMethod::m_phi
protected

Shape function 'phi' as expansion list.

Definition at line 88 of file SmoothedProfileMethod.h.

Referenced by ReadPhi(), SetUpExpansions(), SolveCorrectedVelocity(), UpdateForcing(), UpdatePhiUp(), v_InitObject(), and v_SolveUnsteadyStokesSystem().

◆ m_phiEvaluator

SolverUtils::SessionFunctionSharedPtr Nektar::SmoothedProfileMethod::m_phiEvaluator
protected

Function that evaluates the values of \Phi.

Definition at line 90 of file SmoothedProfileMethod.h.

Referenced by ReadPhi(), and UpdatePhiUp().

◆ m_pressureP

MultiRegions::ExpListSharedPtr Nektar::SmoothedProfileMethod::m_pressureP
protected

◆ m_timeDependentPhi

bool Nektar::SmoothedProfileMethod::m_timeDependentPhi
protected

Flag that is true when phi depends on time.

Definition at line 92 of file SmoothedProfileMethod.h.

Referenced by ReadPhi(), and UpdatePhiUp().

◆ m_timeDependentUp

bool Nektar::SmoothedProfileMethod::m_timeDependentUp
protected

Flag signaling if \(\u_p\) depends on time.

Definition at line 82 of file SmoothedProfileMethod.h.

Referenced by ReadPhi(), and UpdatePhiUp().

◆ m_up

Array<OneD, Array<OneD, NekDouble> > Nektar::SmoothedProfileMethod::m_up
protected

Velocity of the immersed body(ies)

Definition at line 77 of file SmoothedProfileMethod.h.

Referenced by UpdateForcing(), UpdatePhiUp(), and v_InitObject().

◆ m_upPrev

Array<OneD, Array<OneD, NekDouble> > Nektar::SmoothedProfileMethod::m_upPrev
protected

Definition at line 78 of file SmoothedProfileMethod.h.

Referenced by UpdatePhiUp(), and v_SolveUnsteadyStokesSystem().

◆ m_velName

std::vector<std::string> Nektar::SmoothedProfileMethod::m_velName
protected

Vector storing the names of the components of \u_p.

Definition at line 80 of file SmoothedProfileMethod.h.

Referenced by UpdatePhiUp(), and v_InitObject().

◆ solverTypeLookupId

string Nektar::SmoothedProfileMethod::solverTypeLookupId
staticprotected
Initial value:
=
"SolverType", "SmoothedProfileMethod", eSmoothedProfileMethod)
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.
@ eSmoothedProfileMethod

Definition at line 98 of file SmoothedProfileMethod.h.