OpenMS
Loading...
Searching...
No Matches
SavitzkyGolayFilter.h
Go to the documentation of this file.
1// Copyright (c) 2002-present, OpenMS Inc. -- EKU Tuebingen, ETH Zurich, and FU Berlin
2// SPDX-License-Identifier: BSD-3-Clause
3//
4// --------------------------------------------------------------------------
5// $Maintainer: Timo Sachsenberg $
6// $Authors: Eva Lange $
7// --------------------------------------------------------------------------
8
9#pragma once
10
16
17namespace OpenMS
18{
76 class OPENMS_DLLAPI SavitzkyGolayFilter :
77 public ProgressLogger,
79 {
80public:
83
86
87 // low level template to filters spectra and chromatograms
88 // raw data and meta data needs to be copied to the output container before calling this function
89 template<class InputIt, class OutputIt>
90 void filter(InputIt first, InputIt last, OutputIt d_first)
91 {
92 size_t n = std::distance(first, last);
93
94 if (frame_size_ > n) { return; }
95
96 int i;
97 UInt j;
98 int mid = (frame_size_ / 2);
99 double help;
100
101 // compute the transient on
102 OutputIt out_it = d_first;
103
104 for (i = 0; i <= mid; ++i)
105 {
106 InputIt it_forward = (first - i);
107 help = 0;
108 for (j = 0; j < frame_size_; ++j)
109 {
110 help += it_forward->getIntensity() * coeffs_[(i + 1) * frame_size_ - 1 - j];
111 ++it_forward;
112 }
113
114 out_it->setPosition(first->getPosition());
115 out_it->setIntensity(std::max(0.0, help));
116 ++out_it;
117 ++first;
118 }
119
120 // compute the steady state output
121 InputIt it_help = (last - mid);
122
123 while (first != it_help)
124 {
125 InputIt it_forward = (first - mid);
126 help = 0;
127
128 for (j = 0; j < frame_size_; ++j)
129 {
130 help += it_forward->getIntensity() * coeffs_[mid * frame_size_ + j];
131 ++it_forward;
132 }
133
134 out_it->setPosition(first->getPosition());
135 out_it->setIntensity(std::max(0.0, help));
136 ++out_it;
137 ++first;
138 }
139
140 // compute the transient off
141 for (i = (mid - 1); i >= 0; --i)
142 {
143 InputIt it_forward = (first - (frame_size_ - i - 1));
144 help = 0;
145
146 for (j = 0; j < frame_size_; ++j)
147 {
148 help += it_forward->getIntensity() * coeffs_[i * frame_size_ + j];
149 ++it_forward;
150 }
151
152 out_it->setPosition(first->getPosition());
153 out_it->setIntensity(std::max(0.0, help));
154 ++out_it;
155 ++first;
156 }
157
158 }
159
163 void filter(MSSpectrum & spectrum)
164 {
165 // copy the data AND META DATA to the output container
166 MSSpectrum output = spectrum;
167 // filter
168 filter(spectrum.begin(), spectrum.end(), output.begin());
169 // swap back
170 std::swap(spectrum, output);
171 }
172
176 void filter(MSChromatogram & chromatogram)
177 {
178 // copy the data AND META DATA to the output container
179 MSChromatogram output = chromatogram;
180 // filter
181 filter(chromatogram.begin(), chromatogram.end(), output.begin());
182 // swap back
183 std::swap(chromatogram, output);
184 }
185
189 void filter(Mobilogram & mobilogram)
190 {
191 // copy the data AND META DATA to the output container
192 Mobilogram output = mobilogram;
193 // filter
194 filter(mobilogram.begin(), mobilogram.end(), output.begin());
195 // swap back
196 std::swap(mobilogram, output);
197 }
198
203 {
204 Size progress = 0;
205 startProgress(0, map.size() + map.getChromatograms().size(), "smoothing data");
206 for (Size i = 0; i < map.size(); ++i)
207 {
208 filter(map[i]);
209 setProgress(++progress);
210 }
211 for (Size i = 0; i < map.getChromatograms().size(); ++i)
212 {
213 filter(map.getChromatogram(i));
214 setProgress(++progress);
215 }
216 endProgress();
217 }
218
219protected:
221 std::vector<double> coeffs_;
222
225
228
229 // Docu in base class
230 void updateMembers_() override;
231 };
232
233} // namespace OpenMS
A base class for all classes handling default parameters.
Definition DefaultParamHandler.h:66
The representation of a chromatogram.
Definition MSChromatogram.h:30
In-Memory representation of a mass spectrometry run.
Definition MSExperiment.h:49
MSChromatogram & getChromatogram(Size id)
returns a single chromatogram
Size size() const noexcept
The number of spectra.
const std::vector< MSChromatogram > & getChromatograms() const
returns the chromatogram list
The representation of a 1D spectrum.
Definition MSSpectrum.h:44
The representation of a 1D ion mobilogram.
Definition Mobilogram.h:32
Iterator begin() noexcept
Definition Mobilogram.h:144
Iterator end() noexcept
Definition Mobilogram.h:157
Base class for all classes that want to report their progress.
Definition ProgressLogger.h:27
Computes the Savitzky-Golay filter coefficients using QR decomposition.
Definition SavitzkyGolayFilter.h:79
std::vector< double > coeffs_
Coefficients.
Definition SavitzkyGolayFilter.h:221
void filter(MSChromatogram &chromatogram)
Removed the noise from an MSChromatogram.
Definition SavitzkyGolayFilter.h:176
void filter(InputIt first, InputIt last, OutputIt d_first)
Definition SavitzkyGolayFilter.h:90
void filter(MSSpectrum &spectrum)
Removed the noise from an MSSpectrum containing profile data.
Definition SavitzkyGolayFilter.h:163
SavitzkyGolayFilter()
Constructor.
void filter(Mobilogram &mobilogram)
Removed the noise from an Mobilogram.
Definition SavitzkyGolayFilter.h:189
~SavitzkyGolayFilter() override
Destructor.
void filterExperiment(PeakMap &map)
Removed the noise from an MSExperiment containing profile data.
Definition SavitzkyGolayFilter.h:202
void updateMembers_() override
This method is used to update extra member variables at the end of the setParameters() method.
UInt frame_size_
UInt of the filter kernel (number of pre-tabulated coefficients)
Definition SavitzkyGolayFilter.h:224
UInt order_
The order of the smoothing polynomial.
Definition SavitzkyGolayFilter.h:227
unsigned int UInt
Unsigned integer type.
Definition Types.h:64
size_t Size
Size type e.g. used as variable which can hold result of size()
Definition Types.h:97
Main OpenMS namespace.
Definition openswathalgo/include/OpenMS/OPENSWATHALGO/DATAACCESS/ISpectrumAccess.h:19