Nektar++
FilterAeroForcesSPM.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: FilterAeroForcesSPM.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: Outputs values at specific points during time-stepping.
33//
34///////////////////////////////////////////////////////////////////////////////
35
36#ifndef NEKTAR_INCNAVIERSTOKES_FILTERS_FILTERAEROFORCESSPM_H
37#define NEKTAR_INCNAVIERSTOKES_FILTERS_FILTERAEROFORCESSPM_H
38
40
41namespace Nektar
42{
44{
45public:
46 /// Creates an instance of this class
49 const std::weak_ptr<SolverUtils::EquationSystem> &pEquation,
50 const std::map<std::string, std::string> &pParams)
51 {
54 pSession, pEquation, pParams);
55 return p;
56 }
57
58 /// Name of the class
59 static std::string className;
60
63 const std::weak_ptr<SolverUtils::EquationSystem> &pEquation,
64 const std::map<std::string, std::string> &pParams);
65
66 ~FilterAeroForcesSPM() override;
67
68 // Calculates the forces and fills the array 'm_Forces' up
70 const Array<OneD, Array<OneD, NekDouble>> &pUpPrev,
72 NekDouble time, NekDouble dt);
73
74protected:
75 unsigned int m_index;
76 unsigned int m_outputFrequency;
77 std::string m_outputFile;
78 std::ofstream m_outputStream;
79 // Time when we start calculating the forces
81 /// STL vector containing the names of the different directions
82 std::vector<std::string> m_dirNames;
83 /// Dimension of the fluid domain
85 /// Array storing the last value of the aerodynamic forces
87
88 void v_Initialise(
90 const NekDouble &time) override;
91
92 void v_Update(
94 const NekDouble &time) override;
95
96 void v_Finalise(
98 const NekDouble &time) override;
99
100 bool v_IsTimeDependent() override;
101
102private:
103};
104
105typedef std::shared_ptr<FilterAeroForcesSPM> FilterAeroForcesSPMSharedPtr;
106} // namespace Nektar
107
108#endif /* NEKTAR_INCNAVIERSTOKES_FILTERS_FILTERAEROFORCESSPM_H */
void v_Initialise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
NekDouble m_spaceDim
Dimension of the fluid domain.
void v_Update(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
void CalculateForces(const Array< OneD, Array< OneD, NekDouble > > &pIntVel, const Array< OneD, Array< OneD, NekDouble > > &pUpPrev, const MultiRegions::ExpListSharedPtr &pPhi, NekDouble time, NekDouble dt)
Determine the total force on the body defined by (note that if the shape function represents more th...
static std::string className
Name of the class.
FilterAeroForcesSPM(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< SolverUtils::EquationSystem > &pEquation, const std::map< std::string, std::string > &pParams)
std::vector< std::string > m_dirNames
STL vector containing the names of the different directions.
Array< OneD, NekDouble > m_Forces
Array storing the last value of the aerodynamic forces.
static SolverUtils::FilterSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< SolverUtils::EquationSystem > &pEquation, const std::map< std::string, std::string > &pParams)
Creates an instance of this class.
void v_Finalise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< Filter > FilterSharedPtr
A shared pointer to a Driver object.
Definition: Filter.h:51
std::shared_ptr< FilterAeroForcesSPM > FilterAeroForcesSPMSharedPtr
double NekDouble