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

InternalCalibration.h (Maintainer: Alexandra Zerck)

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: Alexandra Zerck $
00025 // --------------------------------------------------------------------------
00026 
00027 
00028 #ifndef OPENMS_FILTERING_CALIBRATION_INTERNALCALIBRATION_H
00029 #define OPENMS_FILTERING_CALIBRATION_INTERNALCALIBRATION_H
00030 
00031 #include <OpenMS/KERNEL/MSExperiment.h>
00032 #include <OpenMS/KERNEL/RawDataPoint1D.h>
00033 #include <OpenMS/TRANSFORMATIONS/RAW2PEAK/PeakPickerCWT.h>
00034 #include <OpenMS/DATASTRUCTURES/DefaultParamHandler.h>
00035 
00036 #include <gsl/gsl_spline.h>
00037 
00038 
00039 
00040 namespace OpenMS
00041 {
00042   
00052   class InternalCalibration 
00053     : public DefaultParamHandler, 
00054       public ProgressLogger
00055   {
00056   public:
00058     typedef RawDataPoint1D RawDataPointType;
00059    
00061     typedef PickedPeak1D PickedPeakType;
00062 
00063 
00065     InternalCalibration();
00066 
00068     ~InternalCalibration(){}
00069 
00071     InternalCalibration(InternalCalibration& obj);
00072 
00074     InternalCalibration& operator=(const InternalCalibration& obj);
00075 
00076 
00082     template<typename InputPeakType>
00083     void calibrate(MSExperiment<InputPeakType>& exp, std::vector<double>& ref_masses,bool peak_data=false);
00084 
00086     inline const DoubleReal getWindowLength() const {return window_length_;}
00088     inline void setWindowLength(const DoubleReal window_length) 
00089     {
00090       window_length_ = window_length;
00091       param_.setValue("window_length",window_length);
00092     }
00093 
00095     inline const MSExperiment<PickedPeakType>& getPeaks() const {return exp_peaks_;}
00097     inline void setPeaks(const MSExperiment<PickedPeakType>& exp_peaks) {exp_peaks_ = exp_peaks;}
00098 
00100     inline const std::vector<std::vector<UInt> >& getMonoisotopicPeaks() const {return monoiso_peaks_;}
00102     inline void setMonoisotopicPeaks(const std::vector<std::vector<UInt> >& monoiso_peaks) {monoiso_peaks_ = monoiso_peaks;}
00103 
00104   protected:
00105 
00106     DoubleReal window_length_;
00107 
00108     MSExperiment<PickedPeakType> exp_peaks_;
00109 
00110     std::vector<std::vector<UInt> > monoiso_peaks_;
00111     
00113     void getMonoisotopicPeaks_();
00114 
00115     // The actual calibration function
00116     template<typename InputPeakType>
00117     void calibrate_(MSExperiment<InputPeakType>& exp, std::vector<double>& ref_masses);
00118     
00119     void updateMembers_();  
00120 
00121   };// class InternalCalibration
00122 
00123 
00124   template<typename InputPeakType>
00125   void InternalCalibration::calibrate(MSExperiment<InputPeakType>& exp, std::vector<double>& ref_masses,bool peak_data)
00126   {
00127 #ifdef DEBUG_CALIBRATION
00128     std::cout.precision(12);
00129 #endif
00130   
00131     if(peak_data)
00132       {
00133         exp_peaks_ = exp;
00134       }
00135     else
00136       {
00137         exp_peaks_.clear();
00138         monoiso_peaks_.clear();
00139         
00140         // sort data
00141         exp.sortSpectra(true);
00142         
00143         // pick peaks (only in a certain distance to the reference masses)
00144         PeakPickerCWT pp;
00145         pp.setParameters(param_.copy("PeakPicker:",true));
00146         typename MSExperiment<InputPeakType>::ConstIterator exp_iter = exp.begin();
00147         typename MSExperiment<InputPeakType>::SpectrumType::ConstIterator spec_iter_l,spec_iter_r; 
00148         for(;exp_iter != exp.end();++exp_iter)
00149         {
00150             MSSpectrum<PickedPeakType> spec;
00151             // pick region around each reference mass
00152             std::vector<double>::iterator vec_iter = ref_masses.begin();
00153             for(;vec_iter != ref_masses.end();++vec_iter)
00154             {
00155                 MSSpectrum<PickedPeakType> tmp_spec;    
00156                 // determine region
00157                 spec_iter_l =  (exp_iter->MZBegin(*vec_iter-window_length_));
00158                 // check borders (avoid )
00159                 spec_iter_r =  (exp_iter->MZBegin(*vec_iter+window_length_));
00160                 if((spec_iter_l >= exp_iter->end()) || (spec_iter_r >= exp_iter->end())) continue;
00161 
00162                 // pick region
00163                 pp.pick(spec_iter_l,spec_iter_r,tmp_spec);
00164                 typename MSSpectrum<PickedPeakType>::Iterator spec_iter = tmp_spec.begin();
00165                 // store
00166                 for(;spec_iter != tmp_spec.end(); ++spec_iter)
00167                 {
00168                     spec.push_back(*spec_iter);
00169                 }
00170             }
00171             if(!spec.empty()) exp_peaks_.push_back(spec);
00172         }
00173       }
00174     
00175     calibrate_(exp,ref_masses);
00176   }
00177 
00178   template<typename InputPeakType>
00179   void InternalCalibration::calibrate_(MSExperiment<InputPeakType>& exp, std::vector<double>& ref_masses)
00180   {
00181     // get monoisotopic peaks
00182     getMonoisotopicPeaks_();
00183 
00184     
00185     size_t num_ref_peaks = ref_masses.size();
00186     std::vector<double> corr_masses,rel_errors;
00187     corr_masses.resize(num_ref_peaks,0.);
00188     rel_errors.resize(num_ref_peaks,0.);
00189     startProgress(0,monoiso_peaks_.size(),"calibrate spectra");    
00190     // for each spectrum
00191     for(size_t spec=0;spec <  monoiso_peaks_.size(); ++spec)
00192       {
00193         UInt corr_peaks=0;
00194         for(size_t peak=0;peak <  monoiso_peaks_[spec].size(); ++peak)
00195           {
00196             for(size_t ref_peak=0; ref_peak < num_ref_peaks;++ref_peak)
00197               {
00198                 if( fabs(exp_peaks_[spec][monoiso_peaks_[spec][peak]].getMZ() - ref_masses[ref_peak]) < 1 )
00199                   {
00200                     corr_masses[ref_peak] = exp_peaks_[spec][monoiso_peaks_[spec][peak]].getMZ();
00201                     ++corr_peaks;
00202                     break;
00203                   }
00204                 
00205               }
00206           }
00207         if(corr_peaks < 2)
00208           {
00209             std::cout << "spec: "<<spec
00210                       << " less than 2 reference masses were detected within a reasonable error range\n";
00211             std::cout << "This spectrum cannot be calibrated!\n";
00212             continue;
00213           }
00214         
00215         double* x = new double[corr_peaks];
00216         double* y = new double[corr_peaks];
00217         UInt p =0;
00218         // determine rel error in ppm for the two reference masses
00219         for(size_t ref_peak=0; ref_peak < num_ref_peaks;++ref_peak)
00220           {
00221             if(corr_masses[ref_peak] != 0.)
00222               {
00223                 rel_errors[ref_peak] = (ref_masses[ref_peak]-corr_masses[ref_peak])/corr_masses[ref_peak] * 1e6;
00224                 x[p] =corr_masses[ref_peak];
00225                 y[p] =rel_errors[ref_peak];
00226                 
00227                 ++p;
00228                 
00229               }
00230           }
00231         
00232         
00233         // linear interpolation
00234         gsl_interp* interp = gsl_interp_alloc(gsl_interp_linear,corr_peaks);
00235         gsl_interp_init(interp, x, y, corr_peaks);
00236         gsl_interp_accel* acc = gsl_interp_accel_alloc();
00237 
00238         
00239         // use interp to internally calibrate the whole spectrum
00240         for(unsigned int peak=0;peak <  exp[spec].size(); ++peak)
00241           {
00242             exp[spec][peak].setMZ(exp[spec][peak].getMZ() + gsl_interp_eval(interp,x,y,
00243                                                                             exp[spec][peak].getMZ(),
00244                                                                             acc)/1e6*exp[spec][peak].getMZ());
00245 #ifdef DEBUG_CALIBRATION
00246             std::cout << exp[spec][peak].getMZ()<< "\t"
00247                       << exp[spec][peak].getMZ() + gsl_interp_eval(interp,x,y,
00248                                                                    exp[spec][peak].getMZ(),
00249                                                                    acc)/1e6*exp[spec][peak].getMZ()
00250                       << std::endl;
00251 #endif
00252           }
00253         delete[] x;
00254         delete[] y;
00255         setProgress(spec);
00256       }// for(size_t spec=0;spec <  monoiso_peaks.size(); ++spec)
00257     endProgress();
00258     
00259   }// calibrate(MSExperiment<InputPeakType> exp, std::vector<Real> ref_masses)
00260 
00261   
00262 } // namespace OpenMS
00263 
00264 #endif // OPENMS_FILTERING_CALIBRATION_INTERNALCALIBRATION_H
00265 

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