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  int dim = pFields.size() - 1;
149  int nExtraFields = dim == 2 ? 3 : 6;
150  int origFields = pFields.size();
151 
152  // Allocate storage
153  m_fields.resize(origFields + nExtraFields);
154  m_delta.resize(dim);
155 
156  for (int n = 0; n < m_fields.size(); ++n)
157  {
158  m_fields[n] = Array<OneD, NekDouble>(pFields[0]->GetTotPoints(), 0.0);
159  }
160  for (int n = 0; n < m_delta.size(); ++n)
161  {
162  m_delta[n] = Array<OneD, NekDouble>(pFields[0]->GetTotPoints(), 0.0);
163  }
164 
165  // Initialise output arrays
166  FilterFieldConvert::v_Initialise(pFields, time);
167 
168  // Update m_fields if using restart file
169  if (m_numSamples)
170  {
171  for (int j = 0; j < m_fields.size(); ++j)
172  {
173  pFields[0]->BwdTrans(m_outFields[j], m_fields[j]);
174  if (pFields[0]->GetWaveSpace())
175  {
176  pFields[0]->HomogeneousBwdTrans(m_fields[j], m_fields[j]);
177  }
178  }
179  }
180 }
181 
184 {
185  int dim = pFields.size() - 1;
186  int origFields = pFields.size();
187 
188  // Fill name of variables
189  for (int n = 0; n < origFields; ++n)
190  {
191  m_variables.push_back(pFields[n]->GetSession()->GetVariable(n));
192  }
193  if (dim == 2)
194  {
195  m_variables.push_back("uu");
196  m_variables.push_back("uv");
197  m_variables.push_back("vv");
198  }
199  else if (dim == 3)
200  {
201  m_variables.push_back("uu");
202  m_variables.push_back("uv");
203  m_variables.push_back("uw");
204  m_variables.push_back("vv");
205  m_variables.push_back("vw");
206  m_variables.push_back("ww");
207  }
208  else
209  {
210  ASSERTL0(false, "Unsupported dimension");
211  }
212 }
213 
216  std::vector<Array<OneD, NekDouble> > &fieldcoeffs,
217  const NekDouble &time)
218 {
219  int i, j, n;
220  int nq = pFields[0]->GetTotPoints();
221  int dim = pFields.size() - 1;
222  bool waveSpace = pFields[0]->GetWaveSpace();
223  NekDouble nSamples = (NekDouble)m_numSamples;
224 
225  // For moving average, take first sample as initial vector
226  NekDouble alpha = m_alpha;
227  if (m_numSamples == 1)
228  {
229  alpha = 1.0;
230  }
231 
232  // Define auxiliary constants for averages
233  NekDouble facOld, facAvg, facStress, facDelta;
234  if (m_movAvg)
235  {
236  facOld = 1.0 - alpha;
237  facAvg = alpha;
238  facStress = alpha;
239  facDelta = 1.0;
240  }
241  else
242  {
243  facOld = 1.0;
244  facAvg = 1.0;
245  facStress = nSamples / (nSamples - 1);
246  facDelta = 1.0 / nSamples;
247  }
248 
249  Array<OneD, NekDouble> vel(nq);
250  Array<OneD, NekDouble> tmp(nq);
251 
252  // Update original velocities in phys space and calculate (\bar{u} - u_n)
253  for (n = 0; n < dim; ++n)
254  {
255  if (waveSpace)
256  {
257  pFields[n]->HomogeneousBwdTrans(pFields[n]->GetPhys(), vel);
258  }
259  else
260  {
261  vel = pFields[n]->GetPhys();
262  }
264  nq, facAvg, vel, 1, facOld, m_fields[n], 1, m_fields[n], 1);
265  Vmath::Svtvm(nq, facDelta, m_fields[n], 1, vel, 1, m_delta[n], 1);
266  }
267  // Update pressure (directly to outFields)
268  Vmath::Svtsvtp(m_outFields[dim].size(),
269  facAvg,
270  pFields[dim]->GetCoeffs(),
271  1,
272  facOld,
273  m_outFields[dim],
274  1,
275  m_outFields[dim],
276  1);
277 
278  // Ignore Reynolds stress for first sample (its contribution is zero)
279  if (m_numSamples == 1)
280  {
281  return;
282  }
283 
284  // Calculate C_{n} = facOld * C_{n-1} + facStress * deltaI * deltaJ
285  for (i = 0, n = dim + 1; i < dim; ++i)
286  {
287  for (j = i; j < dim; ++j, ++n)
288  {
289  Vmath::Vmul(nq, m_delta[i], 1, m_delta[j], 1, tmp, 1);
291  nq, facStress, tmp, 1, facOld, m_fields[n], 1, m_fields[n], 1);
292  }
293  }
294 }
295 
298  const NekDouble &time)
299 {
300  int dim = pFields.size() - 1;
301 
302  m_fieldMetaData["NumberOfFieldDumps"] =
303  boost::lexical_cast<std::string>(m_numSamples);
304 
305  // Set wavespace to false, as calculations were performed in physical space
306  bool waveSpace = pFields[0]->GetWaveSpace();
307  pFields[0]->SetWaveSpace(false);
308 
309  // Forward transform and put into m_outFields (except pressure)
310  for (int i = 0; i < m_fields.size(); ++i)
311  {
312  if (i != dim)
313  {
314  pFields[0]->FwdTrans_IterPerExp(m_fields[i], m_outFields[i]);
315  }
316  }
317 
318  // Restore waveSpace
319  pFields[0]->SetWaveSpace(waveSpace);
320 }
321 
323 {
324  if (m_movAvg)
325  {
326  return 1.0;
327  }
328  else
329  {
330  return 1.0 / m_numSamples;
331  }
332 }
333 
334 }
335 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
NekDouble Evaluate() const
Definition: Equation.cpp:95
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:200
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)
LibUtilities::SessionReaderSharedPtr m_session
Definition: Filter.h:86
virtual void v_ProcessSample(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, const NekDouble &time)
virtual void v_Initialise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
SOLVER_UTILS_EXPORT FilterReynoldsStresses(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< SolverUtils::EquationSystem > &pEquation, const std::map< std::string, std::string > &pParams)
virtual void v_FillVariablesName(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields)
virtual void v_PrepareOutput(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
static std::string className
Name of the class.
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
std::vector< Array< OneD, NekDouble > > m_fields
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:1
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)
vvtvvtp (scalar times vector plus scalar times vector):
Definition: Vmath.cpp:691
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:192
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 plus vector): z = alpha*x - y
Definition: Vmath.cpp:602