36 #ifndef NEKTAR_SOLVERS_SMOOTHEDPROFILEMETHOD_H
37 #define NEKTAR_SOLVERS_SMOOTHEDPROFILEMETHOD_H
139 std::string attrName);
147 template <
typename T>
154 *std::dynamic_pointer_cast<T>(
m_fields[iVel]));
157 for (
int i = 0; i < nvel; ++i)
160 *std::dynamic_pointer_cast<T>(
m_fields[iVel]));
167 m_phi->SetWaveSpace(
true);
168 for (
int i = 0; i < nvel; ++i)
170 m_fs[i]->SetWaveSpace(
true);
178 typedef std::shared_ptr<SmoothedProfileMethod>
MultiRegions::ExpListSharedPtr m_pressure
Pointer to field holding pressure field.
Array< OneD, int > m_velocity
int which identifies which components of m_fields contains the velocity (u,v,w);
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
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:
bool m_filePhi
Flag indicating that phi was defined in a file.
void SetCorrectionPressureBCs()
Updates the BCs for boundaries with Neumann or periodic BCs in the pressure:
Array< OneD, Array< OneD, NekDouble > > m_up
Velocity of the immersed body(ies)
virtual void v_InitObject()
Init object for UnsteadySystem class.
std::vector< std::string > m_velName
Vector storing the names of the components of \u_p.
NekDouble m_gamma0
Stiffly-stable scheme coefficient.
int m_forcesFilter
Position of "AeroForcesSPM" filter in 'm_session->GetFilters()'.
virtual ~SmoothedProfileMethod()
Destroy the Smoothed Profile Method object.
bool m_timeDependentPhi
Flag that is true when phi depends on time.
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
SmoothedProfileMethod(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Construct a new Smoothed Profile Method object.
void UpdatePhiUp(NekDouble time)
Calculates the values of the shape function.
SolverUtils::SessionFunctionSharedPtr m_phiEvaluator
Function that evaluates the values of \Phi.
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:
TiXmlElement * GetFunctionHdl(std::string functionName)
Returns a handle to the requested function. Returns NULL if it does not exist.
virtual void v_GenerateSummary(SolverUtils::SummaryList &s)
Generates the summary of the current simulation.
static std::string className
Name of class.
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...
bool m_timeDependentUp
Flag signaling if depends on time.
MultiRegions::ExpListSharedPtr m_pressureP
Correction pressure field for SPM.
void SetUpExpansions(int nvel)
Initialises the expansions for the intermediate pressure, the 'Phi' function and the IB forcing.
bool GetVarTimeDependence(std::string funcName, std::string attrName)
True if the function is timedependent, false otherwise.
Array< OneD, Array< OneD, NekDouble > > m_upPrev
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, MultiRegions::ExpListSharedPtr > m_fs
Forcing function 'f_s'.
MultiRegions::ExpListSharedPtr m_phi
Shape function 'phi' as expansion list.
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 :
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
enum HomogeneousType m_HomogeneousType
Defines a forcing term to be explicitly applied.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::vector< std::pair< std::string, std::string > > SummaryList
std::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.
std::shared_ptr< SessionFunction > SessionFunctionSharedPtr
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
The above copyright notice and this permission notice shall be included.
std::shared_ptr< SmoothedProfileMethod > SmoothedProfileMethodSharedPtr