88 template<
class InputIt,
class OutputIt>
89 void filter(InputIt first, InputIt last, OutputIt d_first)
91 size_t n = std::distance(first, last);
93 if (frame_size_ > n) {
return; }
97 int mid = (frame_size_ / 2);
101 OutputIt out_it = d_first;
103 for (i = 0; i <= mid; ++i)
105 InputIt it_forward = (first - i);
107 for (j = 0; j < frame_size_; ++j)
109 help += it_forward->getIntensity() * coeffs_[(i + 1) * frame_size_ - 1 - j];
113 out_it->setPosition(first->getPosition());
114 out_it->setIntensity(std::max(0.0, help));
120 InputIt it_help = (last - mid);
122 while (first != it_help)
124 InputIt it_forward = (first - mid);
127 for (j = 0; j < frame_size_; ++j)
129 help += it_forward->getIntensity() * coeffs_[mid * frame_size_ + j];
133 out_it->setPosition(first->getPosition());
134 out_it->setIntensity(std::max(0.0, help));
140 for (i = (mid - 1); i >= 0; --i)
142 InputIt it_forward = (first - (frame_size_ - i - 1));
145 for (j = 0; j < frame_size_; ++j)
147 help += it_forward->getIntensity() * coeffs_[i * frame_size_ + j];
151 out_it->setPosition(first->getPosition());
152 out_it->setIntensity(std::max(0.0, help));
167 filter(spectrum.begin(), spectrum.end(), output.begin());
169 std::swap(spectrum, output);
180 filter(chromatogram.begin(), chromatogram.end(), output.begin());
182 std::swap(chromatogram, output);
192 for (
Size i = 0; i < map.
size(); ++i)
195 setProgress(++progress);
200 setProgress(++progress);
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:46
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:121
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