OpenMS
SavitzkyGolayFilter.h
Go to the documentation of this file.
1 // Copyright (c) 2002-present, The OpenMS Team -- 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 
15 
16 namespace OpenMS
17 {
75  class OPENMS_DLLAPI SavitzkyGolayFilter :
76  public ProgressLogger,
77  public DefaultParamHandler
78  {
79 public:
82 
85 
86  // low level template to filters spectra and chromatograms
87  // raw data and meta data needs to be copied to the output container before calling this function
88  template<class InputIt, class OutputIt>
89  void filter(InputIt first, InputIt last, OutputIt d_first)
90  {
91  size_t n = std::distance(first, last);
92 
93  if (frame_size_ > n) { return; }
94 
95  int i;
96  UInt j;
97  int mid = (frame_size_ / 2);
98  double help;
99 
100  // compute the transient on
101  OutputIt out_it = d_first;
102 
103  for (i = 0; i <= mid; ++i)
104  {
105  InputIt it_forward = (first - i);
106  help = 0;
107  for (j = 0; j < frame_size_; ++j)
108  {
109  help += it_forward->getIntensity() * coeffs_[(i + 1) * frame_size_ - 1 - j];
110  ++it_forward;
111  }
112 
113  out_it->setPosition(first->getPosition());
114  out_it->setIntensity(std::max(0.0, help));
115  ++out_it;
116  ++first;
117  }
118 
119  // compute the steady state output
120  InputIt it_help = (last - mid);
121 
122  while (first != it_help)
123  {
124  InputIt it_forward = (first - mid);
125  help = 0;
126 
127  for (j = 0; j < frame_size_; ++j)
128  {
129  help += it_forward->getIntensity() * coeffs_[mid * frame_size_ + j];
130  ++it_forward;
131  }
132 
133  out_it->setPosition(first->getPosition());
134  out_it->setIntensity(std::max(0.0, help));
135  ++out_it;
136  ++first;
137  }
138 
139  // compute the transient off
140  for (i = (mid - 1); i >= 0; --i)
141  {
142  InputIt it_forward = (first - (frame_size_ - i - 1));
143  help = 0;
144 
145  for (j = 0; j < frame_size_; ++j)
146  {
147  help += it_forward->getIntensity() * coeffs_[i * frame_size_ + j];
148  ++it_forward;
149  }
150 
151  out_it->setPosition(first->getPosition());
152  out_it->setIntensity(std::max(0.0, help));
153  ++out_it;
154  ++first;
155  }
156 
157  }
158 
162  void filter(MSSpectrum & spectrum)
163  {
164  // copy the data AND META DATA to the output container
165  MSSpectrum output = spectrum;
166  // filter
167  filter(spectrum.begin(), spectrum.end(), output.begin());
168  // swap back
169  std::swap(spectrum, output);
170  }
171 
175  void filter(MSChromatogram & chromatogram)
176  {
177  // copy the data AND META DATA to the output container
178  MSChromatogram output = chromatogram;
179  // filter
180  filter(chromatogram.begin(), chromatogram.end(), output.begin());
181  // swap back
182  std::swap(chromatogram, output);
183  }
184 
189  {
190  Size progress = 0;
191  startProgress(0, map.size() + map.getChromatograms().size(), "smoothing data");
192  for (Size i = 0; i < map.size(); ++i)
193  {
194  filter(map[i]);
195  setProgress(++progress);
196  }
197  for (Size i = 0; i < map.getChromatograms().size(); ++i)
198  {
199  filter(map.getChromatogram(i));
200  setProgress(++progress);
201  }
202  endProgress();
203  }
204 
205 protected:
207  std::vector<double> coeffs_;
208 
211 
214 
215  // Docu in base class
216  void updateMembers_() override;
217  };
218 
219 } // 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:45
const std::vector< MSChromatogram > & getChromatograms() const
returns the chromatogram list
MSChromatogram & getChromatogram(Size id)
returns a single chromatogram
Size size() const noexcept
The number of spectra.
Definition: MSExperiment.h:120
The representation of a 1D spectrum.
Definition: MSSpectrum.h:44
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:78
std::vector< double > coeffs_
Coefficients.
Definition: SavitzkyGolayFilter.h:207
void filter(MSChromatogram &chromatogram)
Removed the noise from an MSChromatogram.
Definition: SavitzkyGolayFilter.h:175
void filter(InputIt first, InputIt last, OutputIt d_first)
Definition: SavitzkyGolayFilter.h:89
void filter(MSSpectrum &spectrum)
Removed the noise from an MSSpectrum containing profile data.
Definition: SavitzkyGolayFilter.h:162
SavitzkyGolayFilter()
Constructor.
~SavitzkyGolayFilter() override
Destructor.
void filterExperiment(PeakMap &map)
Removed the noise from an MSExperiment containing profile data.
Definition: SavitzkyGolayFilter.h:188
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:210
UInt order_
The order of the smoothing polynomial.
Definition: SavitzkyGolayFilter.h:213
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