Nektar++
SmoothedProfileMethod.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: SmoothedProfileMethod.h
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// License for the specific language governing rights and limitations under
14// Permission is hereby granted, free of charge, to any person obtaining a
15// copy of this software and associated documentation files (the "Software"),
16// to deal in the Software without restriction, including without limitation
17// the rights to use, copy, modify, merge, publish, distribute, sublicense,
18// and/or sell copies of the Software, and to permit persons to whom the
19// Software is furnished to do so, subject to the following conditions:
20//
21// The above copyright notice and this permission notice shall be included
22// in all copies or substantial portions of the Software.
23//
24// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30// DEALINGS IN THE SOFTWARE.
31//
32// Description: Smoothed Profile Method header
33//
34///////////////////////////////////////////////////////////////////////////////
35
36#ifndef NEKTAR_SOLVERS_SMOOTHEDPROFILEMETHOD_H
37#define NEKTAR_SOLVERS_SMOOTHEDPROFILEMETHOD_H
38
40#include <fstream>
41
42namespace Nektar
43{
45{
46public:
47 /// Creates an instance of this class
51 {
54 pGraph);
55 p->InitObject();
56 return p;
57 }
58
59 /// Name of class
60 static std::string className;
61
62 // Constructor
65
66 // Destructor
67 ~SmoothedProfileMethod() override;
68
69 void v_InitObject(bool DeclareField = true) override;
70
72
73protected:
74 /// Correction pressure field for SPM
76 /// Velocity of the immersed body(ies)
79 /// Vector storing the names of the components of \u_p
80 std::vector<std::string> m_velName;
81 /// Flag signaling if \f$\u_p\f$ depends on time
83 /// Stiffly-stable scheme \f$\gamma_0\f$ coefficient
85 /// Forcing function 'f_s'
87 /// Shape function 'phi' as expansion list
89 /// Function that evaluates the values of \Phi
91 /// Flag that is true when phi depends on time
93 /// Flag indicating that phi was defined in a file
95 /// Position of "AeroForcesSPM" filter in 'm_session->GetFilters()'
97
98 static std::string solverTypeLookupId;
99
100 // Interface for 'v_SolveUnsteadyStokesSystem'
102 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
103 Array<OneD, Array<OneD, NekDouble>> &outarray, NekDouble time,
104 NekDouble a_iixDt) override;
105 // Sets the parameters and BCs for the Poisson equation
107 const Array<OneD, const Array<OneD, NekDouble>> &fields,
109 // Solves the Poisson equation for the correction pressure
111 // Explicitly corrects the velocity by using the force 'fs'
114 NekDouble dt);
115 // Set proper BCs for the corrected pressure 'p_p'
117 // Calculates the shape function values
118 // (only for non-moving boundaries)
119 void UpdatePhiUp(NekDouble time);
120 // Calculates the virtual force 'fs'
121 void UpdateForcing(const Array<OneD, const Array<OneD, NekDouble>> &fields,
122 NekDouble dt);
123 // Gets time-dependence information from the first elment of 'name'
124 bool GetVarTimeDependence(std::string funcName, std::string attrName);
125 // Returns a Tinyxml handle of the requested function element
126 TiXmlElement *GetFunctionHdl(std::string functionName);
127 // Reads and set the values of Phi from an analyitic func. or a file
128 void ReadPhi();
129
130 /// Initialises the expansions for the intermediate pressure, the 'Phi'
131 /// function and the IB forcing
132 template <typename T> void SetUpExpansions(int nvel)
133 {
134 int iVel = m_velocity[0];
136 *std::dynamic_pointer_cast<T>(m_pressure));
138 *std::dynamic_pointer_cast<T>(m_fields[iVel]));
139
141 for (int i = 0; i < nvel; ++i)
142 {
144 *std::dynamic_pointer_cast<T>(m_fields[iVel]));
145 }
146
147 // Set to wave space if homogeneous case
149 {
150 m_pressureP->SetWaveSpace(true);
151 m_phi->SetWaveSpace(true);
152 for (int i = 0; i < nvel; ++i)
153 {
154 m_fs[i]->SetWaveSpace(true);
155 }
156 }
157 }
158
159private:
160};
161
162typedef std::shared_ptr<SmoothedProfileMethod> SmoothedProfileMethodSharedPtr;
163
164} // namespace Nektar
165
166#endif // SMOOTHEDPROFILEMETHOD_H
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 v_GenerateSummary(SolverUtils::SummaryList &s) override
Generates the summary of the current simulation.
void SetCorrectionPressureBCs()
Updates the BCs for boundaries with Neumann or periodic BCs in the pressure:
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...
Array< OneD, Array< OneD, NekDouble > > m_up
Velocity of the immersed body(ies)
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()'.
void v_InitObject(bool DeclareField=true) override
Initialisation object for EquationSystem.
bool m_timeDependentPhi
Flag that is true when phi depends on time.
~SmoothedProfileMethod() override
Destroy the Smoothed Profile Method object.
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.
static std::string className
Name of class.
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
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.
Definition: Forcing.h:71
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
Definition: Misc.h:46
std::shared_ptr< EquationSystem > EquationSystemSharedPtr
A shared pointer to an EquationSystem object.
std::shared_ptr< SessionFunction > SessionFunctionSharedPtr
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::shared_ptr< SmoothedProfileMethod > SmoothedProfileMethodSharedPtr
double NekDouble