|  | OpenMS
    2.6.0
    | 
 
 
  
  
 
Go to the documentation of this file.
  114     template<
class InputIt, 
class OutputIt>
 
  115     void filter(InputIt first, InputIt last, OutputIt d_first)
 
  117       size_t n = std::distance(first, last);
 
  119       if (frame_size_ > n) { 
return; }
 
  123       int mid = (frame_size_ / 2);
 
  127       OutputIt out_it = d_first;
 
  129       for (i = 0; i <= mid; ++i)
 
  131         InputIt it_forward = (first - i);
 
  133         for (j = 0; j < frame_size_; ++j)
 
  135           help += it_forward->getIntensity() * coeffs_[(i + 1) * frame_size_ - 1 - j];
 
  139         out_it->setPosition(first->getPosition());
 
  140         out_it->setIntensity(std::max(0.0, help));
 
  146       InputIt it_help = (last - mid);
 
  148       while (first != it_help)
 
  150         InputIt it_forward = (first - mid);
 
  153         for (j = 0; j < frame_size_; ++j)
 
  155           help += it_forward->getIntensity() * coeffs_[mid * frame_size_ + j];
 
  159         out_it->setPosition(first->getPosition());
 
  160         out_it->setIntensity(std::max(0.0, help));
 
  166       for (i = (mid - 1); i >= 0; --i)
 
  168         InputIt it_forward = (first - (frame_size_ - i - 1));
 
  171         for (j = 0; j < frame_size_; ++j)
 
  173           help += it_forward->getIntensity() * coeffs_[i * frame_size_ + j];
 
  177         out_it->setPosition(first->getPosition());
 
  178         out_it->setIntensity(std::max(0.0, help));
 
  193       filter(spectrum.begin(), spectrum.end(), output.begin());
 
  195       std::swap(spectrum, output);
 
  206       filter(chromatogram.begin(), chromatogram.end(), output.begin());
 
  208       std::swap(chromatogram, output);
 
  218       for (
Size i = 0; i < map.
size(); ++i)
 
  221         setProgress(++progress);
 
  226         setProgress(++progress);
 
  242     void updateMembers_() 
override;
 
  
void filterExperiment(PeakMap &map)
Removed the noise from an MSExperiment containing profile data.
Definition: SavitzkyGolayFilter.h:214
UInt frame_size_
UInt of the filter kernel (number of pre-tabulated coefficients)
Definition: SavitzkyGolayFilter.h:236
In-Memory representation of a mass spectrometry experiment.
Definition: MSExperiment.h:77
size_t Size
Size type e.g. used as variable which can hold result of size()
Definition: Types.h:127
Size size() const
Definition: MSExperiment.h:127
A base class for all classes handling default parameters.
Definition: DefaultParamHandler.h:92
Main OpenMS namespace.
Definition: FeatureDeconvolution.h:46
void filter(MSSpectrum &spectrum)
Removed the noise from an MSSpectrum containing profile data.
Definition: SavitzkyGolayFilter.h:188
Base class for all classes that want to report their progress.
Definition: ProgressLogger.h:54
const std::vector< MSChromatogram > & getChromatograms() const
returns the chromatogram list
unsigned int UInt
Unsigned integer type.
Definition: Types.h:94
void filter(InputIt first, InputIt last, OutputIt d_first)
Definition: SavitzkyGolayFilter.h:115
UInt order_
The order of the smoothing polynomial.
Definition: SavitzkyGolayFilter.h:239
MSChromatogram & getChromatogram(Size id)
returns a single chromatogram
std::vector< double > coeffs_
Coefficients.
Definition: SavitzkyGolayFilter.h:233
The representation of a chromatogram.
Definition: MSChromatogram.h:54
Computes the Savitzky-Golay filter coefficients using QR decomposition.
Definition: SavitzkyGolayFilter.h:101
The representation of a 1D spectrum.
Definition: MSSpectrum.h:67
void filter(MSChromatogram &chromatogram)
Removed the noise from an MSChromatogram.
Definition: SavitzkyGolayFilter.h:201