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

ContinuousWaveletTransformNumIntegration.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 
00028 #ifndef OPENMS_TRANSFORMATIONS_RAW2PEAK_CONTINUOUSWAVELETTRANSFORMNUMINTEGRATION_H
00029 #define OPENMS_TRANSFORMATIONS_RAW2PEAK_CONTINUOUSWAVELETTRANSFORMNUMINTEGRATION_H
00030 
00031 #include <OpenMS/TRANSFORMATIONS/RAW2PEAK/ContinuousWaveletTransform.h>
00032 
00033 #include<math.h>
00034 
00035 #ifdef DEBUG_PEAK_PICKING
00036 #include<iostream>
00037 #include<fstream>
00038 #endif
00039 
00040 namespace OpenMS
00041 {
00049   class ContinuousWaveletTransformNumIntegration : public ContinuousWaveletTransform
00050   {
00051   public:
00053     typedef ContinuousWaveletTransform::RawDataPointConstIterator RawDataPointConstIterator;
00054 
00055     using ContinuousWaveletTransform::signal_;
00056     using ContinuousWaveletTransform::wavelet_;
00057     using ContinuousWaveletTransform::scale_;
00058     using ContinuousWaveletTransform::spacing_;
00059     using ContinuousWaveletTransform::end_left_padding_;
00060     using ContinuousWaveletTransform::begin_right_padding_;
00061     using ContinuousWaveletTransform::signal_length_;
00062 
00063 
00065     ContinuousWaveletTransformNumIntegration()
00066         : ContinuousWaveletTransform()
00067     {}
00068 
00070     virtual ~ContinuousWaveletTransformNumIntegration() {}
00071 
00081     template < typename InputPeakIterator >
00082     void transform(InputPeakIterator begin_input,
00083                    InputPeakIterator end_input,
00084                    float resolution,
00085                    unsigned int zeros=0)
00086     {
00087 
00088 #ifdef DEBUG_PEAK_PICKING
00089       std::cout << "ContinuousWaveletTransformNumIntegration::transform: start " << begin_input->getMZ() << " until " << (end_input-1)->getMZ() << std::endl;
00090 #endif
00091       if (fabs(resolution-1) < 0.0001)
00092       {
00093         // resolution = 1 corresponds to the cwt at supporting points which have a distance corresponding to the minimal spacing in [begin_input,end_input)
00094         unsigned int n = distance(begin_input,end_input);
00095         signal_length_ = n;
00096 
00097         unsigned int i;
00098 
00099         signal_.clear();
00100         signal_.resize(n);
00101 
00102         // TODO avoid to compute the cwt for the zeros in signal
00103 #ifdef DEBUG_PEAK_PICKING
00104         std::cout << "---------START TRANSFORM---------- \n";
00105 #endif
00106         InputPeakIterator help = begin_input;
00107         for (i=0; i < n; ++i)
00108         {
00109           signal_[i].setMZ(help->getMZ());
00110           signal_[i].setIntensity(integrate_(help,begin_input,end_input));
00111           ++help;
00112         }
00113 #ifdef DEBUG_PEAK_PICKING
00114         std::cout << "---------END TRANSFORM----------" << std::endl;
00115 #endif
00116         // no zeropadding
00117         begin_right_padding_=n;
00118         end_left_padding_=-1;
00119       }
00120       else
00121       {
00122         unsigned int n = (unsigned int) resolution * distance(begin_input, end_input);
00123         double origin  = begin_input->getMZ();
00124         double spacing = ((end_input-1)->getMZ()-origin)/(n-1);
00125         
00126         // zero-padding at the ends?
00127         if(zeros > 0)
00128           {
00129             n += (2*zeros);
00130           }
00131         
00132 
00133         std::vector<double> processed_input(n);
00134         signal_.clear();
00135         signal_.resize(n);
00136 
00137         InputPeakIterator it_help = begin_input;
00138         if(zeros >0)
00139           {
00140             processed_input[0]=it_help->getMZ() - zeros*spacing;
00141             for(unsigned int i = 0; i < zeros; ++i) processed_input[i]=0;
00142           }
00143         else processed_input[0]=it_help->getIntensity();
00144         
00145         double x;
00146         for (unsigned int k=1; k < n-zeros; ++k)
00147         {
00148           x = origin + k*spacing;
00149           // go to the real data point next to x
00150           while (((it_help+1) < end_input) && ((it_help+1)->getMZ() < x))
00151           {
00152             ++it_help;
00153           }
00154           processed_input[k] = getInterpolatedValue_(x,it_help);
00155         }
00156         if(zeros >0)
00157           {
00158             for(unsigned int i = 0; i < zeros; ++i) processed_input[n-zeros+i]=0;
00159           }
00160 
00161         // TODO avoid to compute the cwt for the zeros in signal
00162         for (unsigned int i=0; i < n; ++i)
00163         {
00164           signal_[i].setMZ(origin + i*spacing);
00165           signal_[i].setIntensity(integrate_(processed_input,spacing,i));
00166         }
00167 
00168         if(zeros == 0)
00169           {
00170             begin_right_padding_=n;
00171             end_left_padding_=-1;
00172           }
00173         else
00174           {
00175             begin_right_padding_=n-zeros;
00176             end_left_padding_=zeros-1;
00177           }
00178       }
00179     }
00180 
00181 
00192     virtual void init(double scale, double spacing);
00193 
00194   protected:
00195 
00197     template < typename InputPeakIterator >
00198     double integrate_(InputPeakIterator x, InputPeakIterator first, InputPeakIterator last)
00199     {
00200 #ifdef DEBUG_PEAK_PICKING
00201       std::cout << "integrate_" << std::endl;
00202 #endif
00203 
00204       double v=0.;
00205       int middle = wavelet_.size();
00206 
00207       double start_pos = ((x->getMZ()-(middle*spacing_)) > first->getMZ()) ? (x->getMZ()-(middle*spacing_))
00208                          : first->getMZ();
00209       double end_pos = ((x->getMZ()+(middle*spacing_)) < (last-1)->getMZ()) ? (x->getMZ()+(middle*spacing_))
00210                        : (last-1)->getMZ();
00211 
00212       InputPeakIterator help = x;
00213 
00214 #ifdef DEBUG_PEAK_PICKING
00215       std::cout << "integrate from middle to start_pos "<< help->getMZ() << " until " << start_pos << std::endl;
00216 #endif
00217 
00218       //integrate from middle to start_pos
00219       while ((help != first) && ((help-1)->getMZ() > start_pos))
00220       {
00221         // search for the corresponding datapoint of help in the wavelet (take the left most adjacent point)
00222         double distance = fabs(x->getMZ() - help->getMZ());
00223         int index_w_r = (int)round(distance / spacing_);
00224         double wavelet_right =  wavelet_[index_w_r];
00225 
00226 #ifdef DEBUG_PEAK_PICKING
00227         std::cout << "distance x help " << distance<< std::endl;
00228         std::cout << "distance in wavelet_ " << index_w_r*spacing_ << std::endl;
00229         std::cout << "wavelet_right "  <<  wavelet_right << std::endl;
00230 #endif
00231 
00232         // search for the corresponding datapoint for (help-1) in the wavelet (take the left most adjacent point)
00233         distance = fabs(x->getMZ() - (help-1)->getMZ());
00234         int index_w_l = (int)round(distance / spacing_);
00235         double wavelet_left =  wavelet_[index_w_l];
00236 
00237         // start the interpolation for the true value in the wavelet
00238 
00239 #ifdef DEBUG_PEAK_PICKING
00240         std::cout << " help-1 " << (help-1)->getMZ() << " distance x, help-1" << distance << std::endl;
00241         std::cout << "distance in wavelet_ " << index_w_l*spacing_ << std::endl;
00242         std::cout << "wavelet_ at left " <<   wavelet_left << std::endl;
00243 
00244         std::cout << " intensity " << fabs((help-1)->getMZ()-help->getMZ()) / 2. << " * " << (help-1)->getIntensity() << " * " << wavelet_left <<" + " << (help)->getIntensity()<< "* " << wavelet_right
00245         << std::endl;
00246 #endif
00247 
00248         v+= fabs((help-1)->getMZ()-help->getMZ()) / 2. * ((help-1)->getIntensity()*wavelet_left + help->getIntensity()*wavelet_right);
00249         --help;
00250       }
00251 
00252 
00253       //integrate from middle to end_pos
00254       help = x;
00255 #ifdef DEBUG_PEAK_PICKING
00256       std::cout << "integrate from middle to endpos "<< (help)->getMZ() << " until " << end_pos << std::endl;
00257 #endif
00258       while ((help != (last-1)) && ((help+1)->getMZ() < end_pos))
00259       {
00260         // search for the corresponding datapoint for help in the wavelet (take the left most adjacent point)
00261         double distance = fabs(x->getMZ() - help->getMZ());
00262         int index_w_l = (int)round(distance / spacing_);
00263         double wavelet_left =  wavelet_[index_w_l];
00264 
00265 #ifdef DEBUG_PEAK_PICKING
00266         std::cout << " help " << (help)->getMZ() << " distance x, help" << distance << std::endl;
00267         std::cout << "distance in wavelet_ " << index_w_l*spacing_ << std::endl;
00268         std::cout << "wavelet_ at left " <<   wavelet_left << std::endl;
00269 #endif
00270 
00271         // search for the corresponding datapoint for (help+1) in the wavelet (take the left most adjacent point)
00272         distance = fabs(x->getMZ() - (help+1)->getMZ());
00273         int index_w_r = (int)round(distance / spacing_);
00274         double wavelet_right =  wavelet_[index_w_r];
00275 
00276 #ifdef DEBUG_PEAK_PICKING
00277         std::cout << " help+1 " << (help+1)->getMZ() << " distance x, help+1" << distance << std::endl;
00278         std::cout << "distance in wavelet_ " << index_w_r*spacing_ << std::endl;
00279         std::cout << "wavelet_ at right " <<   wavelet_right << std::endl;
00280 #endif
00281 
00282         v+= fabs(help->getMZ() - (help+1)->getMZ()) / 2. * (help->getIntensity()*wavelet_left + (help+1)->getIntensity()*wavelet_right);
00283         ++help;
00284       }
00285 
00286 
00287 #ifdef DEBUG_PEAK_PICKING
00288       std::cout << "return" << (v/sqrt(scale_)) << std::endl;
00289 #endif
00290       return v / sqrt(scale_);
00291     }
00292 
00294     double integrate_(const std::vector<double>& processed_input, double spacing_data, int index);
00295 
00297     inline double marr_(double x)
00298     {
00299       return (1-x*x)*exp(-x*x/2);
00300     }
00301   };
00302 } //namespace OpenMS
00303 #endif

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