Home  · Classes  · Annotated Classes  · Modules  · Members  · Namespaces  · Related Pages

SavitzkyGolayFilter.h (Maintainer: Eva Lange)

Go to the documentation of this file.
00001 // -*- Mode: C++; tab-width: 2; -*-
00002 // vi: set ts=2:
00003 //
00004 // --------------------------------------------------------------------------
00005 //                   OpenMS Mass Spectrometry Framework
00006 // --------------------------------------------------------------------------
00007 //  Copyright (C) 2003-2008 -- Oliver Kohlbacher, Knut Reinert
00008 //
00009 //  This library is free software; you can redistribute it and/or
00010 //  modify it under the terms of the GNU Lesser General Public
00011 //  License as published by the Free Software Foundation; either
00012 //  version 2.1 of the License, or (at your option) any later version.
00013 //
00014 //  This library is distributed in the hope that it will be useful,
00015 //  but WITHOUT ANY WARRANTY; without even the implied warranty of
00016 //  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00017 //  Lesser General Public License for more details.
00018 //
00019 //  You should have received a copy of the GNU Lesser General Public
00020 //  License along with this library; if not, write to the Free Software
00021 //  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00022 //
00023 // --------------------------------------------------------------------------
00024 // $Maintainer: Eva Lange $
00025 // --------------------------------------------------------------------------
00026 
00027 #ifndef OPENMS_FILTERING_SMOOTHING_SAVITZKYGOLAYFILTER_H
00028 #define OPENMS_FILTERING_SMOOTHING_SAVITZKYGOLAYFILTER_H
00029 
00030 #include <OpenMS/FILTERING/SMOOTHING/SmoothFilter.h>
00031 #include <OpenMS/DATASTRUCTURES/DefaultParamHandler.h>
00032 
00033 #include <gsl/gsl_vector.h>
00034 #include <gsl/gsl_matrix.h>
00035 #include <gsl/gsl_linalg.h>
00036 #include <gsl/gsl_permutation.h>
00037 #include <gsl/gsl_pow_int.h>
00038 
00039 namespace OpenMS
00040 {
00096   class SavitzkyGolayFilter 
00097     : public SmoothFilter, 
00098       public DefaultParamHandler
00099   {
00100     public:
00101       using SmoothFilter::coeffs_;
00102 
00104       SavitzkyGolayFilter();
00105 
00107       virtual ~SavitzkyGolayFilter()
00108       {
00109       }
00110 
00123       template < typename InputPeakIterator, typename OutputPeakContainer  >
00124       void filter(InputPeakIterator first, InputPeakIterator last, OutputPeakContainer& smoothed_data_container) 
00125       {
00126         UInt n = distance(first,last);
00127         if (n <= frame_size_)
00128         {
00129           smoothed_data_container.resize(n);
00130           for (UInt i = 0; i < n; ++i)
00131           {
00132             smoothed_data_container[i] = *first;
00133           }
00134           return;
00135         }
00136 
00137         smoothed_data_container.resize(n);
00138 
00139         int i;
00140         UInt j;
00141         int mid=(frame_size_/2);
00142         double help;
00143 
00144         InputPeakIterator it_forward;
00145         InputPeakIterator it_help;
00146         typename OutputPeakContainer::iterator out_it = smoothed_data_container.begin();
00147 
00148         // compute the transient on
00149         for (i=0; i <= mid; ++i)
00150         {
00151           it_forward=(first-i);
00152           help=0;
00153 
00154           for (j=0; j < frame_size_; ++j)
00155           {
00156             help+=it_forward->getIntensity()*coeffs_[(i+1)*frame_size_-1-j];
00157             ++it_forward;
00158           }
00159 
00160 
00161           out_it->setPosition(first->getPosition());
00162           out_it->setIntensity(std::max(0.0,help));
00163           ++out_it;
00164           ++first;
00165         }
00166 
00167         // compute the steady state output
00168         it_help=(last-mid);
00169         while (first!=it_help)
00170         {
00171           it_forward=(first-mid);
00172           help=0;
00173 
00174           for (j=0; j < frame_size_; ++j)
00175           {
00176             help+=it_forward->getIntensity()*coeffs_[mid*frame_size_+j];
00177             ++it_forward;
00178           }
00179 
00180 
00181           out_it->setPosition(first->getPosition());
00182           out_it->setIntensity(std::max(0.0,help));
00183           ++out_it;
00184           ++first;
00185         }
00186 
00187         // compute the transient off
00188         for (i=(mid-1); i >= 0; --i)
00189         {
00190           it_forward=(first-(frame_size_-i-1));
00191           help=0;
00192 
00193           for (j=0; j < frame_size_; ++j)
00194           {
00195             help+=it_forward->getIntensity()*coeffs_[i*frame_size_+j];
00196             ++it_forward;
00197           }
00198 
00199           out_it->setPosition(first->getPosition());
00200           out_it->setIntensity(std::max(0.0,help));
00201           ++out_it;
00202           ++first;
00203         }
00204       }
00205 
00206 
00219       template <typename InputPeakContainer, typename OutputPeakContainer >
00220       void filter(const InputPeakContainer& input_peak_container, OutputPeakContainer& baseline_filtered_container)
00221       {
00222         // copy the spectrum settings
00223         static_cast<SpectrumSettings&>(baseline_filtered_container) = input_peak_container;
00224         
00225         filter(input_peak_container.begin(), input_peak_container.end(), baseline_filtered_container);
00226       }
00227 
00239       template <typename InputSpectrumIterator, typename OutputPeakType >
00240       void filterExperiment(InputSpectrumIterator first, InputSpectrumIterator last, MSExperiment<OutputPeakType>& ms_exp_filtered)
00241       {
00242         UInt n = distance(first,last);
00243         ms_exp_filtered.reserve(n);
00244         startProgress(0,n,"smoothing data");
00245         
00246         // smooth each scan
00247         for (UInt i = 0; i < n; ++i)
00248         {
00249           MSSpectrum< OutputPeakType > spectrum;
00250           InputSpectrumIterator input_it(first+i);
00251 
00252           // smooth scan i
00253           filter(*input_it,spectrum);
00254 
00255           // copy the spectrum settings
00256           static_cast<SpectrumSettings&>(spectrum) = *input_it;
00257           spectrum.setType(SpectrumSettings::RAWDATA);
00258 
00259           // copy the spectrum information
00260           spectrum.getPrecursorPeak() = input_it->getPrecursorPeak();
00261           spectrum.setRT(input_it->getRT());
00262           spectrum.setMSLevel(input_it->getMSLevel());
00263           spectrum.getName() = input_it->getName();
00264 
00265           ms_exp_filtered.push_back(spectrum);
00266         }
00267         endProgress();
00268       }
00269 
00270 
00280       template <typename InputPeakType, typename OutputPeakType >
00281       void filterExperiment(const MSExperiment< InputPeakType >& ms_exp_raw, MSExperiment<OutputPeakType>& ms_exp_filtered)
00282       {
00283         // copy the experimental settings
00284         static_cast<ExperimentalSettings&>(ms_exp_filtered) = ms_exp_raw;
00285 
00286         filterExperiment(ms_exp_raw.begin(), ms_exp_raw.end(), ms_exp_filtered);
00287       }
00288     
00289     protected:
00291       UInt frame_size_;
00293       UInt order_;
00294       // Docu in base class
00295       virtual void updateMembers_();
00296       
00297   };
00298 
00299 } // namespace OpenMS
00300 #endif // OPENMS_FILTERING_SMOOTHING_SAVITZKYGOLAYFILTER_H
00301 

Generated Tue Apr 1 15:36:37 2008 -- using doxygen 1.5.4 OpenMS / TOPP 1.1