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

#include <SmoothedProfileMethod.h>

Inheritance diagram for Nektar::SmoothedProfileMethod:
[legend]

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

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

Protected Attributes

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 []
 

Friends

class MemoryManager< SmoothedProfileMethod >
 

Additional Inherited Members

- Public Member Functions inherited from Nektar::VelocityCorrectionScheme
void SetUpPressureForcing (const Array< OneD, const Array< OneD, NekDouble > > &fields, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt)
 
void SetUpViscousForcing (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &Forcing, const NekDouble aii_Dt)
 
void SolvePressure (const Array< OneD, NekDouble > &Forcing)
 
void SolveViscous (const Array< OneD, const Array< OneD, NekDouble > > &Forcing, const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble aii_Dt)
 
void SolveUnsteadyStokesSystem (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time, const NekDouble a_iixDt)
 
void EvaluateAdvection_SetPressureBCs (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
 
- Public Member Functions inherited from Nektar::IncNavierStokes
int GetNConvectiveFields (void)
 
void AddForcing (const SolverUtils::ForcingSharedPtr &pForce)
 
bool DefinedForcing (const std::string &sForce)
 
- Public Member Functions inherited from Nektar::SolverUtils::AdvectionSystem
SOLVER_UTILS_EXPORT AdvectionSystem (const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
 
SOLVER_UTILS_EXPORT ~AdvectionSystem () override=default
 
SOLVER_UTILS_EXPORT AdvectionSharedPtr GetAdvObject ()
 Returns the advection object held by this instance. More...
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleGetElmtCFLVals (const bool FlagAcousticCFL=true)
 
SOLVER_UTILS_EXPORT NekDouble GetCFLEstimate (int &elmtid)
 
- Public Member Functions inherited from Nektar::SolverUtils::UnsteadySystem
SOLVER_UTILS_EXPORT ~UnsteadySystem () override=default
 Destructor. More...
 
SOLVER_UTILS_EXPORT NekDouble GetTimeStep (const Array< OneD, const Array< OneD, NekDouble > > &inarray)
 Calculate the larger time-step mantaining the problem stable. More...
 
SOLVER_UTILS_EXPORT NekDouble GetTimeStep ()
 
SOLVER_UTILS_EXPORT void SetTimeStep (const NekDouble timestep)
 
SOLVER_UTILS_EXPORT void SteadyStateResidual (int step, Array< OneD, NekDouble > &L2)
 
SOLVER_UTILS_EXPORT LibUtilities::TimeIntegrationSchemeSharedPtrGetTimeIntegrationScheme ()
 Returns the time integration scheme. More...
 
SOLVER_UTILS_EXPORT LibUtilities::TimeIntegrationSchemeOperatorsGetTimeIntegrationSchemeOperators ()
 Returns the time integration scheme operators. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::EquationSystem
virtual SOLVER_UTILS_EXPORT ~EquationSystem ()
 Destructor. More...
 
SOLVER_UTILS_EXPORT void InitObject (bool DeclareField=true)
 Initialises the members of this object. More...
 
SOLVER_UTILS_EXPORT void DoInitialise (bool dumpInitialConditions=true)
 Perform any initialisation necessary before solving the problem. More...
 
SOLVER_UTILS_EXPORT void DoSolve ()
 Solve the problem. More...
 
SOLVER_UTILS_EXPORT void TransCoeffToPhys ()
 Transform from coefficient to physical space. More...
 
SOLVER_UTILS_EXPORT void TransPhysToCoeff ()
 Transform from physical to coefficient space. More...
 
SOLVER_UTILS_EXPORT void Output ()
 Perform output operations after solve. More...
 
SOLVER_UTILS_EXPORT std::string GetSessionName ()
 Get Session name. More...
 
template<class T >
std::shared_ptr< T > as ()
 
SOLVER_UTILS_EXPORT void ResetSessionName (std::string newname)
 Reset Session name. More...
 
SOLVER_UTILS_EXPORT LibUtilities::SessionReaderSharedPtr GetSession ()
 Get Session name. More...
 
SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr GetPressure ()
 Get pressure field if available. More...
 
SOLVER_UTILS_EXPORT void ExtraFldOutput (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 
SOLVER_UTILS_EXPORT void PrintSummary (std::ostream &out)
 Print a summary of parameters and solver characteristics. More...
 
SOLVER_UTILS_EXPORT void SetLambda (NekDouble lambda)
 Set parameter m_lambda. More...
 
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction (std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
 Get a SessionFunction by name. More...
 
SOLVER_UTILS_EXPORT void SetInitialConditions (NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
 Initialise the data in the dependent fields. More...
 
SOLVER_UTILS_EXPORT void EvaluateExactSolution (int field, Array< OneD, NekDouble > &outfield, const NekDouble time)
 Evaluates an exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln, bool Normalised=false)
 Compute the L2 error between fields and a given exact solution. More...
 
SOLVER_UTILS_EXPORT NekDouble L2Error (unsigned int field, bool Normalised=false)
 Compute the L2 error of the fields. More...
 
SOLVER_UTILS_EXPORT NekDouble LinfError (unsigned int field, const Array< OneD, NekDouble > &exactsoln=NullNekDouble1DArray)
 Linf error computation. More...
 
SOLVER_UTILS_EXPORT NekDouble H1Error (unsigned int field, const Array< OneD, NekDouble > &exactsoln, bool Normalised=false)
 Compute the H1 error between fields and a given exact solution. More...
 
SOLVER_UTILS_EXPORT Array< OneD, NekDoubleErrorExtraPoints (unsigned int field)
 Compute error (L2 and L_inf) over an larger set of quadrature points return [L2 Linf]. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n)
 Write checkpoint file of m_fields. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_Output (const int n, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 Write checkpoint file of custom data fields. More...
 
SOLVER_UTILS_EXPORT void Checkpoint_BaseFlow (const int n)
 Write base flow file of m_fields. More...
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname)
 Write field data to the given filename. More...
 
SOLVER_UTILS_EXPORT void WriteFld (const std::string &outname, MultiRegions::ExpListSharedPtr &field, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 Write input fields to the given filename. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields)
 Input field data from the given file. More...
 
SOLVER_UTILS_EXPORT void ImportFldToMultiDomains (const std::string &infile, Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const int ndomains)
 Input field data from the given file to multiple domains. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, std::vector< std::string > &fieldStr, Array< OneD, Array< OneD, NekDouble > > &coeffs)
 Output a field. Input field data into array from the given file. More...
 
SOLVER_UTILS_EXPORT void ImportFld (const std::string &infile, MultiRegions::ExpListSharedPtr &pField, std::string &pFieldName)
 Output a field. Input field data into ExpList from the given file. More...
 
SOLVER_UTILS_EXPORT void SessionSummary (SummaryList &vSummary)
 Write out a session summary. More...
 
SOLVER_UTILS_EXPORT Array< OneD, MultiRegions::ExpListSharedPtr > & UpdateFields ()
 
SOLVER_UTILS_EXPORT LibUtilities::FieldMetaDataMapUpdateFieldMetaDataMap ()
 Get hold of FieldInfoMap so it can be updated. More...
 
SOLVER_UTILS_EXPORT NekDouble GetTime ()
 Return final time. More...
 
SOLVER_UTILS_EXPORT int GetNcoeffs ()
 
SOLVER_UTILS_EXPORT int GetNcoeffs (const int eid)
 
SOLVER_UTILS_EXPORT int GetNumExpModes ()
 
SOLVER_UTILS_EXPORT const Array< OneD, int > GetNumExpModesPerExp ()
 
SOLVER_UTILS_EXPORT int GetNvariables ()
 
SOLVER_UTILS_EXPORT const std::string GetVariable (unsigned int i)
 
SOLVER_UTILS_EXPORT int GetTraceTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTraceNpoints ()
 
SOLVER_UTILS_EXPORT int GetExpSize ()
 
SOLVER_UTILS_EXPORT int GetPhys_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetCoeff_Offset (int n)
 
SOLVER_UTILS_EXPORT int GetTotPoints ()
 
SOLVER_UTILS_EXPORT int GetTotPoints (int n)
 
SOLVER_UTILS_EXPORT int GetNpoints ()
 
SOLVER_UTILS_EXPORT int GetSteps ()
 
SOLVER_UTILS_EXPORT void SetSteps (const int steps)
 
SOLVER_UTILS_EXPORT NekDouble GetTimeStep ()
 
SOLVER_UTILS_EXPORT void CopyFromPhysField (const int i, Array< OneD, NekDouble > &output)
 
SOLVER_UTILS_EXPORT void CopyToPhysField (const int i, const Array< OneD, const NekDouble > &input)
 
SOLVER_UTILS_EXPORT Array< OneD, NekDouble > & UpdatePhysField (const int i)
 
SOLVER_UTILS_EXPORT void ZeroPhysFields ()
 
SOLVER_UTILS_EXPORT void FwdTransFields ()
 
SOLVER_UTILS_EXPORT void SetModifiedBasis (const bool modbasis)
 
SOLVER_UTILS_EXPORT int GetCheckpointNumber ()
 
SOLVER_UTILS_EXPORT void SetCheckpointNumber (int num)
 
SOLVER_UTILS_EXPORT int GetCheckpointSteps ()
 
SOLVER_UTILS_EXPORT void SetCheckpointSteps (int num)
 
SOLVER_UTILS_EXPORT int GetInfoSteps ()
 
SOLVER_UTILS_EXPORT void SetInfoSteps (int num)
 
SOLVER_UTILS_EXPORT void SetIterationNumberPIT (int num)
 
SOLVER_UTILS_EXPORT void SetWindowNumberPIT (int num)
 
SOLVER_UTILS_EXPORT Array< OneD, const Array< OneD, NekDouble > > GetTraceNormals ()
 
SOLVER_UTILS_EXPORT void SetTime (const NekDouble time)
 
SOLVER_UTILS_EXPORT void SetTimeStep (const NekDouble timestep)
 
SOLVER_UTILS_EXPORT void SetInitialStep (const int step)
 
SOLVER_UTILS_EXPORT void SetBoundaryConditions (NekDouble time)
 Evaluates the boundary conditions at the given time. More...
 
SOLVER_UTILS_EXPORT bool NegatedOp ()
 Identify if operator is negated in DoSolve. More...
 
- Public Member Functions inherited from Nektar::SolverUtils::ALEHelper
virtual ~ALEHelper ()=default
 
virtual SOLVER_UTILS_EXPORT void v_ALEInitObject (int spaceDim, Array< OneD, MultiRegions::ExpListSharedPtr > &fields)
 
SOLVER_UTILS_EXPORT void InitObject (int spaceDim, Array< OneD, MultiRegions::ExpListSharedPtr > &fields)
 
virtual SOLVER_UTILS_EXPORT void v_UpdateGridVelocity (const NekDouble &time)
 
virtual SOLVER_UTILS_EXPORT void v_ALEPreMultiplyMass (Array< OneD, Array< OneD, NekDouble > > &fields)
 
SOLVER_UTILS_EXPORT void ALEDoElmtInvMass (Array< OneD, Array< OneD, NekDouble > > &traceNormals, Array< OneD, Array< OneD, NekDouble > > &fields, NekDouble time)
 Update m_fields with u^n by multiplying by inverse mass matrix. That's then used in e.g. checkpoint output and L^2 error calculation. More...
 
SOLVER_UTILS_EXPORT void ALEDoElmtInvMassBwdTrans (const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
 
SOLVER_UTILS_EXPORT void MoveMesh (const NekDouble &time, Array< OneD, Array< OneD, NekDouble > > &traceNormals)
 
const Array< OneD, const Array< OneD, NekDouble > > & GetGridVelocity ()
 
SOLVER_UTILS_EXPORT const Array< OneD, const Array< OneD, NekDouble > > & GetGridVelocityTrace ()
 
SOLVER_UTILS_EXPORT void ExtraFldOutputGridVelocity (std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
 
- Public Member Functions inherited from Nektar::SolverUtils::FluidInterface
virtual ~FluidInterface ()=default
 
SOLVER_UTILS_EXPORT void GetVelocity (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &velocity)
 Extract array with velocity from physfield. More...
 
SOLVER_UTILS_EXPORT bool HasConstantDensity ()
 
SOLVER_UTILS_EXPORT void GetDensity (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &density)
 Extract array with density from physfield. More...
 
SOLVER_UTILS_EXPORT void GetPressure (const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &pressure)
 Extract array with pressure from physfield. More...
 
SOLVER_UTILS_EXPORT void SetMovingFrameVelocities (const Array< OneD, NekDouble > &vFrameVels, const int step)
 
SOLVER_UTILS_EXPORT bool GetMovingFrameVelocities (Array< OneD, NekDouble > &vFrameVels, const int step)
 
SOLVER_UTILS_EXPORT void SetMovingFrameDisp (const Array< OneD, NekDouble > &vFrameDisp, const int step)
 
SOLVER_UTILS_EXPORT void SetMovingFramePivot (const Array< OneD, NekDouble > &vFramePivot)
 
SOLVER_UTILS_EXPORT bool GetMovingFrameDisp (Array< OneD, NekDouble > &vFrameDisp, const int step)
 
SOLVER_UTILS_EXPORT void SetAeroForce (Array< OneD, NekDouble > forces)
 Set aerodynamic force and moment. More...
 
SOLVER_UTILS_EXPORT void GetAeroForce (Array< OneD, NekDouble > forces)
 Get aerodynamic force and moment. More...
 
- Protected Types inherited from Nektar::SolverUtils::EquationSystem
enum  HomogeneousType { eHomogeneous1D , eHomogeneous2D , eHomogeneous3D , eNotHomogeneous }
 Parameter for homogeneous expansions. More...
 

Detailed Description

Definition at line 45 of file SmoothedProfileMethod.h.

Constructor & Destructor Documentation

◆ SmoothedProfileMethod()

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

Construct a new Smoothed Profile Method object.

Parameters
pSession
pGraph

Definition at line 61 of file SmoothedProfileMethod.cpp.

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

◆ ~SmoothedProfileMethod()

Nektar::SmoothedProfileMethod::~SmoothedProfileMethod ( )
overrideprotecteddefault

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

54 {
57 pGraph);
58 p->InitObject();
59 return p;
60 }
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 563 of file SmoothedProfileMethod.cpp.

564{
565 // Get the handler of first function block
566 TiXmlElement *conds = m_session->GetElement("Nektar/Conditions");
567 TiXmlElement *function = conds->FirstChildElement("FUNCTION");
568
569 // Loop over functions until the block 'name' is found
570 std::string functionType = function->Attribute("NAME");
571 while (function && !boost::iequals(functionType, functionName))
572 {
573 function = function->NextSiblingElement("FUNCTION");
574 functionType = function->Attribute("NAME");
575 }
576
577 return function;
578}
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  elemName 
)
protected

True if the function is timedependent, false otherwise.

Parameters
name
type
attribute
Returns
string

Definition at line 523 of file SmoothedProfileMethod.cpp.

525{
526 // Get the handler of the function
527 TiXmlElement *function = GetFunctionHdl(funcName);
528
529 // Go to the first element
530 TiXmlElement *functionDef = function->FirstChildElement();
531 ASSERTL0(functionDef, "At least one element must be defined in " + funcName)
532
533 // And search the element with name 'elemName'
534 std::string varName = functionDef->Attribute("VAR");
535 while (functionDef && !boost::iequals(varName, elemName))
536 {
537 functionDef = functionDef->NextSiblingElement();
538 varName = functionDef->Attribute("VAR");
539 }
540
541 ASSERTL0(functionDef,
542 "Variable " + elemName + " must be defined in " + funcName + ".");
543
544 // And return the value of USERDEFINEDTYPE
545 std::string attr;
546 int err = functionDef->QueryStringAttribute("USERDEFINEDTYPE", &attr);
547 bool output = boost::iequals(attr, "TimeDependent");
548
549 ASSERTL0((err == TIXML_NO_ATTRIBUTE) || (err == TIXML_SUCCESS && output),
550 "USERDEFINEDTYPE in " + elemName +
551 " must be TimeDependent if defined");
552
553 return output;
554}
#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 580 of file SmoothedProfileMethod.cpp.

581{
582 // Function evaluator for Phi and Up
583 m_phiEvaluator = GetFunction("ShapeFunction");
584
585 TiXmlElement *function = GetFunctionHdl("ShapeFunction");
586 TiXmlElement *child = function->FirstChildElement();
587 m_filePhi = false;
588
589 // If defined by using a file
590 if (boost::iequals(child->ValueStr(), "F"))
591 {
592 // Get name of STL file
593 std::string fileName;
594 int status = child->QueryStringAttribute("FILE", &fileName);
595 ASSERTL0(status == TIXML_SUCCESS,
596 "An FLD file with the values "
597 "of the phi function has to be supplied.")
598 ASSERTL0(boost::iequals(fileName.substr(fileName.length() - 4), ".fld"),
599 "A valid FLD file must be supplied in the "
600 "'ShapeFunction' field.")
601
602 // Get phi values from XML file (after "FieldConvert" the STL file)
603 // First, load the data
604 std::vector<LibUtilities::FieldDefinitionsSharedPtr> fieldDef;
605 std::vector<std::vector<NekDouble>> fieldData;
606 LibUtilities::FieldMetaDataMap fieldMetaData;
607 LibUtilities::FieldIOSharedPtr phiFile =
608 LibUtilities::FieldIO::CreateForFile(m_session, fileName);
609 phiFile->Import(fileName, fieldDef, fieldData, fieldMetaData);
610
611 // Only Phi field should be defined in the file
612 ASSERTL0(fieldData.size() == 1, "Only one field (phi) must be "
613 "defined in the FLD file.")
614
615 // Extract Phi field to output
616 std::string tmp("phi");
617 m_phi->ExtractDataToCoeffs(fieldDef[0], fieldData[0], tmp,
618 m_phi->UpdateCoeffs());
619 m_phi->BwdTrans(m_phi->GetCoeffs(), m_phi->UpdatePhys());
620 m_filePhi = true;
621 m_timeDependentPhi = false;
622 m_timeDependentUp = false;
623 }
624 else
625 {
626 // Check if Phi is timedependent
627 m_timeDependentPhi = GetVarTimeDependence("ShapeFunction", "Phi");
628
629 // If so, check if its velocity changes as well
630 m_timeDependentUp = GetVarTimeDependence("ShapeFunction", "Up");
631 switch (m_velocity.size())
632 {
633 case 2:
635 GetVarTimeDependence("ShapeFunction", "Vp");
636 break;
637 case 3:
639 GetVarTimeDependence("ShapeFunction", "Vp");
641 GetVarTimeDependence("ShapeFunction", "Wp");
642 break;
643 }
644 }
645}
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:287
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 395 of file SmoothedProfileMethod.cpp.

396{
397 size_t nvel = m_velocity.size();
398 Array<OneD, ExpListSharedPtr> BndExp;
399 Array<OneD, SpatialDomains::BoundaryConditionShPtr> BndCond;
400
401 // Get the BC expansions
402 BndExp = m_pressureP->GetBndCondExpansions();
403 BndCond = m_pressureP->GetBndConditions();
404
405 // For each boundary...
406 for (size_t b = 0; b < BndExp.size(); ++b)
407 {
408 // Only for BCs based on the derivative
409 if (BndCond[b]->GetBoundaryConditionType() ==
411 BndCond[b]->GetBoundaryConditionType() == SpatialDomains::ePeriodic)
412 {
413 // Calculate f_s values
414 Array<OneD, Array<OneD, NekDouble>> f_s(nvel);
415 for (size_t i = 0; i < nvel; ++i)
416 {
417 f_s[i] = m_fs[0]->GetBndCondExpansions()[b]->GetPhys();
418 }
419
420 // BC is f_s * n
421 BndExp[b]->NormVectorIProductWRTBase(f_s, BndExp[b]->UpdatePhys());
422 }
423 }
424}
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 292 of file SmoothedProfileMethod.cpp.

295{
296 size_t physTot = m_fs[0]->GetNpoints();
297 size_t nvel = m_velocity.size();
298
299 // Set boundary conditions
301
302 // Divergence of 'fs'
303 m_fields[m_velocity[0]]->PhysDeriv(eX, m_fs[0]->GetPhys(), Forcing[0]);
304
305 // Using 'Forcing[1]' as storage
306 for (size_t i = 1; i < nvel; ++i)
307 {
308 size_t ind = m_velocity[i];
309 m_fields[ind]->PhysDeriv(DirCartesianMap[i], m_fs[i]->GetPhys(),
310 Forcing[1]);
311 Vmath::Vadd(physTot, Forcing[1], 1, Forcing[0], 1, Forcing[0], 1);
312 }
313}
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 133 of file SmoothedProfileMethod.h.

134 {
135 int iVel = m_velocity[0];
137 *std::dynamic_pointer_cast<T>(m_pressure));
139 *std::dynamic_pointer_cast<T>(m_fields[iVel]));
140
141 m_fs = Array<OneD, MultiRegions::ExpListSharedPtr>(nvel);
142 for (int i = 0; i < nvel; ++i)
143 {
145 *std::dynamic_pointer_cast<T>(m_fields[iVel]));
146 }
147
148 // Set to wave space if homogeneous case
150 {
151 m_pressureP->SetWaveSpace(true);
152 m_phi->SetWaveSpace(true);
153 for (int i = 0; i < nvel; ++i)
154 {
155 m_fs[i]->SetWaveSpace(true);
156 }
157 }
158 }
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 348 of file SmoothedProfileMethod.cpp.

351{
352 size_t physTot = m_phi->GetNpoints();
353
354 // Gradient of p_p
355 size_t nvel = m_velocity.size();
356 if (nvel == 2)
357 {
358 m_pressureP->PhysDeriv(m_pressureP->GetPhys(), Forcing[0], Forcing[1]);
359 }
360 else
361 {
362 m_pressureP->PhysDeriv(m_pressureP->GetPhys(), Forcing[0], Forcing[1],
363 Forcing[2]);
364 }
365
366 // Velocity correction
367 for (size_t i = 0; i < nvel; ++i)
368 {
369 // Adding -(1-m_phi)*grad(p_p) instead of -grad(p_p) reduces the
370 // flux through the walls, but the flow is not incompressible
371 if (m_session->DefinesSolverInfo("ForceBoundary") &&
372 boost::iequals(m_session->GetSolverInfo("ForceBoundary"), "True"))
373 {
374 Vmath::Vvtvm(physTot, m_phi->GetPhys(), 1, Forcing[i], 1,
375 Forcing[i], 1, Forcing[i], 1);
376 Vmath::Vadd(physTot, m_fs[i]->GetPhys(), 1, Forcing[i], 1,
377 Forcing[i], 1);
378 }
379 else
380 {
381 Vmath::Vsub(physTot, m_fs[i]->GetPhys(), 1, Forcing[i], 1,
382 Forcing[i], 1);
383 }
384 Blas::Daxpy(physTot, dt / m_gamma0, Forcing[i], 1, fields[i], 1);
385 }
386}
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 323 of file SmoothedProfileMethod.cpp.

325{
327 // Factor 'lambda=0' in Helmholtz equation to get the Poisson form
329
330 // Solve the Poisson equation
331 m_pressureP->HelmSolve(Forcing, m_pressureP->UpdateCoeffs(), factors);
332
333 // Update node values from coefficients
334 m_pressureP->BwdTrans(m_pressureP->GetCoeffs(), m_pressureP->UpdatePhys());
335}
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 475 of file SmoothedProfileMethod.cpp.

477{
478 size_t nvel = m_velocity.size();
479 size_t nq = m_phi->GetNpoints();
480
481 for (size_t i = 0; i < nvel; ++i)
482 {
483 // In homogeneous cases, switch out of wave space
484 Array<OneD, NekDouble> tmpField(nq);
485 size_t ind = m_velocity[i];
486
488 m_fields[ind]->GetWaveSpace())
489 {
490 m_fields[ind]->HomogeneousBwdTrans(nq, fields[i], tmpField);
491 m_fs[i]->HomogeneousBwdTrans(nq, m_fs[i]->GetPhys(),
492 m_fs[i]->UpdatePhys());
493 }
494 else
495 {
496 tmpField = fields[i];
497 }
498
499 Vmath::Vsub(nq, m_up[i], 1, tmpField, 1, m_fs[i]->UpdatePhys(), 1);
500 Vmath::Vmul(nq, m_phi->GetPhys(), 1, m_fs[i]->GetPhys(), 1,
501 m_fs[i]->UpdatePhys(), 1);
502 Vmath::Smul(nq, m_gamma0 / dt, m_fs[i]->GetPhys(), 1,
503 m_fs[i]->UpdatePhys(), 1);
504
505 // And go back to wave space if the 'if' was executed
507 m_fields[ind]->GetWaveSpace())
508 {
509 m_fs[i]->HomogeneousFwdTrans(nq, m_fs[i]->GetPhys(),
510 m_fs[i]->UpdatePhys());
511 }
512 }
513}
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 431 of file SmoothedProfileMethod.cpp.

432{
433 // Initialise 'm_up' and 'm_phi' during first step
434 if (t <= 0.0)
435 {
436 if (!m_filePhi)
437 {
438 // Update 'm_phi' only if it was provided as a function
439 m_phiEvaluator->Evaluate("Phi", m_phi->UpdatePhys(), t);
440 }
441
442 // Initialize both variables for the first step
443 m_phiEvaluator->Evaluate(m_velName, m_up, t);
444
445 // Initialise 'm_upPrev' in all cases
446 m_upPrev = m_up;
447 }
448 // If timedependent 'm_phi'
449 // Phi functions from files are not timedependent
450 else if (m_timeDependentPhi)
451 {
452 m_phiEvaluator->Evaluate("Phi", m_phi->UpdatePhys(), t);
453
454 // And if velocities are timedependent as well
456 {
457 // Store previous value of u_p during simulation
458 m_upPrev = m_up;
459 m_phiEvaluator->Evaluate(m_velName, m_up, t);
460 }
461 }
462}
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)
overrideprotectedvirtual

Generates the summary of the current simulation.

Parameters
s

Reimplemented from Nektar::SolverUtils::UnsteadySystem.

Definition at line 224 of file SmoothedProfileMethod.cpp.

225{
227 SolverUtils::AddSummaryItem(s, "IB formulation",
228 "Smoothed Profile Method (SPM)");
229}
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)
overrideprotectedvirtual

Initialisation object for EquationSystem.

Continuous field

Setting up the normals

Setting up the normals

Reimplemented from Nektar::IncNavierStokes.

Definition at line 69 of file SmoothedProfileMethod.cpp.

70{
72
73 // Update implicit time-intregration class operators
76
77 // Number of dims as number of velocity vectors
78 size_t nvel = m_velocity.size();
79
80 // Initialization of correction pressure and shape function
81 switch (nvel)
82 {
83 case 1:
85 {
86 SetUpExpansions<ContField>(nvel);
87 }
89 {
90 SetUpExpansions<DisContField>(nvel);
91 }
92 break;
93
94 case 2:
96 {
97 SetUpExpansions<ContField>(nvel);
98 }
100 {
101 SetUpExpansions<DisContField>(nvel);
102 }
103 break;
104
105 case 3:
107 {
109 {
110 SetUpExpansions<ContField>(nvel);
111 }
113 {
114 SetUpExpansions<ContField3DHomogeneous1D>(nvel);
115 }
118 {
119 SetUpExpansions<ContField3DHomogeneous2D>(nvel);
120 }
121 }
123 {
125 {
126 SetUpExpansions<DisContField>(nvel);
127 }
129 {
130 SetUpExpansions<DisContField3DHomogeneous1D>(nvel);
131 }
134 {
135 SetUpExpansions<DisContField3DHomogeneous2D>(nvel);
136 }
137 }
138 break;
139 }
140
141 // Read 'm_phi' and its velocity
142 ASSERTL0(m_session->DefinesFunction("ShapeFunction"),
143 "ShapeFunction must be defined in the session file.")
144 ReadPhi();
145
146 // Allocate the vector 'm_up'
147 size_t physTot = m_phi->GetTotPoints();
148 m_velName.push_back("Up");
149 if (nvel > 1)
150 {
151 m_velName.push_back("Vp");
152 }
153 if (nvel == 3)
154 {
155 m_velName.push_back("Wp");
156 }
157
158 m_up = Array<OneD, Array<OneD, NekDouble>>(nvel);
159 for (size_t i = 0; i < nvel; ++i)
160 {
161 m_up[i] = Array<OneD, NekDouble>(physTot, 0.0);
162 }
163
164 // Make sure that m_phi and m_up are defined
165 UpdatePhiUp(0.0);
166
167 // Get the time integration scheme.
168 LibUtilities::TimeIntScheme timeInt;
169 if (m_session->DefinesTimeIntScheme())
170 {
171 timeInt = m_session->GetTimeIntScheme();
172 }
173 else
174 {
175 timeInt.method = m_session->GetSolverInfo("TimeIntegrationMethod");
176 timeInt.order = timeInt.method.back() - '0';
177
178 // Remove everything past the IMEX.
179 timeInt.method = timeInt.method.substr(0, 4);
180 }
181
182 // Select 'm_gamma0' depending on IMEX order
183 ASSERTL0(
184 boost::iequals(timeInt.method, "IMEX") && 1 <= timeInt.order &&
185 timeInt.order <= 4,
186 "The TimeIntegrationMethod scheme must be IMEX with order '1' to '4'.")
187
188 switch (timeInt.order)
189 {
190 case 1:
191 m_gamma0 = 1.0;
192 break;
193
194 case 2:
195 m_gamma0 = 3.0 / 2.0;
196 break;
197
198 case 3:
199 m_gamma0 = 11.0 / 6.0;
200 break;
201
202 case 4:
203 m_gamma0 = 25.0 / 12.0;
204 break;
205 }
206
207 // Check if the aeroforces filter is active, negative if inactive
208 m_forcesFilter = -1;
209 for (size_t i = 0; i < m_session->GetFilters().size(); ++i)
210 {
211 if (boost::iequals(m_session->GetFilters()[i].name, "AeroForcesSPM"))
212 {
213 m_forcesFilter = i;
214 break;
215 }
216 }
217}
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.
SOLVER_UTILS_EXPORT int GetTotPoints()
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
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 245 of file SmoothedProfileMethod.cpp.

249{
251 time, a_iixDt);
252
253 size_t physTot = m_pressureP->GetNpoints();
254
255 /* SPM correction of velocity */
256 // Update 'm_phi' and 'm_up' if needed (evaluated at next time step)
257 UpdatePhiUp(time + a_iixDt);
258 // Update calculation of IB forcing 'm_fs'
259 UpdateForcing(outarray, a_iixDt);
260 // Estimate forces only if requested
261 if (m_forcesFilter >= 0)
262 {
263 std::static_pointer_cast<FilterAeroForcesSPM>(
265 ->CalculateForces(outarray, m_upPrev, m_phi, time, a_iixDt);
266 }
267 // Set BC conditions for pressure p_p
268 SetUpCorrectionPressure(outarray, m_F);
269 // Solve Poisson equation for pressure p_p
271 // Solve velocity in the next step with IB
272 SolveCorrectedVelocity(m_F, outarray, a_iixDt);
273
274 // Add pressures to get final value
275 Vmath::Vadd(physTot, m_pressure->GetPhys(), 1, m_pressureP->GetPhys(), 1,
276 m_pressure->UpdatePhys(), 1);
277 m_pressure->FwdTrans(m_pressure->GetPhys(), m_pressure->UpdateCoeffs());
278
279 // Add presure to outflow bc if using convective like BCs
280 m_extrapolation->AddPressureToOutflowBCs(m_kinvis);
281}
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().

Friends And Related Function Documentation

◆ MemoryManager< SmoothedProfileMethod >

friend class MemoryManager< SmoothedProfileMethod >
friend

Definition at line 163 of file SmoothedProfileMethod.h.

Member Data Documentation

◆ className

std::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 63 of file SmoothedProfileMethod.h.

◆ m_filePhi

bool Nektar::SmoothedProfileMethod::m_filePhi
protected

Flag indicating that phi was defined in a file.

Definition at line 86 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 88 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 76 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 80 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 82 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 84 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 74 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 69 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 70 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 72 of file SmoothedProfileMethod.h.

Referenced by UpdatePhiUp(), and v_InitObject().

◆ solverTypeLookupId

std::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 90 of file SmoothedProfileMethod.h.