OpenMS  2.5.0
SavitzkyGolayFilter.h
Go to the documentation of this file.
1 // --------------------------------------------------------------------------
2 // OpenMS -- Open-Source Mass Spectrometry
3 // --------------------------------------------------------------------------
4 // Copyright The OpenMS Team -- Eberhard Karls University Tuebingen,
5 // ETH Zurich, and Freie Universitaet Berlin 2002-2020.
6 //
7 // This software is released under a three-clause BSD license:
8 // * Redistributions of source code must retain the above copyright
9 // notice, this list of conditions and the following disclaimer.
10 // * Redistributions in binary form must reproduce the above copyright
11 // notice, this list of conditions and the following disclaimer in the
12 // documentation and/or other materials provided with the distribution.
13 // * Neither the name of any author or any participating institution
14 // may be used to endorse or promote products derived from this software
15 // without specific prior written permission.
16 // For a full list of authors, refer to the file AUTHORS.
17 // --------------------------------------------------------------------------
18 // THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
19 // AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
20 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
21 // ARE DISCLAIMED. IN NO EVENT SHALL ANY OF THE AUTHORS OR THE CONTRIBUTING
22 // INSTITUTIONS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
23 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
24 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS;
25 // OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY,
26 // WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR
27 // OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF
28 // ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
29 //
30 // --------------------------------------------------------------------------
31 // $Maintainer: Timo Sachsenberg $
32 // $Authors: Eva Lange $
33 // --------------------------------------------------------------------------
34 
35 #pragma once
36 
41 
42 namespace OpenMS
43 {
101  class OPENMS_DLLAPI SavitzkyGolayFilter :
102  public ProgressLogger,
103  public DefaultParamHandler
104  {
105 public:
108 
110  ~SavitzkyGolayFilter() override;
111 
112  // low level template to filters spectra and chromatograms
113  // raw data and meta data needs to be copied to the output container before calling this function
114  template<class InputIt, class OutputIt>
115  void filter(InputIt first, InputIt last, OutputIt d_first)
116  {
117  size_t n = std::distance(first, last);
118 
119  if (frame_size_ > n) { return; }
120 
121  int i;
122  UInt j;
123  int mid = (frame_size_ / 2);
124  double help;
125 
126  // compute the transient on
127  OutputIt out_it = d_first;
128 
129  for (i = 0; i <= mid; ++i)
130  {
131  InputIt it_forward = (first - i);
132  help = 0;
133  for (j = 0; j < frame_size_; ++j)
134  {
135  help += it_forward->getIntensity() * coeffs_[(i + 1) * frame_size_ - 1 - j];
136  ++it_forward;
137  }
138 
139  out_it->setPosition(first->getPosition());
140  out_it->setIntensity(std::max(0.0, help));
141  ++out_it;
142  ++first;
143  }
144 
145  // compute the steady state output
146  InputIt it_help = (last - mid);
147 
148  while (first != it_help)
149  {
150  InputIt it_forward = (first - mid);
151  help = 0;
152 
153  for (j = 0; j < frame_size_; ++j)
154  {
155  help += it_forward->getIntensity() * coeffs_[mid * frame_size_ + j];
156  ++it_forward;
157  }
158 
159  out_it->setPosition(first->getPosition());
160  out_it->setIntensity(std::max(0.0, help));
161  ++out_it;
162  ++first;
163  }
164 
165  // compute the transient off
166  for (i = (mid - 1); i >= 0; --i)
167  {
168  InputIt it_forward = (first - (frame_size_ - i - 1));
169  help = 0;
170 
171  for (j = 0; j < frame_size_; ++j)
172  {
173  help += it_forward->getIntensity() * coeffs_[i * frame_size_ + j];
174  ++it_forward;
175  }
176 
177  out_it->setPosition(first->getPosition());
178  out_it->setIntensity(std::max(0.0, help));
179  ++out_it;
180  ++first;
181  }
182 
183  }
184 
188  void filter(MSSpectrum & spectrum)
189  {
190  // copy the data AND META DATA to the output container
191  MSSpectrum output = spectrum;
192  // filter
193  filter(spectrum.begin(), spectrum.end(), output.begin());
194  // swap back
195  std::swap(spectrum, output);
196  }
197 
201  void filter(MSChromatogram & chromatogram)
202  {
203  // copy the data AND META DATA to the output container
204  MSChromatogram output = chromatogram;
205  // filter
206  filter(chromatogram.begin(), chromatogram.end(), output.begin());
207  // swap back
208  std::swap(chromatogram, output);
209  }
210 
215  {
216  Size progress = 0;
217  startProgress(0, map.size() + map.getChromatograms().size(), "smoothing data");
218  for (Size i = 0; i < map.size(); ++i)
219  {
220  filter(map[i]);
221  setProgress(++progress);
222  }
223  for (Size i = 0; i < map.getChromatograms().size(); ++i)
224  {
225  filter(map.getChromatogram(i));
226  setProgress(++progress);
227  }
228  endProgress();
229  }
230 
231 protected:
233  std::vector<double> coeffs_;
234 
237 
240 
241  // Docu in base class
242  void updateMembers_() override;
243  };
244 
245 } // namespace OpenMS
DefaultParamHandler.h
OpenMS::SavitzkyGolayFilter::filterExperiment
void filterExperiment(PeakMap &map)
Removed the noise from an MSExperiment containing profile data.
Definition: SavitzkyGolayFilter.h:214
OpenMS::SavitzkyGolayFilter::frame_size_
UInt frame_size_
UInt of the filter kernel (number of pre-tabulated coefficients)
Definition: SavitzkyGolayFilter.h:236
OpenMS::MSExperiment
In-Memory representation of a mass spectrometry experiment.
Definition: MSExperiment.h:77
OpenMS::Size
size_t Size
Size type e.g. used as variable which can hold result of size()
Definition: Types.h:127
OpenMS::MSExperiment::size
Size size() const
Definition: MSExperiment.h:127
OpenMS::DefaultParamHandler
A base class for all classes handling default parameters.
Definition: DefaultParamHandler.h:91
OpenMS
Main OpenMS namespace.
Definition: FeatureDeconvolution.h:46
OpenMS::SavitzkyGolayFilter::filter
void filter(MSSpectrum &spectrum)
Removed the noise from an MSSpectrum containing profile data.
Definition: SavitzkyGolayFilter.h:188
OpenMS::ProgressLogger
Base class for all classes that want to report their progress.
Definition: ProgressLogger.h:54
OpenMS::MSExperiment::getChromatograms
const std::vector< MSChromatogram > & getChromatograms() const
returns the chromatogram list
ProgressLogger.h
OpenMS::UInt
unsigned int UInt
Unsigned integer type.
Definition: Types.h:94
MSExperiment.h
OpenMS::SavitzkyGolayFilter::filter
void filter(InputIt first, InputIt last, OutputIt d_first)
Definition: SavitzkyGolayFilter.h:115
OpenMS::SavitzkyGolayFilter::order_
UInt order_
The order of the smoothing polynomial.
Definition: SavitzkyGolayFilter.h:239
OpenMS::MSExperiment::getChromatogram
MSChromatogram & getChromatogram(Size id)
returns a single chromatogram
OpenMS::SavitzkyGolayFilter::coeffs_
std::vector< double > coeffs_
Coefficients.
Definition: SavitzkyGolayFilter.h:233
OpenMS::MSChromatogram
The representation of a chromatogram.
Definition: MSChromatogram.h:54
OpenMS::SavitzkyGolayFilter
Computes the Savitzky-Golay filter coefficients using QR decomposition.
Definition: SavitzkyGolayFilter.h:101
OpenMS::MSSpectrum
The representation of a 1D spectrum.
Definition: MSSpectrum.h:67
StandardTypes.h
OpenMS::SavitzkyGolayFilter::filter
void filter(MSChromatogram &chromatogram)
Removed the noise from an MSChromatogram.
Definition: SavitzkyGolayFilter.h:201