Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 // 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: Calculates exponential moving average of solution fields
33 // during time-stepping.
34 //
35 ///////////////////////////////////////////////////////////////////////////////
36 
38 
39 namespace Nektar
40 {
41 namespace SolverUtils
42 {
44  GetFilterFactory().RegisterCreatorFunction("MovingAverage",
46 
49  const ParamMap &pParams)
50  : FilterSampler(pSession, pParams)
51 {
52  ParamMap::const_iterator it;
53 
54  // Load filter parameter
55  it = pParams.find("alpha");
56  if (it == pParams.end())
57  {
58  it = pParams.find("tau");
59  if (it == pParams.end())
60  {
61  ASSERTL0(false, "FilterMovingAverage needs either alpha or tau.");
62  }
63  else
64  {
65  // Load time constant
66  LibUtilities::Equation equ(m_session, it->second);
67  NekDouble tau = equ.Evaluate();
68  // Load delta T between samples
69  NekDouble dT;
70  m_session->LoadParameter("TimeStep", dT);
71  dT = dT * m_sampleFrequency;
72  // Calculate alpha
73  m_alpha = dT / (tau + dT);
74  }
75  }
76  else
77  {
78  LibUtilities::Equation equ(m_session, it->second);
79  m_alpha = equ.Evaluate();
80  // Check if tau was also defined
81  it = pParams.find("tau");
82  if (it != pParams.end())
83  {
84  ASSERTL0(false,
85  "Cannot define both alpha and tau in MovingAverage.");
86  }
87  }
88  // Check bounds of m_alpha
89  ASSERTL0(m_alpha > 0 && m_alpha < 1, "Alpha out of bounds.");
90 }
91 
93 {
94 }
95 
98  const NekDouble &time)
99 {
100  int nfield = pFields.num_elements();
101  m_variables.resize(pFields.num_elements());
102  // Fill name of variables
103  for (int n = 0; n < nfield; ++n)
104  {
105  m_variables[n] = pFields[n]->GetSession()->GetVariable(n);
106  }
107  // Now let FilterSampler initialise the output
108  FilterSampler::v_Initialise(pFields, time);
109 }
110 
113  const NekDouble &time)
114 {
115  // Take first sample as initial vector
116  NekDouble alpha = m_alpha;
117  if (m_numSamples == 1)
118  {
119  alpha = 1.0;
120  }
121 
122  // \bar{u}_n = alpha * u_n + (1-alpha) * \bar{u}_{n-1}
123  for (int n = 0; n < pFields.num_elements(); ++n)
124  {
125  Vmath::Svtsvtp(m_outFields[n].num_elements(),
126  alpha,
127  pFields[n]->GetCoeffs(),
128  1,
129  (1.0 - alpha),
130  m_outFields[n],
131  1,
132  m_outFields[n],
133  1);
134  }
135 }
136 
138 {
139  return true;
140 }
141 }
142 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
std::vector< Array< OneD, NekDouble > > m_outFields
Definition: FilterSampler.h:90
virtual SOLVER_UTILS_EXPORT ~FilterMovingAverage()
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
static std::string className
Name of the class.
std::vector< std::string > m_variables
Definition: FilterSampler.h:91
static FilterSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const std::map< std::string, std::string > &pParams)
Creates an instance of this class.
NekDouble Evaluate() const
Definition: Equation.h:102
virtual SOLVER_UTILS_EXPORT void v_Initialise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
virtual void v_Initialise(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
double NekDouble
std::map< std::string, std::string > ParamMap
Definition: Filter.h:67
LibUtilities::SessionReaderSharedPtr m_session
Definition: Filter.h:84
SOLVER_UTILS_EXPORT FilterMovingAverage(const LibUtilities::SessionReaderSharedPtr &pSession, const ParamMap &pParams)
FilterFactory & GetFilterFactory()
Definition: Filter.cpp:42
virtual void v_ProcessSample(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, const NekDouble &time)
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:577
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:215