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 
37 namespace Nektar
38 {
39 namespace SolverUtils
40 {
42  GetFilterFactory().RegisterCreatorFunction("ReynoldsStresses",
44 
45 /**
46  * @class FilterReynoldsStresses
47  *
48  * @brief Append Reynolds stresses to the average fields
49  *
50  * This class appends the average fields with the Reynolds stresses of the form
51  * \f$ \overline{u' v'} \f$.
52  *
53  * For the default case, this is achieved by calculating
54  * \f$ C_{n} = \Sigma_{i=1}^{n} (u_i - \bar{u}_n)(v_i - \bar{v}_n)\f$
55  * using the recursive relation:
56  *
57  * \f[ C_{n} = C_{n-1} + \frac{n}{n-1} (u_n - \bar{u}_n)(v_n - \bar{v}_n) \f]
58  *
59  * The FilterSampler base class then divides the result by n, leading
60  * to the Reynolds stress.
61  *
62  * It is also possible to perform the averages using an exponential moving
63  * average, in which case either the moving average parameter \f$ \alpha \f$
64  * or the time constant \f$ \tau \f$ must be prescribed.
65  */
68  const std::weak_ptr<SolverUtils::EquationSystem> &pEquation,
69  const std::map<std::string, std::string> &pParams)
70  : FilterFieldConvert(pSession, pEquation, pParams)
71 {
72  // Load sampling frequency
73  auto it = pParams.find("SampleFrequency");
74  if (it == pParams.end())
75  {
77  }
78  else
79  {
80  LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
81  m_sampleFrequency = round(equ.Evaluate());
82  }
83 
84  // Check if should use moving average
85  it = pParams.find("MovingAverage");
86  if (it == pParams.end())
87  {
88  m_movAvg = false;
89  }
90  else
91  {
92  std::string sOption = it->second.c_str();
93  m_movAvg = (boost::iequals(sOption, "true")) ||
94  (boost::iequals(sOption, "yes"));
95  }
96 
97  if (!m_movAvg)
98  {
99  return;
100  }
101 
102  // Load alpha parameter for moving average
103  it = pParams.find("alpha");
104  if (it == pParams.end())
105  {
106  it = pParams.find("tau");
107  if (it == pParams.end())
108  {
109  ASSERTL0(false, "MovingAverage needs either alpha or tau.");
110  }
111  else
112  {
113  // Load time constant
114  LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
115  NekDouble tau = equ.Evaluate();
116  // Load delta T between samples
117  NekDouble dT;
118  m_session->LoadParameter("TimeStep", dT);
119  dT = dT * m_sampleFrequency;
120  // Calculate alpha
121  m_alpha = dT / (tau + dT);
122  }
123  }
124  else
125  {
126  LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
127  m_alpha = equ.Evaluate();
128  // Check if tau was also defined
129  it = pParams.find("tau");
130  if (it != pParams.end())
131  {
132  ASSERTL0(false,
133  "Cannot define both alpha and tau in MovingAverage.");
134  }
135  }
136  // Check bounds of m_alpha
137  ASSERTL0(m_alpha > 0 && m_alpha < 1, "Alpha out of bounds.");
138 }
139 
141 {
142 }
143 
146  const NekDouble &time)
147 {
148  size_t dim = pFields.size() - 1;
149  size_t nExtraFields = (dim + 1) * dim / 2;
150  size_t origFields = pFields.size();
151  size_t nqtot = pFields[0]->GetTotPoints();
152 
153  // Allocate storage
154  m_fields.resize(origFields + nExtraFields);
155  m_delta.resize(dim);
156 
157  for (size_t n = 0; n < m_fields.size(); ++n)
158  {
159  m_fields[n] = Array<OneD, NekDouble>(nqtot, 0.0);
160  }
161  for (size_t n = 0; n < m_delta.size(); ++n)
162  {
163  m_delta[n] = Array<OneD, NekDouble>(nqtot, 0.0);
164  }
165 
166  // Initialise output arrays
167  FilterFieldConvert::v_Initialise(pFields, time);
168 
169  // Update m_fields if using restart file
170  if (m_numSamples)
171  {
172  for (size_t j = 0; j < m_fields.size(); ++j)
173  {
174  pFields[0]->BwdTrans(m_outFields[j], m_fields[j]);
175  if (pFields[0]->GetWaveSpace())
176  {
177  pFields[0]->HomogeneousBwdTrans(nqtot, m_fields[j],
178  m_fields[j]);
179  }
180  }
181  }
182 }
183 
186 {
187  size_t dim = pFields.size() - 1;
188  size_t origFields = pFields.size();
189 
190  // Fill name of variables
191  for (size_t n = 0; n < origFields; ++n)
192  {
193  m_variables.push_back(pFields[n]->GetSession()->GetVariable(n));
194  }
195  for (int i = 0; i < dim; ++i)
196  {
197  for (int j = i; j < dim; ++j)
198  {
199  std::string var = pFields[i]->GetSession()->GetVariable(i) +
200  pFields[j]->GetSession()->GetVariable(j);
201  m_variables.push_back(var);
202  }
203  }
204 }
205 
208  std::vector<Array<OneD, NekDouble>> &fieldcoeffs, const NekDouble &time)
209 {
210  boost::ignore_unused(fieldcoeffs, time);
211 
212  size_t i, j, n;
213  size_t nq = pFields[0]->GetTotPoints();
214  size_t dim = pFields.size() - 1;
215  bool waveSpace = pFields[0]->GetWaveSpace();
216  NekDouble nSamples = (NekDouble)m_numSamples;
217 
218  // For moving average, take first sample as initial vector
219  NekDouble alpha = m_alpha;
220  if (m_numSamples == 1)
221  {
222  alpha = 1.0;
223  }
224 
225  // Define auxiliary constants for averages
226  NekDouble facOld, facAvg, facStress, facDelta;
227  if (m_movAvg)
228  {
229  facOld = 1.0 - alpha;
230  facAvg = alpha;
231  facStress = alpha;
232  facDelta = 1.0;
233  }
234  else
235  {
236  facOld = 1.0;
237  facAvg = 1.0;
238  facStress = nSamples / (nSamples - 1);
239  facDelta = 1.0 / nSamples;
240  }
241 
242  Array<OneD, NekDouble> vel(nq);
243  Array<OneD, NekDouble> tmp(nq);
244 
245  // Update original velocities in phys space and calculate (\bar{u} - u_n)
246  for (n = 0; n < dim; ++n)
247  {
248  if (waveSpace)
249  {
250  pFields[n]->HomogeneousBwdTrans(nq, pFields[n]->GetPhys(), vel);
251  }
252  else
253  {
254  vel = pFields[n]->GetPhys();
255  }
256  Vmath::Svtsvtp(nq, facAvg, vel, 1, facOld, m_fields[n], 1, m_fields[n],
257  1);
258  Vmath::Svtvm(nq, facDelta, m_fields[n], 1, vel, 1, m_delta[n], 1);
259  }
260  // Update pressure (directly to outFields)
261  Vmath::Svtsvtp(m_outFields[dim].size(), facAvg, pFields[dim]->GetCoeffs(),
262  1, facOld, m_outFields[dim], 1, m_outFields[dim], 1);
263 
264  // Ignore Reynolds stress for first sample (its contribution is zero)
265  if (m_numSamples == 1)
266  {
267  return;
268  }
269 
270  // Calculate C_{n} = facOld * C_{n-1} + facStress * deltaI * deltaJ
271  for (i = 0, n = dim + 1; i < dim; ++i)
272  {
273  for (j = i; j < dim; ++j, ++n)
274  {
275  Vmath::Vmul(nq, m_delta[i], 1, m_delta[j], 1, tmp, 1);
276  Vmath::Svtsvtp(nq, facStress, tmp, 1, facOld, m_fields[n], 1,
277  m_fields[n], 1);
278  }
279  }
280 }
281 
284  const NekDouble &time)
285 {
286  boost::ignore_unused(time);
287 
288  size_t dim = pFields.size() - 1;
289 
290  m_fieldMetaData["NumberOfFieldDumps"] =
291  boost::lexical_cast<std::string>(m_numSamples);
292 
293  // Set wavespace to false, as calculations were performed in physical space
294  bool waveSpace = pFields[0]->GetWaveSpace();
295  pFields[0]->SetWaveSpace(false);
296 
297  // Forward transform and put into m_outFields (except pressure)
298  for (size_t i = 0; i < m_fields.size(); ++i)
299  {
300  if (i != dim)
301  {
302  pFields[0]->FwdTransLocalElmt(m_fields[i], m_outFields[i]);
303  }
304  }
305 
306  // Restore waveSpace
307  pFields[0]->SetWaveSpace(waveSpace);
308 }
309 
311 {
312  if (m_movAvg)
313  {
314  return 1.0;
315  }
316  else
317  {
318  return 1.0 / m_numSamples;
319  }
320 }
321 
322 } // namespace SolverUtils
323 } // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
LibUtilities::FieldMetaDataMap m_fieldMetaData
std::vector< Array< OneD, NekDouble > > m_outFields
virtual SOLVER_UTILS_EXPORT void v_Initialise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
LibUtilities::SessionReaderSharedPtr m_session
Definition: Filter.h:85
SOLVER_UTILS_EXPORT FilterReynoldsStresses(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< SolverUtils::EquationSystem > &pEquation, const std::map< std::string, std::string > &pParams)
static std::string className
Name of the class.
virtual void v_ProcessSample(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, std::vector< Array< OneD, NekDouble >> &fieldcoeffs, const NekDouble &time) override
static 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.
std::vector< Array< OneD, NekDouble > > m_delta
virtual void v_Initialise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
virtual void v_FillVariablesName(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields) override
std::vector< Array< OneD, NekDouble > > m_fields
virtual void v_PrepareOutput(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time) override
std::shared_ptr< SessionReader > SessionReaderSharedPtr
FilterFactory & GetFilterFactory()
Definition: Filter.cpp:41
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
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)
svtvvtp (scalar times vector plus scalar times vector):
Definition: Vmath.cpp:751
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.cpp:209
void Svtvm(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
svtvp (scalar times vector minus vector): z = alpha*x - y
Definition: Vmath.cpp:664