Nektar++
FilterReynoldsStresses.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: FilterReynoldsStresses.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// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: Append Reynolds stresses to the average fields
32//
33///////////////////////////////////////////////////////////////////////////////
34
36
37namespace Nektar::SolverUtils
38{
40 GetFilterFactory().RegisterCreatorFunction("ReynoldsStresses",
42
43/**
44 * @class FilterReynoldsStresses
45 *
46 * @brief Append Reynolds stresses to the average fields
47 *
48 * This class appends the average fields with the Reynolds stresses of the form
49 * \f$ \overline{u' v'} \f$.
50 *
51 * For the default case, this is achieved by calculating
52 * \f$ C_{n} = \Sigma_{i=1}^{n} (u_i - \bar{u}_n)(v_i - \bar{v}_n)\f$
53 * using the recursive relation:
54 *
55 * \f[ C_{n} = C_{n-1} + \frac{n}{n-1} (u_n - \bar{u}_n)(v_n - \bar{v}_n) \f]
56 *
57 * The FilterSampler base class then divides the result by n, leading
58 * to the Reynolds stress.
59 *
60 * It is also possible to perform the averages using an exponential moving
61 * average, in which case either the moving average parameter \f$ \alpha \f$
62 * or the time constant \f$ \tau \f$ must be prescribed.
63 */
66 const std::shared_ptr<SolverUtils::EquationSystem> &pEquation,
67 const std::map<std::string, std::string> &pParams)
68 : FilterFieldConvert(pSession, pEquation, pParams)
69{
70 // Load sampling frequency
71 auto it = pParams.find("SampleFrequency");
72 if (it == pParams.end())
73 {
75 }
76 else
77 {
78 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
79 m_sampleFrequency = round(equ.Evaluate());
80 }
81
82 // Check if should use moving average
83 it = pParams.find("MovingAverage");
84 if (it == pParams.end())
85 {
86 m_movAvg = false;
87 }
88 else
89 {
90 std::string sOption = it->second.c_str();
91 m_movAvg = (boost::iequals(sOption, "true")) ||
92 (boost::iequals(sOption, "yes"));
93 }
94
95 if (!m_movAvg)
96 {
97 return;
98 }
99
100 // Load alpha parameter for moving average
101 it = pParams.find("alpha");
102 if (it == pParams.end())
103 {
104 it = pParams.find("tau");
105 if (it == pParams.end())
106 {
107 ASSERTL0(false, "MovingAverage needs either alpha or tau.");
108 }
109 else
110 {
111 // Load time constant
112 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
113 NekDouble tau = equ.Evaluate();
114 // Load delta T between samples
115 NekDouble dT;
116 m_session->LoadParameter("TimeStep", dT);
117 dT = dT * m_sampleFrequency;
118 // Calculate alpha
119 m_alpha = dT / (tau + dT);
120 }
121 }
122 else
123 {
124 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
125 m_alpha = equ.Evaluate();
126 // Check if tau was also defined
127 it = pParams.find("tau");
128 if (it != pParams.end())
129 {
130 ASSERTL0(false,
131 "Cannot define both alpha and tau in MovingAverage.");
132 }
133 }
134 // Check bounds of m_alpha
135 ASSERTL0(m_alpha > 0 && m_alpha < 1, "Alpha out of bounds.");
136}
137
139{
140}
141
144 const NekDouble &time)
145{
146 size_t dim = pFields.size() - 1;
147 size_t nExtraFields = (dim + 1) * dim / 2;
148 size_t origFields = pFields.size();
149 size_t nqtot = pFields[0]->GetTotPoints();
150
151 // Allocate storage
152 m_fields.resize(origFields + nExtraFields);
153 m_delta.resize(dim);
154
155 for (size_t n = 0; n < m_fields.size(); ++n)
156 {
157 m_fields[n] = Array<OneD, NekDouble>(nqtot, 0.0);
158 }
159 for (size_t n = 0; n < m_delta.size(); ++n)
160 {
161 m_delta[n] = Array<OneD, NekDouble>(nqtot, 0.0);
162 }
163
164 // Initialise output arrays
166
167 // Update m_fields if using restart file
168 if (m_numSamples)
169 {
170 for (size_t j = 0; j < m_fields.size(); ++j)
171 {
172 pFields[0]->BwdTrans(m_outFields[j], m_fields[j]);
173 if (pFields[0]->GetWaveSpace())
174 {
175 pFields[0]->HomogeneousBwdTrans(nqtot, m_fields[j],
176 m_fields[j]);
177 }
178 }
179 }
180}
181
184{
185 size_t dim = pFields.size() - 1;
186 size_t origFields = pFields.size();
187
188 // Fill name of variables
189 for (size_t n = 0; n < origFields; ++n)
190 {
191 m_variables.push_back(pFields[n]->GetSession()->GetVariable(n));
192 }
193 for (int i = 0; i < dim; ++i)
194 {
195 for (int j = i; j < dim; ++j)
196 {
197 std::string var = pFields[i]->GetSession()->GetVariable(i) +
198 pFields[j]->GetSession()->GetVariable(j);
199 m_variables.push_back(var);
200 }
201 }
202}
203
206 [[maybe_unused]] std::vector<Array<OneD, NekDouble>> &fieldcoeffs,
207 [[maybe_unused]] const NekDouble &time)
208{
209 size_t i, j, n;
210 size_t nq = pFields[0]->GetTotPoints();
211 size_t dim = pFields.size() - 1;
212 bool waveSpace = pFields[0]->GetWaveSpace();
213 NekDouble nSamples = (NekDouble)m_numSamples;
214
215 // For moving average, take first sample as initial vector
216 NekDouble alpha = m_alpha;
217 if (m_numSamples == 1)
218 {
219 alpha = 1.0;
220 }
221
222 // Define auxiliary constants for averages
223 NekDouble facOld, facAvg, facStress, facDelta;
224 if (m_movAvg)
225 {
226 facOld = 1.0 - alpha;
227 facAvg = alpha;
228 facStress = alpha;
229 facDelta = 1.0;
230 }
231 else
232 {
233 facOld = 1.0;
234 facAvg = 1.0;
235 facStress = nSamples / (nSamples - 1);
236 facDelta = 1.0 / nSamples;
237 }
238
241
242 // Update original velocities in phys space and calculate (\bar{u} - u_n)
243 for (n = 0; n < dim; ++n)
244 {
245 if (waveSpace)
246 {
247 pFields[n]->HomogeneousBwdTrans(nq, pFields[n]->GetPhys(), vel);
248 }
249 else
250 {
251 vel = pFields[n]->GetPhys();
252 }
253 Vmath::Svtsvtp(nq, facAvg, vel, 1, facOld, m_fields[n], 1, m_fields[n],
254 1);
255 Vmath::Svtvm(nq, facDelta, m_fields[n], 1, vel, 1, m_delta[n], 1);
256 }
257 // Update pressure (directly to outFields)
258 Vmath::Svtsvtp(m_outFields[dim].size(), facAvg, pFields[dim]->GetCoeffs(),
259 1, facOld, m_outFields[dim], 1, m_outFields[dim], 1);
260
261 // Ignore Reynolds stress for first sample (its contribution is zero)
262 if (m_numSamples == 1)
263 {
264 return;
265 }
266
267 // Calculate C_{n} = facOld * C_{n-1} + facStress * deltaI * deltaJ
268 for (i = 0, n = dim + 1; i < dim; ++i)
269 {
270 for (j = i; j < dim; ++j, ++n)
271 {
272 Vmath::Vmul(nq, m_delta[i], 1, m_delta[j], 1, tmp, 1);
273 Vmath::Svtsvtp(nq, facStress, tmp, 1, facOld, m_fields[n], 1,
274 m_fields[n], 1);
275 }
276 }
277}
278
281 [[maybe_unused]] const NekDouble &time)
282{
283 size_t dim = pFields.size() - 1;
284
285 m_fieldMetaData["NumberOfFieldDumps"] = std::to_string(m_numSamples);
286
287 // Set wavespace to false, as calculations were performed in physical space
288 bool waveSpace = pFields[0]->GetWaveSpace();
289 pFields[0]->SetWaveSpace(false);
290
291 // Forward transform and put into m_outFields (except pressure)
292 for (size_t i = 0; i < m_fields.size(); ++i)
293 {
294 if (i != dim)
295 {
296 pFields[0]->FwdTransLocalElmt(m_fields[i], m_outFields[i]);
297 }
298 }
299
300 // Restore waveSpace
301 pFields[0]->SetWaveSpace(waveSpace);
302}
303
305{
306 if (m_movAvg)
307 {
308 return 1.0;
309 }
310 else
311 {
312 return 1.0 / m_numSamples;
313 }
314}
315
316} // namespace Nektar::SolverUtils
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
LibUtilities::FieldMetaDataMap m_fieldMetaData
std::vector< Array< OneD, NekDouble > > m_outFields
SOLVER_UTILS_EXPORT void v_Initialise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
LibUtilities::SessionReaderSharedPtr m_session
Definition: Filter.h:83
SOLVER_UTILS_EXPORT FilterReynoldsStresses(const LibUtilities::SessionReaderSharedPtr &pSession, const std::shared_ptr< SolverUtils::EquationSystem > &pEquation, const std::map< std::string, std::string > &pParams)
static std::string className
Name of the class.
std::vector< Array< OneD, NekDouble > > m_delta
void v_Initialise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
void v_FillVariablesName(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields) override
static 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< Array< OneD, NekDouble > > m_fields
void v_PrepareOutput(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
SOLVER_UTILS_EXPORT ~FilterReynoldsStresses() override
void v_ProcessSample(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, const NekDouble &time) override
std::shared_ptr< SessionReader > SessionReaderSharedPtr
FilterFactory & GetFilterFactory()
double NekDouble
void Svtsvtp(int n, const T alpha, const T *x, int incx, const T beta, const T *y, int incy, T *z, int incz)
Svtsvtp (scalar times vector plus scalar times vector):
Definition: Vmath.hpp:473
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 Svtvm(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Svtvm (scalar times vector minus vector): z = alpha*x - y.
Definition: Vmath.hpp:424