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
38using namespace std;
39
40namespace Nektar
41{
44 "AeroForcesSPM", FilterAeroForcesSPM::create);
45
46/**
47 *
48 */
51 const std::weak_ptr<SolverUtils::EquationSystem> &pEquation,
52 const ParamMap &pParams)
53 : Filter(pSession, pEquation)
54{
55 // OutputFile
56 auto it = pParams.find("OutputFile");
57 if (it == pParams.end())
58 {
59 m_outputFile = m_session->GetSessionName();
60 }
61 else
62 {
63 ASSERTL0(it->second.length() > 0, "Missing parameter 'OutputFile'.");
64 m_outputFile = it->second;
65 }
66 if (!(m_outputFile.length() >= 4 &&
67 m_outputFile.substr(m_outputFile.length() - 4) == ".fce"))
68 {
69 m_outputFile += ".fce";
70 }
71
72 // OutputFrequency
73 it = pParams.find("OutputFrequency");
74 if (it == pParams.end())
75 {
77 }
78 else
79 {
80 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
81 m_outputFrequency = round(equ.Evaluate());
82 }
83
84 // Time after which we need to calculate the forces
85 it = pParams.find("StartTime");
86 if (it == pParams.end())
87 {
88 m_startTime = 0;
89 }
90 else
91 {
92 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
93 m_startTime = equ.Evaluate();
94 }
95}
96
97/**
98 *
99 */
101{
102}
103
104/**
105 *
106 */
109 [[maybe_unused]] const NekDouble &time)
110{
111 // Save space dimension
112 m_spaceDim = pFields[0]->GetGraph()->GetMeshDimension();
113
114 // Fill in the directions vector
115 m_dirNames.push_back("X");
116 if (m_spaceDim > 1)
117 {
118 m_dirNames.push_back("Y");
119 }
120 if (m_spaceDim > 2)
121 {
122 m_dirNames.push_back("Z");
123 }
124
125 // Allocate the aero-forces vector
127
128 // Write header
129 LibUtilities::CommSharedPtr vComm = pFields[0]->GetComm();
130 if (vComm->GetRank() == 0)
131 {
132 // Open output stream
133 bool adaptive;
134 m_session->MatchSolverInfo("Driver", "Adaptive", adaptive, false);
135 if (adaptive)
136 {
137 m_outputStream.open(m_outputFile.c_str(), ofstream::app);
138 }
139 else
140 {
141 m_outputStream.open(m_outputFile.c_str());
142 }
143 m_outputStream << "# Forces acting on bodies\n";
144 m_outputStream << "#";
145 m_outputStream.width(7);
146 m_outputStream << "Time";
147 for (int i = 0; i < m_spaceDim; ++i)
148 {
149 m_outputStream.width(14);
150 m_outputStream << "F_" << m_dirNames[i];
151 }
152
153 m_outputStream << endl;
154 }
155
156 m_index = 0;
157}
158
159/**
160 *
161 */
164 const NekDouble &time)
165{
166 // Only output every m_outputFrequency.
167 if ((m_index++) % m_outputFrequency || (time < m_startTime))
168 {
169 return;
170 }
171
172 // Communicators to exchange results
173 LibUtilities::CommSharedPtr vComm = pFields[0]->GetComm();
174
175 // Write Results
176 if (vComm->GetRank() == 0)
177 {
178 // Write time
179 m_outputStream.width(8);
180 m_outputStream << setprecision(6) << time;
181 // Write forces
182 for (int i = 0; i < m_spaceDim; ++i)
183 {
184 m_outputStream.width(15);
185 m_outputStream << setprecision(8) << m_Forces[i];
186 }
187 m_outputStream.width(10);
188 m_outputStream << endl;
189 }
190}
191
192/**
193 *
194 */
197 [[maybe_unused]] const NekDouble &time)
198{
199 if (pFields[0]->GetComm()->GetRank() == 0)
200 {
201 m_outputStream.close();
202 }
203}
204
205/**
206 *
207 */
209{
210 return true;
211}
212
213/**
214 * @brief Determine the total force on the body defined by \f$\Phi\f$
215 * (note that if the shape function represents more than one
216 * body, this function calculates the value of the final force after adding
217 * up the values for each body). This value must be scaled with the
218 * density to get the real force vector.
219 *
220 * @param pIntVel
221 * @param pUpPrev
222 * @param pPhi
223 * @param time
224 * @param dt
225 */
227 const Array<OneD, Array<OneD, NekDouble>> &pIntVel,
228 const Array<OneD, Array<OneD, NekDouble>> &pUpPrev,
230{
231 // Only output every m_outputFrequency.
232 if ((m_index % m_outputFrequency) || (time < m_startTime))
233 {
234 return;
235 }
236
237 int nq = pIntVel[0].size();
239
240 // Aerodynamic forces are computed according eq. 18a in Luo et al. (2009).
241 // Smoothed profile method for particulate flows: Error analysis and
242 // simulations. Journal of Computational Physics, 228(5)
243 for (int i = 0; i < m_spaceDim; ++i)
244 {
245 // "Scalar" field to be integrated
246 Vmath::Vsub(nq, pIntVel[i], 1, pUpPrev[i], 1, tmp, 1);
247 Vmath::Vmul(nq, pPhi->GetPhys(), 1, tmp, 1, tmp, 1);
248
249 // Integration of force throughout the domain
250 m_Forces[i] = pPhi->Integral(tmp) / dt;
251 }
252}
253
254} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
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
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
LibUtilities::SessionReaderSharedPtr m_session
Definition: Filter.h:83
std::map< std::string, std::string > ParamMap
Definition: Filter.h:65
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()
Definition: Filter.cpp:39
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