Nektar++
FilterMovingAverage.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File FilterMovingAverage.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: Calculates exponential moving average of solution fields
32 // during time-stepping.
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
36 #include <boost/core/ignore_unused.hpp>
37 
39 
40 namespace Nektar
41 {
42 namespace SolverUtils
43 {
45  GetFilterFactory().RegisterCreatorFunction("MovingAverage",
47 
50  const std::weak_ptr<EquationSystem> &pEquation,
51  const ParamMap &pParams)
52  : FilterFieldConvert(pSession, pEquation, pParams)
53 {
54  // Load sampling frequency
55  auto it = pParams.find("SampleFrequency");
56  if (it == pParams.end())
57  {
59  }
60  else
61  {
63  m_session->GetInterpreter(), it->second);
64  m_sampleFrequency = round(equ.Evaluate());
65  }
66 
67  // Load filter parameter
68  it = pParams.find("alpha");
69  if (it == pParams.end())
70  {
71  it = pParams.find("tau");
72  if (it == pParams.end())
73  {
74  ASSERTL0(false, "FilterMovingAverage needs either alpha or tau.");
75  }
76  else
77  {
78  // Load time constant
80  m_session->GetInterpreter(), it->second);
81  NekDouble tau = equ.Evaluate();
82  // Load delta T between samples
83  NekDouble dT;
84  m_session->LoadParameter("TimeStep", dT);
85  dT = dT * m_sampleFrequency;
86  // Calculate alpha
87  m_alpha = dT / (tau + dT);
88  }
89  }
90  else
91  {
93  m_session->GetInterpreter(), it->second);
94  m_alpha = equ.Evaluate();
95  // Check if tau was also defined
96  it = pParams.find("tau");
97  if (it != pParams.end())
98  {
99  ASSERTL0(false,
100  "Cannot define both alpha and tau in MovingAverage.");
101  }
102  }
103  // Check bounds of m_alpha
104  ASSERTL0(m_alpha > 0 && m_alpha < 1, "Alpha out of bounds.");
105 }
106 
108 {
109 }
110 
113  std::vector<Array<OneD, NekDouble> > &fieldcoeffs,
114  const NekDouble &time)
115 {
116  boost::ignore_unused(pFields, time);
117 
118  // Take first sample as initial vector
119  NekDouble alpha = m_alpha;
120  if (m_numSamples == 1)
121  {
122  alpha = 1.0;
123  }
124 
125  // \bar{u}_n = alpha * u_n + (1-alpha) * \bar{u}_{n-1}
126  for (int n = 0; n < m_outFields.size(); ++n)
127  {
128  Vmath::Svtsvtp(m_outFields[n].num_elements(),
129  alpha,
130  fieldcoeffs[n],
131  1,
132  (1.0 - alpha),
133  m_outFields[n],
134  1,
135  m_outFields[n],
136  1);
137  }
138 }
139 
140 }
141 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
virtual SOLVER_UTILS_EXPORT ~FilterMovingAverage()
std::vector< Array< OneD, NekDouble > > m_outFields
virtual void v_ProcessSample(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, const NekDouble &time)
static std::string className
Name of the class.
NekDouble Evaluate() const
Definition: Equation.cpp:95
SOLVER_UTILS_EXPORT FilterMovingAverage(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const ParamMap &pParams)
double NekDouble
static FilterSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const std::map< std::string, std::string > &pParams)
Creates an instance of this class.
std::map< std::string, std::string > ParamMap
Definition: Filter.h:68
LibUtilities::SessionReaderSharedPtr m_session
Definition: Filter.h:86
FilterFactory & GetFilterFactory()
Definition: Filter.cpp:41
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:594
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:199
std::shared_ptr< SessionReader > SessionReaderSharedPtr