Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | 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...
 
virtual ~SmoothedProfileMethod ()
 Destroy the Smoothed Profile Method object. More...
 
virtual void v_InitObject ()
 Init object for UnsteadySystem class. More...
 
virtual void v_GenerateSummary (SolverUtils::SummaryList &s)
 Generates the summary of the current simulation. More...
 
void SolveUnsteadyStokesSystem (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time, const NekDouble a_iixDt)
 
- Public Member Functions inherited from Nektar::VelocityCorrectionScheme
 VelocityCorrectionScheme (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Constructor. More...
 
virtual ~VelocityCorrectionScheme ()
 
void SetUpPressureForcing (const Array< OneD, const Array< OneD, NekDouble > > &fields, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt)
 
void SetUpViscousForcing (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt)
 
void SolvePressure (const Array< OneD, NekDouble > &Forcing)
 
void SolveViscous (const Array< OneD, const Array< OneD, NekDouble > > &Forcing, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble aii_Dt)
 
void SolveUnsteadyStokesSystem (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time, const NekDouble a_iixDt)
 
void EvaluateAdvection_SetPressureBCs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 
- Public Member Functions inherited from Nektar::IncNavierStokes
virtual ~IncNavierStokes ()
 
int GetNConvectiveFields (void)
 
Array< OneD, int > & GetVelocity (void)
 
void AddForcing (const SolverUtils::ForcingSharedPtr &pForce)
 
virtual void GetPressure (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &pressure)
 Extract array with pressure from physfield. More...
 
virtual void GetDensity (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &density)
 Extract array with density from physfield. More...
 
virtual bool HasConstantDensity ()
 
virtual void GetVelocity (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &velocity)
 Extract array with velocity from physfield. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::AdvectionSystem
SOLVER_UTILS_EXPORT AdvectionSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
virtual SOLVER_UTILS_EXPORT ~AdvectionSystem ()
 
SOLVER_UTILS_EXPORT AdvectionSharedPtr GetAdvObject ()
 Returns the advection object held by this instance. More...
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleGetElmtCFLVals (const bool FlagAcousticCFL=true)
 
SOLVER_UTILS_EXPORT NekDouble GetCFLEstimate (int &elmtid)
 
- Public Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
virtual SOLVER_UTILS_EXPORT ~UnsteadySystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT NekDouble GetTimeStep (const Array< OneD, const Array< OneD, NekDouble >> &inarray)
 Calculate the larger time-step mantaining the problem stable. More...
 
SOLVER_UTILS_EXPORT void SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT void SetUpTraceNormals (void)
 
SOLVER_UTILS_EXPORT void InitObject ()
 Initialises the members of this object. More...
 
SOLVER_UTILS_EXPORT void DoInitialise ()
 Perform any initialisation necessary before solving the problem. More...
 
SOLVER_UTILS_EXPORT void DoSolve ()
 Solve the problem. More...
 
SOLVER_UTILS_EXPORT void TransCoeffToPhys ()
 Transform from coefficient to physical space. More...
 
SOLVER_UTILS_EXPORT void TransPhysToCoeff ()
 Transform from physical to coefficient space. More...
 
SOLVER_UTILS_EXPORT void Output ()
 Perform output operations after solve. More...
 
SOLVER_UTILS_EXPORT NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation. More...
 
SOLVER_UTILS_EXPORT std::string GetSessionName ()
 Get Session name. More...
 
template<class T >
std::shared_ptr< T > as ()
 
SOLVER_UTILS_EXPORT void ResetSessionName (std::string newname)
 Reset Session name. More...
 
SOLVER_UTILS_EXPORT LibUtilities::SessionReaderSharedPtr GetSession ()
 Get Session name. More...
 
SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr GetPressure ()
 Get pressure field if available. More...
 
SOLVER_UTILS_EXPORT void ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 
SOLVER_UTILS_EXPORT void PrintSummary (std::ostream &out)
 Print a summary of parameters and solver characteristics. More...
 
SOLVER_UTILS_EXPORT void SetLambda (NekDouble lambda)
 Set parameter m_lambda. More...
 
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction (std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
 Get a SessionFunction by name. More...
 
SOLVER_UTILS_EXPORT void SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 Initialise the data in the dependent fields. More...
 
SOLVER_UTILS_EXPORT void EvaluateExactSolution (int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 Evaluates an exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln, bool Normalised=false)
 Compute the L2 error between fields and a given exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, bool Normalised=false)
 Compute the L2 error of the fields. More...
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleErrorExtraPoints (unsigned int field)
 Compute error (L2 and L_inf) over an larger set of quadrature points return [L2 Linf]. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n)
 Write checkpoint file of m_fields. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 Write checkpoint file of custom data fields. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_BaseFlow (const int n)
 Write base flow file of m_fields. More...
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname)
 Write field data to the given filename. More...
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 Write input fields to the given filename. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 Input field data from the given file. More...
 
SOLVER_UTILS_EXPORT void ImportFldToMultiDomains (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int ndomains)
 Input field data from the given file to multiple domains. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, std::vector< std::string > &fieldStr, Array< OneD, Array< OneD, NekDouble > > &coeffs)
 Output a field. Input field data into array from the given file. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, MultiRegions::ExpListSharedPtr &pField, std::string &pFieldName)
 Output a field. Input field data into ExpList from the given file. More...
 
SOLVER_UTILS_EXPORT void SessionSummary (SummaryList &vSummary)
 Write out a session summary. More...
 
SOLVER_UTILS_EXPORT Array< OneD, MultiRegions::ExpListSharedPtr > & UpdateFields ()
 
SOLVER_UTILS_EXPORT LibUtilities::FieldMetaDataMapUpdateFieldMetaDataMap ()
 Get hold of FieldInfoMap so it can be updated. More...
 
SOLVER_UTILS_EXPORT NekDouble GetFinalTime ()
 Return final time. More...
 
SOLVER_UTILS_EXPORT int GetNcoeffs ()
 
SOLVER_UTILS_EXPORT int GetNcoeffs (const int eid)
 
SOLVER_UTILS_EXPORT int GetNumExpModes ()
 
SOLVER_UTILS_EXPORT const Array< OneD, int > GetNumExpModesPerExp ()
 
SOLVER_UTILS_EXPORT int GetNvariables ()
 
SOLVER_UTILS_EXPORT const std::string GetVariable (unsigned int i)
 
SOLVER_UTILS_EXPORT int GetTraceTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTraceNpoints ()
 
SOLVER_UTILS_EXPORT int GetExpSize ()
 
SOLVER_UTILS_EXPORT int GetPhys_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetCoeff_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTotPoints (int n)
 
SOLVER_UTILS_EXPORT int GetNpoints ()
 
SOLVER_UTILS_EXPORT int GetSteps ()
 
SOLVER_UTILS_EXPORT NekDouble GetTimeStep ()
 
SOLVER_UTILS_EXPORT void CopyFromPhysField (const int i, Array< OneD, NekDouble > &output)
 
SOLVER_UTILS_EXPORT void CopyToPhysField (const int i, Array< OneD, NekDouble > &output)
 
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 void SetTime (const NekDouble time)
 
SOLVER_UTILS_EXPORT void SetInitialStep (const int step)
 
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time. More...
 
virtual SOLVER_UTILS_EXPORT bool v_NegatedOp ()
 Virtual function to identify if operator is negated in DoSolve. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::FluidInterface
virtual ~FluidInterface ()=default
 

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

Protected Member Functions

virtual void v_SolveUnsteadyStokesSystem (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, NekDouble time, NekDouble a_iixDt)
 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...
 
virtual bool v_PostIntegrate (int step)
 
virtual void v_TransCoeffToPhys (void)
 Virtual function for transformation to physical space. More...
 
virtual void v_TransPhysToCoeff (void)
 Virtual function for transformation to coefficient space. More...
 
virtual void v_DoInitialise (void)
 Sets up initial conditions. More...
 
virtual Array< OneD, bool > v_GetSystemSingularChecks ()
 
virtual int v_GetForceDimension ()
 
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, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble aii_Dt)
 
virtual void v_EvaluateAdvection_SetPressureBCs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 
virtual bool v_RequireFwdTrans ()
 
virtual std::string v_GetExtrapolateStr (void)
 
virtual std::string v_GetSubSteppingExtrapolateStr (const std::string &instr)
 
void SetUpSVV (void)
 
void SetUpExtrapolation (void)
 
void SVVVarDiffCoeff (const NekDouble velmag, Array< OneD, NekDouble > &diffcoeff, const Array< OneD, Array< OneD, NekDouble > > &vel=NullNekDoubleArrayOfArray)
 
void AppendSVVFactors (StdRegions::ConstFactorMap &factors, MultiRegions::VarFactorsMap &varFactorsMap)
 
- Protected Member Functions inherited from Nektar::IncNavierStokes
 IncNavierStokes (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Constructor. More...
 
EquationType GetEquationType (void)
 
void EvaluateAdvectionTerms (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 
void WriteModalEnergy (void)
 
void SetBoundaryConditions (NekDouble time)
 time dependent boundary conditions updating More...
 
void SetRadiationBoundaryForcing (int fieldid)
 Set Radiation forcing term. More...
 
void SetZeroNormalVelocity ()
 Set Normal Velocity Component to Zero. More...
 
void SetWomersleyBoundary (const int fldid, const int bndid)
 Set Womersley Profile if specified. More...
 
void SetUpWomersley (const int fldid, const int bndid, std::string womstr)
 Set Up Womersley details. More...
 
virtual MultiRegions::ExpListSharedPtr v_GetPressure ()
 
virtual Array< OneD, NekDoublev_GetMaxStdVelocity (const NekDouble SpeedSoundFactor)
 
virtual bool v_PreIntegrate (int step)
 
- Protected Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
SOLVER_UTILS_EXPORT UnsteadySystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises UnsteadySystem class members. More...
 
SOLVER_UTILS_EXPORT NekDouble MaxTimeStepEstimator ()
 Get the maximum timestep estimator for cfl control. More...
 
virtual SOLVER_UTILS_EXPORT void v_DoSolve ()
 Solves an unsteady problem. More...
 
virtual SOLVER_UTILS_EXPORT void v_AppendOutput1D (Array< OneD, Array< OneD, NekDouble >> &solution1D)
 Print the solution at each solution point in a txt file. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_GetTimeStep (const Array< OneD, const Array< OneD, NekDouble >> &inarray)
 Return the timestep to be used for the next step in the time-marching loop. More...
 
virtual SOLVER_UTILS_EXPORT void v_SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
SOLVER_UTILS_EXPORT void CheckForRestartTime (NekDouble &time, int &nchk)
 
SOLVER_UTILS_EXPORT void SVVVarDiffCoeff (const Array< OneD, Array< OneD, NekDouble >> vel, StdRegions::VarCoeffMap &varCoeffMap)
 Evaluate the SVV diffusion coefficient according to Moura's paper where it should proportional to h time velocity. More...
 
virtual SOLVER_UTILS_EXPORT bool UpdateTimeStepCheck ()
 
- Protected Member Functions inherited from Nektar::SolverUtils::EquationSystem
SOLVER_UTILS_EXPORT EquationSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 Initialises EquationSystem class members. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Virtual function for the L_inf error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT NekDouble v_L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray, bool Normalised=false)
 Virtual function for the L_2 error computation between fields and a given exact solution. More...
 
virtual SOLVER_UTILS_EXPORT void v_SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 
virtual SOLVER_UTILS_EXPORT void v_EvaluateExactSolution (unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 
virtual SOLVER_UTILS_EXPORT void v_Output (void)
 
virtual SOLVER_UTILS_EXPORT void v_ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 

Protected 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...
 
NekDouble m_sVVCutoffRatio
 cutt off ratio from which to start decayhing modes More...
 
NekDouble m_sVVDiffCoeff
 Diffusion coefficient of SVV modes. More...
 
NekDouble m_sVVCutoffRatioHomo1D
 
NekDouble m_sVVDiffCoeffHomo1D
 Diffusion coefficient of SVV modes in homogeneous 1D Direction. More...
 
Array< OneD, NekDoublem_svvVarDiffCoeff
 Array of coefficient if power kernel is used in SVV. More...
 
bool m_IsSVVPowerKernel
 Identifier for Power Kernel otherwise DG kernel. More...
 
Array< OneD, NekDoublem_diffCoeff
 Diffusion coefficients (will be kinvis for velocities) More...
 
StdRegions::VarCoeffMap m_varCoeffLap
 Variable Coefficient map for the Laplacian which can be activated as part of SVV or otherwise. More...
 
NekDouble m_flowrate
 Desired volumetric flowrate. More...
 
NekDouble m_flowrateArea
 Area of the boundary through which we are measuring the flowrate. More...
 
bool m_homd1DFlowinPlane
 
NekDouble m_greenFlux
 Flux of the Stokes function solution. More...
 
NekDouble m_alpha
 Current flowrate correction. More...
 
int m_flowrateBndID
 Boundary ID of the flowrate reference surface. More...
 
int m_planeID
 Plane ID for cases with homogeneous expansion. More...
 
MultiRegions::ExpListSharedPtr m_flowrateBnd
 Flowrate reference surface. More...
 
Array< OneD, Array< OneD, NekDouble > > m_flowrateStokes
 Stokes solution used to impose flowrate. More...
 
std::ofstream m_flowrateStream
 Output stream to record flowrate. More...
 
int m_flowrateSteps
 Interval at which to record flowrate data. More...
 
NekDouble m_flowrateAiidt
 Value of aii_dt used to compute Stokes flowrate solution. More...
 
Array< OneD, Array< OneD, NekDouble > > m_F
 
- Protected Attributes inherited from Nektar::IncNavierStokes
ExtrapolateSharedPtr m_extrapolation
 
std::ofstream m_mdlFile
 modal energy file More...
 
bool m_SmoothAdvection
 bool to identify if advection term smoothing is requested More...
 
std::vector< SolverUtils::ForcingSharedPtrm_forcing
 Forcing terms. More...
 
int m_nConvectiveFields
 Number of fields to be convected;. More...
 
Array< OneD, int > m_velocity
 int which identifies which components of m_fields contains the velocity (u,v,w); More...
 
MultiRegions::ExpListSharedPtr m_pressure
 Pointer to field holding pressure field. More...
 
NekDouble m_kinvis
 Kinematic viscosity. More...
 
int m_energysteps
 dump energy to file at steps time More...
 
EquationType m_equationType
 equation type; More...
 
Array< OneD, Array< OneD, int > > m_fieldsBCToElmtID
 Mapping from BCs to Elmt IDs. More...
 
Array< OneD, Array< OneD, int > > m_fieldsBCToTraceID
 Mapping from BCs to Elmt Edge IDs. More...
 
Array< OneD, Array< OneD, NekDouble > > m_fieldsRadiationFactor
 RHS Factor for Radiation Condition. More...
 
int m_intSteps
 Number of time integration steps AND Order of extrapolation for pressure boundary conditions. More...
 
std::map< int, std::map< int, WomersleyParamsSharedPtr > > m_womersleyParams
 Womersley parameters if required. More...
 
- Protected Attributes inherited from Nektar::SolverUtils::AdvectionSystem
SolverUtils::AdvectionSharedPtr m_advObject
 Advection term. More...
 
- Protected Attributes inherited from Nektar::SolverUtils::UnsteadySystem
int m_infosteps
 Number of time steps between outputting status information. More...
 
int m_abortSteps
 Number of steps between checks for abort conditions. More...
 
int m_filtersInfosteps
 Number of time steps between outputting filters information. More...
 
int m_nanSteps
 
LibUtilities::TimeIntegrationSchemeSharedPtr m_intScheme
 Wrapper to the time integration scheme. More...
 
LibUtilities::TimeIntegrationSchemeOperators m_ode
 The time integration scheme operators to use. More...
 
NekDouble m_epsilon
 
bool m_explicitDiffusion
 Indicates if explicit or implicit treatment of diffusion is used. More...
 
bool m_explicitAdvection
 Indicates if explicit or implicit treatment of advection is used. More...
 
bool m_explicitReaction
 Indicates if explicit or implicit treatment of reaction is used. More...
 
bool m_homoInitialFwd
 Flag to determine if simulation should start in homogeneous forward transformed state. More...
 
NekDouble m_steadyStateTol
 Tolerance to which steady state should be evaluated at. More...
 
int m_steadyStateSteps
 Check for steady state at step interval. More...
 
NekDouble m_steadyStateRes = 1.0
 
NekDouble m_steadyStateRes0 = 1.0
 
Array< OneD, Array< OneD, NekDouble > > m_previousSolution
 Storage for previous solution for steady-state check. More...
 
std::ofstream m_errFile
 
std::vector< int > m_intVariables
 
std::vector< std::pair< std::string, FilterSharedPtr > > m_filters
 
NekDouble m_filterTimeWarning
 Number of time steps between outputting status information. More...
 
NekDouble m_TimeIntegLambda =0.0
 coefff of spacial derivatives(rhs or m_F in GLM) in calculating the residual of the whole equation(used in unsteady time integrations) More...
 
bool m_flagImplicitItsStatistics
 
bool m_flagImplicitSolver = false
 
Array< OneD, NekDoublem_magnitdEstimat
 estimate the magnitude of each conserved varibles More...
 
Array< OneD, NekDoublem_locTimeStep
 local time step(notice only for jfnk other see m_cflSafetyFactor) More...
 
NekDouble m_inArrayNorm =-1.0
 
int m_TotLinItePerStep =0
 
int m_StagesPerStep =1
 
bool m_flagUpdatePreconMat
 
int m_maxLinItePerNewton
 
int m_TotNewtonIts =0
 
int m_TotLinIts =0
 
int m_TotImpStages =0
 
bool m_CalcPhysicalAV = true
 flag to update artificial viscosity More...
 
- Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
LibUtilities::CommSharedPtr m_comm
 Communicator. More...
 
bool m_verbose
 
bool m_root
 
LibUtilities::SessionReaderSharedPtr m_session
 The session reader. More...
 
std::map< std::string, SolverUtils::SessionFunctionSharedPtrm_sessionFunctions
 Map of known SessionFunctions. More...
 
LibUtilities::FieldIOSharedPtr m_fld
 Field input/output. More...
 
Array< OneD, MultiRegions::ExpListSharedPtrm_fields
 Array holding all dependent variables. More...
 
SpatialDomains::BoundaryConditionsSharedPtr m_boundaryConditions
 Pointer to boundary conditions object. More...
 
SpatialDomains::MeshGraphSharedPtr m_graph
 Pointer to graph defining mesh. More...
 
std::string m_sessionName
 Name of the session. More...
 
NekDouble m_time
 Current time of simulation. More...
 
int m_initialStep
 Number of the step where the simulation should begin. More...
 
NekDouble m_fintime
 Finish time of the simulation. More...
 
NekDouble m_timestep
 Time step size. More...
 
NekDouble m_timestepMax = -1.0
 Time step size. More...
 
NekDouble m_lambda
 Lambda constant in real system if one required. More...
 
NekDouble m_checktime
 Time between checkpoints. More...
 
NekDouble m_lastCheckTime
 
NekDouble m_TimeIncrementFactor
 
int m_nchk
 Number of checkpoints written so far. More...
 
int m_steps
 Number of steps to take. More...
 
int m_checksteps
 Number of steps between checkpoints. More...
 
int m_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...
 
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...
 

Additional Inherited Members

- Public Attributes inherited from Nektar::SolverUtils::UnsteadySystem
NekDouble m_cflSafetyFactor
 CFL safety factor (comprise between 0 to 1). More...
 
NekDouble m_cflNonAcoustic
 
NekDouble m_CFLGrowth
 CFL growth rate. More...
 
NekDouble m_CFLEnd
 maximun cfl in cfl growth More...
 
- Protected Types inherited from Nektar::SolverUtils::EquationSystem
enum  HomogeneousType { eHomogeneous1D , eHomogeneous2D , eHomogeneous3D , eNotHomogeneous }
 Parameter for homogeneous expansions. More...
 
- Static Protected Attributes inherited from Nektar::SolverUtils::EquationSystem
static std::string equationSystemTypeLookupIds []
 

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 59 of file SmoothedProfileMethod.cpp.

62  : UnsteadySystem(pSession, pGraph),
63  VelocityCorrectionScheme(pSession, pGraph)
64  {
65 
66  }
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  )
virtual

Destroy the Smoothed Profile Method object.

Definition at line 72 of file SmoothedProfileMethod.cpp.

73  {
74 
75  }

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  pSession, 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 593 of file SmoothedProfileMethod.cpp.

594  {
595  // Get the handler of first function block
596  TiXmlElement *conds = m_session->GetElement("Nektar/Conditions");
597  TiXmlElement *function = conds->FirstChildElement("FUNCTION");
598 
599  // Loop over functions until the block 'name' is found
600  string functionType = function->Attribute("NAME");
601  while (function && !boost::iequals(functionType, functionName))
602  {
603  function = function->NextSiblingElement("FUNCTION");
604  functionType = function->Attribute("NAME");
605  }
606 
607  return function;
608  }
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 552 of file SmoothedProfileMethod.cpp.

554  {
555  // Get the handler of the function
556  TiXmlElement *function = GetFunctionHdl(funcName);
557 
558  // Go to the first element
559  TiXmlElement *functionDef = function->FirstChildElement();
560  ASSERTL0(functionDef, "At least one element must be defined in " +
561  funcName)
562 
563  // And search the element with name 'elemName'
564  string varName = functionDef->Attribute("VAR");
565  while(functionDef && !boost::iequals(varName, elemName))
566  {
567  functionDef = functionDef->NextSiblingElement();
568  varName = functionDef->Attribute("VAR");
569  }
570 
571  ASSERTL0(functionDef, "Variable " + elemName + " must be defined in " +
572  funcName + ".");
573 
574  // And return the value of USERDEFINEDTYPE
575  string attr;
576  int err = functionDef->QueryStringAttribute("USERDEFINEDTYPE", &attr);
577  bool output = boost::iequals(attr, "TimeDependent");
578 
579  ASSERTL0((err == TIXML_NO_ATTRIBUTE) ||
580  (err == TIXML_SUCCESS && output), "USERDEFINEDTYPE in " +
581  elemName + " must be TimeDependent if defined");
582 
583  return output;
584  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
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 610 of file SmoothedProfileMethod.cpp.

611  {
612  // Function evaluator for Phi and Up
613  m_phiEvaluator = GetFunction("ShapeFunction");
614 
615  TiXmlElement *function = GetFunctionHdl("ShapeFunction");
616  TiXmlElement *child = function->FirstChildElement();
617  m_filePhi = false;
618 
619  // If defined by using a file
620  if (boost::iequals(child->ValueStr(), "F"))
621  {
622  // Get name of STL file
623  string fileName;
624  int status = child->QueryStringAttribute("FILE", &fileName);
625  ASSERTL0(status == TIXML_SUCCESS, "An FLD file with the values "
626  "of the phi function has to be supplied.")
627  ASSERTL0(boost::iequals(fileName.substr(fileName.length()-4),
628  ".fld"), "A valid FLD file must be supplied in the "
629  "'ShapeFunction' field.")
630 
631  // Get phi values from XML file (after "FieldConvert" the STL file)
632  // First, load the data
633  std::vector<LibUtilities::FieldDefinitionsSharedPtr> fieldDef;
634  std::vector<std::vector<NekDouble> > fieldData;
635  LibUtilities::FieldMetaDataMap fieldMetaData;
636  LibUtilities::FieldIOSharedPtr phiFile =
637  LibUtilities::FieldIO::CreateForFile(m_session, fileName);
638  phiFile->Import(fileName, fieldDef, fieldData, fieldMetaData);
639 
640  // Only Phi field should be defined in the file
641  ASSERTL0(fieldData.size() == 1, "Only one field (phi) must be "
642  "defined in the FLD file.")
643 
644  // Extract Phi field to output
645  string tmp("phi");
646  m_phi->ExtractDataToCoeffs(fieldDef[0], fieldData[0],
647  tmp, m_phi->UpdateCoeffs());
648  m_phi->BwdTrans(m_phi->GetCoeffs(), m_phi->UpdatePhys());
649  m_filePhi = true;
650  m_timeDependentPhi = false;
651  m_timeDependentUp = false;
652  }
653  else
654  {
655  // Check if Phi is timedependent
656  m_timeDependentPhi = GetVarTimeDependence("ShapeFunction", "Phi");
657 
658  // If so, check if its velocity changes as well
659  m_timeDependentUp = GetVarTimeDependence("ShapeFunction", "Up");
660  switch (m_velocity.size())
661  {
662  case 3:
664  GetVarTimeDependence("ShapeFunction", "Wp");
665  case 2:
667  GetVarTimeDependence("ShapeFunction", "Vp");
668  }
669  }
670  }
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:293
std::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:306
std::map< std::string, std::string > FieldMetaDataMap
Definition: FieldIO.h:52
std::shared_ptr< FieldDefinitions > FieldDefinitionsSharedPtr
Definition: FieldIO.h:179
double NekDouble

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 421 of file SmoothedProfileMethod.cpp.

422  {
423  int nvel = m_velocity.size();
424  Array<OneD, ExpListSharedPtr> BndExp;
425  Array<OneD, SpatialDomains::BoundaryConditionShPtr> BndCond;
426 
427  // Get the BC expansions
428  BndExp = m_pressureP->GetBndCondExpansions();
429  BndCond = m_pressureP->GetBndConditions();
430 
431  // For each boundary...
432  for (int b = 0; b < BndExp.size(); ++b)
433  {
434  // Only for BCs based on the derivative
435  if (BndCond[b]->GetBoundaryConditionType() ==
437  BndCond[b]->GetBoundaryConditionType() ==
439  {
440  // Calculate f_s values
441  Array<OneD, Array<OneD, NekDouble> > f_s(nvel);
442  for (int i = 0; i < nvel; ++i)
443  {
444  f_s[i] = m_fs[0]->GetBndCondExpansions()[b]->GetPhys();
445  }
446 
447  // BC is f_s * n
448  BndExp[b]->NormVectorIProductWRTBase(f_s,
449  BndExp[b]->UpdatePhys());
450  }
451  }
452  }
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 311 of file SmoothedProfileMethod.cpp.

314  {
315  int physTot = m_fs[0]->GetNpoints();
316  int nvel = m_velocity.size();
317 
318  // Set boundary conditions
320 
321  // Divergence of 'fs'
322  m_fields[m_velocity[0]]->PhysDeriv(eX, m_fs[0]->GetPhys(), Forcing[0]);
323 
324  // Using 'Forcing[1]' as storage
325  for (int i = 1; i < nvel; ++i)
326  {
327  int ind = m_velocity[i];
328  m_fields[ind]->PhysDeriv(DirCartesianMap[i],
329  m_fs[i]->GetPhys(), Forcing[1]);
330  Vmath::Vadd(physTot, Forcing[1], 1, Forcing[0], 1, Forcing[0], 1);
331  }
332  }
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:90
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:322

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 148 of file SmoothedProfileMethod.h.

149  {
150  int iVel = m_velocity[0];
152  *std::dynamic_pointer_cast<T>(m_pressure));
154  *std::dynamic_pointer_cast<T>(m_fields[iVel]));
155 
156  m_fs = Array<OneD, MultiRegions::ExpListSharedPtr>(nvel);
157  for (int i = 0; i < nvel; ++i)
158  {
160  *std::dynamic_pointer_cast<T>(m_fields[iVel]));
161  }
162 
163  // Set to wave space if homogeneous case
165  {
166  m_pressureP->SetWaveSpace(true);
167  m_phi->SetWaveSpace(true);
168  for (int i = 0; i < nvel; ++i)
169  {
170  m_fs[i]->SetWaveSpace(true);
171  }
172  }
173  }
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 368 of file SmoothedProfileMethod.cpp.

372  {
373  int physTot = m_phi->GetNpoints();
374 
375  // Gradient of p_p
376  int nvel = m_velocity.size();
377  if (nvel == 2)
378  {
379  m_pressureP->PhysDeriv(m_pressureP->GetPhys(),
380  Forcing[0],
381  Forcing[1]);
382  }
383  else
384  {
385  m_pressureP->PhysDeriv(m_pressureP->GetPhys(),
386  Forcing[0],
387  Forcing[1],
388  Forcing[2]);
389  }
390 
391  // Velocity correction
392  for (int i = 0; i < nvel; ++i)
393  {
394  // Adding -(1-m_phi)*grad(p_p) instead of -grad(p_p) reduces the
395  // flux through the walls, but the flow is not incompressible
396  if (m_session->DefinesSolverInfo("ForceBoundary") &&
397  boost::iequals(m_session->GetSolverInfo("ForceBoundary"),
398  "True"))
399  {
400  Vmath::Vvtvm(physTot, m_phi->GetPhys(), 1, Forcing[i], 1,
401  Forcing[i], 1, Forcing[i], 1);
402  Vmath::Vadd(physTot, m_fs[i]->GetPhys(), 1, Forcing[i], 1,
403  Forcing[i], 1);
404  }
405  else
406  {
407  Vmath::Vsub(physTot, m_fs[i]->GetPhys(), 1, Forcing[i], 1,
408  Forcing[i], 1);
409  }
410  Blas::Daxpy(physTot, dt/m_gamma0, Forcing[i], 1, fields[i], 1);
411  }
412  }
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:167
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 plus vector): z = w*x - y
Definition: Vmath.cpp:541
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:372

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 342 of file SmoothedProfileMethod.cpp.

344  {
346  // Factor 'lambda=0' in Helmholtz equation to get the Poisson form
347  factors[StdRegions::eFactorLambda] = 0.0;
348 
349  // Solve the Poisson equation
350  m_pressureP->HelmSolve(Forcing, m_pressureP->UpdateCoeffs(), factors);
351 
352  // Update node values from coefficients
353  m_pressureP->BwdTrans(m_pressureP->GetCoeffs(),
354  m_pressureP->UpdatePhys());
355  }
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:314

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

Referenced by v_SolveUnsteadyStokesSystem().

◆ SolveUnsteadyStokesSystem()

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

Definition at line 76 of file SmoothedProfileMethod.h.

81  {
82  v_SolveUnsteadyStokesSystem(inarray, outarray, time, a_iixDt);
83  }
virtual void v_SolveUnsteadyStokesSystem(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, NekDouble time, NekDouble a_iixDt)
Linear terms due to pressure and visosity are calculated here. After solving the velocity filed witho...

References v_SolveUnsteadyStokesSystem().

Referenced by v_InitObject().

◆ 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 503 of file SmoothedProfileMethod.cpp.

506  {
507  int nvel = m_velocity.size();
508  int nq = m_phi->GetNpoints();
509 
510  for (int i = 0; i < nvel; ++i)
511  {
512  // In homogeneous cases, switch out of wave space
513  Array<OneD, NekDouble> tmpField(nq);
514  int ind = m_velocity[i];
515 
517  m_fields[ind]->GetWaveSpace())
518  {
519  m_fields[ind]->HomogeneousBwdTrans(fields[i], tmpField);
520  m_fs[i]->HomogeneousBwdTrans(m_fs[i]->GetPhys(),
521  m_fs[i]->UpdatePhys());
522  }
523  else
524  {
525  tmpField = fields[i];
526  }
527 
528  Vmath::Vsub(nq, m_up[i], 1, tmpField, 1, m_fs[i]->UpdatePhys(), 1);
529  Vmath::Vmul(nq, m_phi->GetPhys(), 1, m_fs[i]->GetPhys(), 1,
530  m_fs[i]->UpdatePhys(), 1);
531  Vmath::Smul(nq, m_gamma0/dt, m_fs[i]->GetPhys(), 1,
532  m_fs[i]->UpdatePhys(), 1);
533 
534  // And go back to wave space if the 'if' was executed
536  m_fields[ind]->GetWaveSpace())
537  {
538  m_fs[i]->HomogeneousFwdTrans(m_fs[i]->GetPhys(),
539  m_fs[i]->UpdatePhys());
540  }
541  }
542  }
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.cpp:192
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:225

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 459 of file SmoothedProfileMethod.cpp.

460  {
461  // Initialise 'm_up' and 'm_phi' during first step
462  if (t <= 0.0)
463  {
464  if (!m_filePhi)
465  {
466  // Update 'm_phi' only if it was provided as a function
467  m_phiEvaluator->Evaluate("Phi", m_phi->UpdatePhys(), t);
468  }
469 
470  // Initialize both variables for the first step
471  m_phiEvaluator->Evaluate(m_velName, m_up, t);
472 
473  // Initialise 'm_upPrev' in all cases
474  m_upPrev = m_up;
475  }
476  // If timedependent 'm_phi'
477  // Phi functions from files are not timedependent
478  else if (m_timeDependentPhi)
479  {
480  m_phiEvaluator->Evaluate("Phi", m_phi->UpdatePhys(), t);
481 
482  // And if velocities are timedependent as well
483  if (m_timeDependentUp)
484  {
485  // Store previous value of u_p during simulation
486  m_upPrev = m_up;
487  m_phiEvaluator->Evaluate(m_velName, m_up, t);
488  }
489  }
490  }
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)
virtual

Generates the summary of the current simulation.

Parameters
s

Reimplemented from Nektar::VelocityCorrectionScheme.

Definition at line 238 of file SmoothedProfileMethod.cpp.

239  {
241  SolverUtils::AddSummaryItem(s, "IB formulation",
242  "Smoothed Profile Method (SPM)");
243  }
virtual void v_GenerateSummary(SolverUtils::SummaryList &s)
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 ( )
virtual

Init object for UnsteadySystem class.

Initialization object for UnsteadySystem class.

Reimplemented from Nektar::VelocityCorrectionScheme.

Definition at line 77 of file SmoothedProfileMethod.cpp.

78  {
80 
81  // Update implicit time-intregration class operators
84 
85  // Number of dims as number of velocity vectors
86  int nvel = m_velocity.size();
87 
88  // Initialization of correction pressure and shape function
89  switch (nvel)
90  {
91  case 1:
93  {
94  SetUpExpansions<ContField>(nvel);
95  }
96  else if (m_projectionType == eDiscontinuous)
97  {
98  SetUpExpansions<DisContField>(nvel);
99  }
100  break;
101 
102  case 2:
104  {
105  SetUpExpansions<ContField>(nvel);
106  }
107  else if (m_projectionType == eDiscontinuous)
108  {
109  SetUpExpansions<DisContField>(nvel);
110  }
111  break;
112 
113  case 3:
115  {
117  {
118  SetUpExpansions<ContField>(nvel);
119  }
120  else if (m_HomogeneousType ==
122  {
123  SetUpExpansions<ContField3DHomogeneous1D>(nvel);
124  }
125  else if (m_HomogeneousType ==
129  {
130  SetUpExpansions<ContField3DHomogeneous2D>(nvel);
131  }
132  }
133  else if (m_projectionType == eDiscontinuous)
134  {
136  {
137  SetUpExpansions<DisContField>(nvel);
138  }
139  else if (m_HomogeneousType ==
141  {
142  SetUpExpansions<DisContField3DHomogeneous1D>(nvel);
143  }
144  else if (m_HomogeneousType ==
148  {
149  SetUpExpansions<DisContField3DHomogeneous2D>(nvel);
150  }
151  }
152  break;
153  }
154 
155  // Read 'm_phi' and its velocity
156  ASSERTL0(m_session->DefinesFunction("ShapeFunction"),
157  "ShapeFunction must be defined in the session file.")
158  ReadPhi();
159 
160  // Allocate the vector 'm_up'
161  int physTot = m_phi->GetTotPoints();
162  m_velName.push_back("Up");
163  if (nvel > 1)
164  {
165  m_velName.push_back("Vp");
166  }
167  if (nvel == 3)
168  {
169  m_velName.push_back("Wp");
170  }
171 
172  m_up = Array<OneD, Array<OneD, NekDouble> >(nvel);
173  for (int i = 0; i < nvel; ++i)
174  {
175  m_up[i] = Array<OneD, NekDouble>(physTot, 0.0);
176  }
177 
178  // Make sure that m_phi and m_up are defined
179  UpdatePhiUp(0.0);
180 
181  // Get the time integration scheme.
182  LibUtilities::TimeIntScheme timeInt;
183  if (m_session->DefinesTimeIntScheme())
184  {
185  timeInt = m_session->GetTimeIntScheme();
186  }
187  else
188  {
189  timeInt.method = m_session->GetSolverInfo("TimeIntegrationMethod");
190  timeInt.order = timeInt.method.back() - '0';
191 
192  // Remove everything past the IMEX.
193  timeInt.method = timeInt.method.substr(0,4);
194  }
195 
196  // Select 'm_gamma0' depending on IMEX order
197  ASSERTL0(boost::iequals(timeInt.method, "IMEX") &&
198  1 <= timeInt.order && timeInt.order <= 4,
199  "The TimeIntegrationMethod scheme must be IMEX with order '1' to '4'.")
200 
201  switch (timeInt.order)
202  {
203  case 1:
204  m_gamma0 = 1.0;
205  break;
206 
207  case 2:
208  m_gamma0 = 3.0/2.0;
209  break;
210 
211  case 3:
212  m_gamma0 = 11.0/6.0;
213  break;
214 
215  case 4:
216  m_gamma0 = 25.0/12.0;
217  break;
218  }
219 
220  // Check if the aeroforces filter is active, negative if inactive
221  m_forcesFilter = -1;
222  for (int i = 0; i < m_session->GetFilters().size(); ++i)
223  {
224  if (boost::iequals(m_session->GetFilters()[i].first,
225  "AeroForcesSPM"))
226  {
227  m_forcesFilter = i;
228  break;
229  }
230  }
231  }
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
int m_forcesFilter
Position of "AeroForcesSPM" filter in 'm_session->GetFilters()'.
void UpdatePhiUp(NekDouble time)
Calculates the values of the shape function.
void SolveUnsteadyStokesSystem(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time, const NekDouble a_iixDt)
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.
virtual void v_InitObject()
Init object for UnsteadySystem class.

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(), SolveUnsteadyStokesSystem(), UpdatePhiUp(), and Nektar::VelocityCorrectionScheme::v_InitObject().

◆ 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 
)
protectedvirtual

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

Definition at line 259 of file SmoothedProfileMethod.cpp.

264  {
266  outarray,
267  time,
268  a_iixDt);
269 
270  int physTot = m_pressureP->GetNpoints();
271 
272  /* SPM correction of velocity */
273  // Update 'm_phi' and 'm_up' if needed (evaluated at next time step)
274  UpdatePhiUp(time + a_iixDt);
275  // Update calculation of IB forcing 'm_fs'
276  UpdateForcing(outarray, a_iixDt);
277  // Estimate forces only if requested
278  if (m_forcesFilter >= 0)
279  {
280  static_pointer_cast<FilterAeroForcesSPM>(
281  m_filters[m_forcesFilter].second)->CalculateForces(
282  outarray, m_upPrev, m_phi, time, a_iixDt);
283  }
284  // Set BC conditions for pressure p_p
285  SetUpCorrectionPressure(outarray, m_F);
286  // Solve Poisson equation for pressure p_p
288  // Solve velocity in the next step with IB
289  SolveCorrectedVelocity(m_F, outarray, a_iixDt);
290 
291  // Add pressures to get final value
292  Vmath::Vadd(physTot, m_pressure->GetPhys(), 1,
293  m_pressureP->GetPhys(), 1,
294  m_pressure->UpdatePhys(), 1);
295  m_pressure->FwdTrans(m_pressure->GetPhys(),
296  m_pressure->UpdateCoeffs());
297 
298  // Add presure to outflow bc if using convective like BCs
299  m_extrapolation->AddPressureToOutflowBCs(m_kinvis);
300  }
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
void SolveUnsteadyStokesSystem(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time, const NekDouble a_iixDt)
Array< OneD, Array< OneD, NekDouble > > m_F

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(), Nektar::VelocityCorrectionScheme::SolveUnsteadyStokesSystem(), UpdateForcing(), UpdatePhiUp(), and Vmath::Vadd().

Referenced by SolveUnsteadyStokesSystem().

Member Data Documentation

◆ className

string Nektar::SmoothedProfileMethod::className
static
Initial value:
=
"SmoothedProfileMethod",
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:200
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 106 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 108 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 96 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 100 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 102 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 104 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 94 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 89 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 90 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 92 of file SmoothedProfileMethod.h.

Referenced by UpdatePhiUp(), and v_InitObject().