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
37
38namespace Nektar::SolverUtils
39{
43
46 const std::shared_ptr<EquationSystem> &pEquation, const ParamMap &pParams)
47 : FilterFieldConvert(pSession, pEquation, pParams)
48{
49 // Load sampling frequency
50 auto it = pParams.find("SampleFrequency");
51 if (it == pParams.end())
52 {
54 }
55 else
56 {
57 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
58 m_sampleFrequency = round(equ.Evaluate());
59 }
60
61 // Load filter parameter
62 it = pParams.find("alpha");
63 if (it == pParams.end())
64 {
65 it = pParams.find("tau");
66 if (it == pParams.end())
67 {
68 ASSERTL0(false, "FilterMovingAverage needs either alpha or tau.");
69 }
70 else
71 {
72 // Load time constant
73 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
74 NekDouble tau = equ.Evaluate();
75 // Load delta T between samples
76 NekDouble dT;
77 m_session->LoadParameter("TimeStep", dT);
78 dT = dT * m_sampleFrequency;
79 // Calculate alpha
80 m_alpha = dT / (tau + dT);
81 }
82 }
83 else
84 {
85 LibUtilities::Equation equ(m_session->GetInterpreter(), it->second);
86 m_alpha = equ.Evaluate();
87 // Check if tau was also defined
88 it = pParams.find("tau");
89 if (it != pParams.end())
90 {
91 ASSERTL0(false,
92 "Cannot define both alpha and tau in MovingAverage.");
93 }
94 }
95 // Check bounds of m_alpha
96 ASSERTL0(m_alpha > 0 && m_alpha < 1, "Alpha out of bounds.");
97}
98
100{
101}
102
105 &pFields,
106 std::vector<Array<OneD, NekDouble>> &fieldcoeffs,
107 [[maybe_unused]] const NekDouble &time)
108{
109 // Take first sample as initial vector
110 NekDouble alpha = m_alpha;
111 if (m_numSamples == 1)
112 {
113 alpha = 1.0;
114 }
115
116 // \bar{u}_n = alpha * u_n + (1-alpha) * \bar{u}_{n-1}
117 for (int n = 0; n < m_outFields.size(); ++n)
118 {
119 Vmath::Svtsvtp(m_outFields[n].size(), alpha, fieldcoeffs[n], 1,
120 (1.0 - alpha), m_outFields[n], 1, m_outFields[n], 1);
121 }
122}
123
124} // namespace Nektar::SolverUtils
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
std::vector< Array< OneD, NekDouble > > m_outFields
LibUtilities::SessionReaderSharedPtr m_session
Definition: Filter.h:83
std::map< std::string, std::string > ParamMap
Definition: Filter.h:65
static std::string className
Name of the class.
static FilterSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const std::shared_ptr< EquationSystem > &pEquation, const std::map< std::string, std::string > &pParams)
Creates an instance of this class.
void v_ProcessSample(const Array< OneD, const MultiRegions::ExpListSharedPtr > &pFields, std::vector< Array< OneD, NekDouble > > &fieldcoeffs, const NekDouble &time) override
SOLVER_UTILS_EXPORT ~FilterMovingAverage() override
SOLVER_UTILS_EXPORT FilterMovingAverage(const LibUtilities::SessionReaderSharedPtr &pSession, const std::shared_ptr< EquationSystem > &pEquation, const ParamMap &pParams)
std::shared_ptr< SessionReader > SessionReaderSharedPtr
FilterFactory & GetFilterFactory()
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)
Svtsvtp (scalar times vector plus scalar times vector):
Definition: Vmath.hpp:473