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
40namespace Nektar
41{
42namespace SolverUtils
43{
47
50 const std::weak_ptr<EquationSystem> &pEquation, const ParamMap &pParams)
51 : FilterFieldConvert(pSession, pEquation, pParams)
52{
53 // Load sampling frequency
54 auto it = pParams.find("SampleFrequency");
55 if (it == pParams.end())
56 {
58 }
59 else
60 {
61 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
62 m_sampleFrequency = round(equ.Evaluate());
63 }
64
65 // Load filter parameter
66 it = pParams.find("alpha");
67 if (it == pParams.end())
68 {
69 it = pParams.find("tau");
70 if (it == pParams.end())
71 {
72 ASSERTL0(false, "FilterMovingAverage needs either alpha or tau.");
73 }
74 else
75 {
76 // Load time constant
77 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
78 NekDouble tau = equ.Evaluate();
79 // Load delta T between samples
80 NekDouble dT;
81 m_session->LoadParameter("TimeStep", dT);
82 dT = dT * m_sampleFrequency;
83 // Calculate alpha
84 m_alpha = dT / (tau + dT);
85 }
86 }
87 else
88 {
89 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
90 m_alpha = equ.Evaluate();
91 // Check if tau was also defined
92 it = pParams.find("tau");
93 if (it != pParams.end())
94 {
95 ASSERTL0(false,
96 "Cannot define both alpha and tau in MovingAverage.");
97 }
98 }
99 // Check bounds of m_alpha
100 ASSERTL0(m_alpha > 0 && m_alpha < 1, "Alpha out of bounds.");
101}
102
104{
105}
106
109 std::vector<Array<OneD, NekDouble>> &fieldcoeffs, const NekDouble &time)
110{
111 boost::ignore_unused(pFields, time);
112
113 // Take first sample as initial vector
114 NekDouble alpha = m_alpha;
115 if (m_numSamples == 1)
116 {
117 alpha = 1.0;
118 }
119
120 // \bar{u}_n = alpha * u_n + (1-alpha) * \bar{u}_{n-1}
121 for (int n = 0; n < m_outFields.size(); ++n)
122 {
123 Vmath::Svtsvtp(m_outFields[n].size(), alpha, fieldcoeffs[n], 1,
124 (1.0 - alpha), m_outFields[n], 1, m_outFields[n], 1);
125 }
126}
127
128} // namespace SolverUtils
129} // 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
std::vector< Array< OneD, NekDouble > > m_outFields
LibUtilities::SessionReaderSharedPtr m_session
Definition: Filter.h:85
std::map< std::string, std::string > ParamMap
Definition: Filter.h:67
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< EquationSystem > &pEquation, const std::map< std::string, std::string > &pParams)
Creates an instance of this class.
virtual SOLVER_UTILS_EXPORT ~FilterMovingAverage()
SOLVER_UTILS_EXPORT FilterMovingAverage(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const ParamMap &pParams)
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:746