Nektar++
FilterAeroForcesSPM.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: FilterAeroForcesSPM.cpp
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: Output values of aerodynamic forces during time-stepping.
33//
34///////////////////////////////////////////////////////////////////////////////
35
37
38namespace Nektar
39{
40
43 "AeroForcesSPM", FilterAeroForcesSPM::create);
44
45/**
46 *
47 */
50 const std::shared_ptr<SolverUtils::EquationSystem> &pEquation,
51 const ParamMap &pParams)
52 : Filter(pSession, pEquation)
53{
54 // OutputFile
55 std::string ext = ".fce";
56 m_outputFile = Filter::SetupOutput(ext, pParams);
57
58 // OutputFrequency
59 auto it = pParams.find("OutputFrequency");
60 if (it == pParams.end())
61 {
63 }
64 else
65 {
66 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
67 m_outputFrequency = round(equ.Evaluate());
68 }
69
70 // Time after which we need to calculate the forces
71 it = pParams.find("StartTime");
72 if (it == pParams.end())
73 {
74 m_startTime = 0;
75 }
76 else
77 {
78 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
79 m_startTime = equ.Evaluate();
80 }
81}
82
83/**
84 *
85 */
88 [[maybe_unused]] const NekDouble &time)
89{
90 // Save space dimension
91 m_spaceDim = pFields[0]->GetGraph()->GetMeshDimension();
92
93 // Fill in the directions vector
94 m_dirNames.push_back("X");
95 if (m_spaceDim > 1)
96 {
97 m_dirNames.push_back("Y");
98 }
99 if (m_spaceDim > 2)
100 {
101 m_dirNames.push_back("Z");
102 }
103
104 // Allocate the aero-forces vector
106
107 // Write header
108 LibUtilities::CommSharedPtr vComm = pFields[0]->GetComm();
109 if (vComm->GetRank() == 0)
110 {
111 // Open output stream
112 bool adaptive;
113 m_session->MatchSolverInfo("Driver", "Adaptive", adaptive, false);
114 if (adaptive)
115 {
116 m_outputStream.open(m_outputFile.c_str(), std::ofstream::app);
117 }
118 else
119 {
120 m_outputStream.open(m_outputFile.c_str());
121 }
122 m_outputStream << "# Forces acting on bodies\n";
123 m_outputStream << "#";
124 m_outputStream.width(7);
125 m_outputStream << "Time";
126 for (int i = 0; i < m_spaceDim; ++i)
127 {
128 m_outputStream.width(14);
129 m_outputStream << "F_" << m_dirNames[i];
130 }
131
132 m_outputStream << std::endl;
133 }
134
135 m_index = 0;
136}
137
138/**
139 *
140 */
143 const NekDouble &time)
144{
145 // Only output every m_outputFrequency.
146 if ((m_index++) % m_outputFrequency || (time < m_startTime))
147 {
148 return;
149 }
150
151 // Communicators to exchange results
152 LibUtilities::CommSharedPtr vComm = pFields[0]->GetComm();
153
154 // Write Results
155 if (vComm->GetRank() == 0)
156 {
157 // Write time
158 m_outputStream.width(8);
159 m_outputStream << std::setprecision(6) << time;
160 // Write forces
161 for (int i = 0; i < m_spaceDim; ++i)
162 {
163 m_outputStream.width(15);
164 m_outputStream << std::setprecision(8) << m_Forces[i];
165 }
166 m_outputStream.width(10);
167 m_outputStream << std::endl;
168 }
169}
170
171/**
172 *
173 */
176 [[maybe_unused]] const NekDouble &time)
177{
178 if (pFields[0]->GetComm()->GetRank() == 0)
179 {
180 m_outputStream.close();
181 }
182}
183
184/**
185 *
186 */
188{
189 return true;
190}
191
192/**
193 * @brief Determine the total force on the body defined by \f$\Phi\f$
194 * (note that if the shape function represents more than one
195 * body, this function calculates the value of the final force after adding
196 * up the values for each body). This value must be scaled with the
197 * density to get the real force vector.
198 *
199 * @param pIntVel
200 * @param pUpPrev
201 * @param pPhi
202 * @param time
203 * @param dt
204 */
206 const Array<OneD, Array<OneD, NekDouble>> &pIntVel,
207 const Array<OneD, Array<OneD, NekDouble>> &pUpPrev,
209{
210 // Only output every m_outputFrequency.
211 if ((m_index % m_outputFrequency) || (time < m_startTime))
212 {
213 return;
214 }
215
216 int nq = pIntVel[0].size();
218
219 // Aerodynamic forces are computed according eq. 18a in Luo et al. (2009).
220 // Smoothed profile method for particulate flows: Error analysis and
221 // simulations. Journal of Computational Physics, 228(5)
222 for (int i = 0; i < m_spaceDim; ++i)
223 {
224 // "Scalar" field to be integrated
225 Vmath::Vsub(nq, pIntVel[i], 1, pUpPrev[i], 1, tmp, 1);
226 Vmath::Vmul(nq, pPhi->GetPhys(), 1, tmp, 1, tmp, 1);
227
228 // Integration of force throughout the domain
229 m_Forces[i] = pPhi->Integral(tmp) / dt;
230 }
231}
232
233} // namespace Nektar
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::shared_ptr< SolverUtils::EquationSystem > &pEquation, const std::map< std::string, std::string > &pParams)
static SolverUtils::FilterSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const std::shared_ptr< SolverUtils::EquationSystem > &pEquation, const std::map< std::string, std::string > &pParams)
Creates an instance of this class.
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.
void v_Finalise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
SOLVER_UTILS_EXPORT std::string SetupOutput(const std::string ext, const ParamMap &pParams)
Definition: Filter.h:139
LibUtilities::SessionReaderSharedPtr m_session
Definition: Filter.h:93
std::map< std::string, std::string > ParamMap
Definition: Filter.h:66
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
FilterFactory & GetFilterFactory()
double NekDouble
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 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