00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028 #ifndef OPENMS_TRANSFORMATIONS_RAW2PEAK_TWODOPTIMIZATION_H
00029 #define OPENMS_TRANSFORMATIONS_RAW2PEAK_TWODOPTIMIZATION_H
00030
00031
00032 #undef DEBUG_2D
00033
00034 #ifdef DEBUG_2D
00035 #include<iostream>
00036 #include<fstream>
00037 #endif
00038
00039 #include <vector>
00040 #include <utility>
00041 #include <cmath>
00042 #include <set>
00043
00044 #include <OpenMS/TRANSFORMATIONS/RAW2PEAK/PeakShape.h>
00045 #include <OpenMS/KERNEL/PickedPeak1D.h>
00046 #include <OpenMS/KERNEL/MSExperiment.h>
00047 #include <OpenMS/KERNEL/MSSpectrum.h>
00048 #include <OpenMS/KERNEL/DPeak.h>
00049 #include <OpenMS/CONCEPT/Exception.h>
00050 #include <OpenMS/DATASTRUCTURES/Param.h>
00051 #include <OpenMS/DATASTRUCTURES/IsotopeCluster.h>
00052 #include <OpenMS/DATASTRUCTURES/DefaultParamHandler.h>
00053
00054 #ifndef OPENMS_SYSTEM_STOPWATCH_H
00055 # include <OpenMS/SYSTEM/StopWatch.h>
00056 #endif
00057
00058
00059 #include <gsl/gsl_vector.h>
00060 #include <gsl/gsl_multifit_nlin.h>
00061 #include <gsl/gsl_blas.h>
00062 #include <OpenMS/TRANSFORMATIONS/RAW2PEAK/OptimizePeakDeconvolution.h>
00063 #include <OpenMS/TRANSFORMATIONS/RAW2PEAK/OptimizePick.h>
00064 #include <OpenMS/TRANSFORMATIONS/RAW2PEAK/PeakShape.h>
00065
00066 namespace OpenMS
00067 {
00068
00069 typedef std::pair<unsigned int,unsigned int> Idx ;
00070 typedef std::set<Idx> IndexSet;
00071
00072 namespace OptimizationFunctions
00073 {
00074
00076 typedef RawDataPoint1D RawDataPointType;
00078 typedef MSExperiment<PickedPeak1D > ExperimentPickedType;
00079 extern std::vector<std::pair<int,int> > signal2D;
00080 extern std::multimap<double,IsotopeCluster>::iterator iso_map_iter;
00081 extern unsigned int total_nr_peaks;
00082 extern std::map<int, std::vector<ExperimentPickedType::SpectrumType::Iterator > > matching_peaks;
00083 extern MSExperiment<PickedPeak1D>::Iterator picked_peaks_iter;
00084 extern MSExperiment<RawDataPointType>::ConstIterator raw_data_first;
00085
00086
00091
00092 int residual2D(const gsl_vector* x, void* params , gsl_vector* f);
00094 int jacobian2D(const gsl_vector* x, void* params, gsl_matrix* J);
00096 int evaluate2D(const gsl_vector* x, void* params, gsl_vector* f, gsl_matrix* J);
00097
00098 }
00099
00114 class TwoDOptimization : public DefaultParamHandler
00115 {
00116 public:
00117
00119 struct IndexLess
00120 : public std::binary_function <Idx,Idx, bool>
00121 {
00122 inline bool operator () (const Idx& a, const Idx& b) const
00123 {
00124 return (a.first < b.first);
00125 }
00126 };
00127
00129 TwoDOptimization();
00130
00132 TwoDOptimization(const TwoDOptimization& opt);
00133
00135 virtual ~TwoDOptimization(){}
00136
00138 TwoDOptimization& operator=(const TwoDOptimization& opt);
00139
00140
00142 inline DoubleReal getMZTolerance() const {return tolerance_mz_;}
00144 inline void setMZTolerance(double tolerance_mz)
00145 {
00146 tolerance_mz_ = tolerance_mz;
00147 param_.setValue("2d:tolerance_mz",tolerance_mz);
00148 }
00149
00151 inline DoubleReal getMaxPeakDistance() const {return max_peak_distance_;}
00153 inline void setMaxPeakDistance(double max_peak_distance)
00154 {
00155 max_peak_distance_ = max_peak_distance;
00156 param_.setValue("2d:max_peak_distance",max_peak_distance);
00157 }
00158
00160 inline DoubleReal getMaxAbsError() const {return eps_abs_;}
00162 inline void setMaxAbsError(double eps_abs)
00163 {
00164 eps_abs_ = eps_abs;
00165 param_.setValue("delta_abs_error",eps_abs);
00166 }
00167
00169 inline DoubleReal getMaxRelError() const {return eps_rel_;}
00171 inline void setMaxRelError(double eps_rel)
00172 {
00173 eps_rel_ = eps_rel;
00174 param_.setValue("delta_rel_error",eps_rel);
00175 }
00176
00178 inline Int getMaxIterations() const {return max_iteration_;}
00180 inline void setMaxIterations(int max_iteration)
00181 {
00182 max_iteration_ = max_iteration;
00183 param_.setValue("iterations",max_iteration);
00184 }
00185
00187 inline const OptimizationFunctions::PenaltyFactorsIntensity& getPenalties() const {return penalties_;}
00189 inline void setPenalties(OptimizationFunctions::PenaltyFactorsIntensity& penalties)
00190 {
00191 penalties_ = penalties;
00192 param_.setValue("penalties:position",penalties.pos);
00193 param_.setValue("penalties:height",penalties.height);
00194 param_.setValue("penalties:left_width",penalties.lWidth);
00195 param_.setValue("penalties:right_width",penalties.rWidth);
00196 }
00197
00198
00199
00201 template <typename InputSpectrumIterator,typename OutputPeakType>
00202 void optimize(InputSpectrumIterator& first,
00203 InputSpectrumIterator& last,
00204 MSExperiment< OutputPeakType >& ms_exp,bool real2D=true);
00205
00206
00207 protected:
00208
00210 std::multimap<double,IsotopeCluster> iso_map_;
00211
00213 std::multimap<double,IsotopeCluster>::const_iterator curr_region_;
00214
00216 double max_peak_distance_;
00217
00219 double tolerance_mz_;
00220
00222 std::map<int, std::vector<MSExperiment<PickedPeak1D >::SpectrumType::Iterator > > matching_peaks_;
00223
00225 double eps_abs_;
00226
00228 double eps_rel_;
00229
00231 int max_iteration_;
00232
00234 bool real_2D_;
00235
00236
00238 OptimizationFunctions::PenaltyFactorsIntensity penalties_;
00239
00244 std::vector<double>::iterator searchInScan_(std::vector<double>::iterator scan_begin,
00245 std::vector<double>::iterator scan_end ,
00246 double current_mz);
00247
00249 template <typename InputSpectrumIterator,typename OutputPeakType>
00250 void optimizeRegions_(InputSpectrumIterator& first,
00251 InputSpectrumIterator& last,
00252 MSExperiment<OutputPeakType>& ms_exp);
00253
00255 template <typename InputSpectrumIterator,typename OutputPeakType>
00256 void optimizeRegionsScanwise_(InputSpectrumIterator& first,
00257 InputSpectrumIterator& last,
00258 MSExperiment<OutputPeakType>& ms_exp);
00259
00260
00262 template<typename InputSpectrumIterator,typename OutputPeakType>
00263 void getRegionEndpoints_(MSExperiment<OutputPeakType>& exp,
00264 InputSpectrumIterator& first,
00265 InputSpectrumIterator& last,
00266 unsigned int iso_map_idx,
00267 double noise_level);
00268
00270 void findMatchingPeaks_(std::multimap<double, IsotopeCluster>::iterator& it,
00271 MSExperiment<PickedPeak1D>& ms_exp);
00272
00274
00276 void updateMembers_();
00277 };
00278
00279
00280 template <typename InputSpectrumIterator,typename OutputPeakType>
00281 void TwoDOptimization::optimize(InputSpectrumIterator& first,
00282 InputSpectrumIterator& last,
00283 MSExperiment< OutputPeakType >& ms_exp,bool real2D)
00284 {
00285 real_2D_ = real2D;
00286 typedef typename InputSpectrumIterator::value_type InputSpectrumType;
00287 typedef typename InputSpectrumType::value_type RawDataPointType;
00288 typedef MSSpectrum<RawDataPointType> SpectrumType;
00289
00290 typename MSExperiment<OutputPeakType>::Iterator ms_exp_it = ms_exp.begin();
00291 typename MSExperiment<OutputPeakType>::Iterator ms_exp_it_end = ms_exp.end();
00292 if(ms_exp.size() == 0)
00293 {
00294 std::cout << "empty experiment"<<std::endl;
00295 return;
00296 }
00297
00298 std::vector<double> iso_last_scan;
00299 std::vector<double> iso_curr_scan;
00300 std::vector<std::multimap<double,IsotopeCluster>::iterator> clusters_last_scan;
00301 std::vector<std::multimap<double,IsotopeCluster>::iterator> clusters_curr_scan;
00302 std::multimap<double,IsotopeCluster>::iterator cluster_iter;
00303 double current_rt=ms_exp_it->getRT(),last_rt = 0;
00304
00305
00306 max_peak_distance_ = param_.getValue("2d:max_peak_distance");
00307 double tolerance_mz = param_.getValue("2d:tolerance_mz");
00308
00309 UInt current_charge = 0;
00310 double mz_in_hash = 0;
00311
00312
00313 for (unsigned int curr_scan =0; ms_exp_it+curr_scan != ms_exp_it_end;++curr_scan)
00314 {
00315 unsigned int nr_peaks_in_scan = (ms_exp_it +curr_scan)->size();
00316
00317 current_rt = (ms_exp_it+curr_scan)->getRT();
00318 typename MSExperiment<OutputPeakType>::SpectrumType::Iterator peak_it = (ms_exp_it+curr_scan)->begin();
00319 typename MSExperiment<OutputPeakType>::SpectrumType::Iterator peak_it_last = (ms_exp_it+curr_scan)->end();
00320
00321
00322 iso_last_scan = iso_curr_scan;
00323 iso_curr_scan.clear();
00324 clusters_last_scan = clusters_curr_scan;
00325 clusters_curr_scan.clear();
00326
00327 #ifdef DEBUG_2D
00328 std::cout << "Next scan with rt: " << current_rt << std::endl;
00329 std::cout << "Next scan, rt = "<<current_rt<<" last_rt: "<<last_rt<< std::endl;
00330 std::cout << "---------------------------------------------------------------------------" << std::endl;
00331 #endif
00332 MSSpectrum<RawDataPointType> s;
00333 s.setRT(current_rt);
00334
00335 if ( last_rt == 0||
00336 ( (lower_bound(first,last,s, typename SpectrumType::RTLess())-1)->getRT() == last_rt))
00337 {
00338
00339
00340 for(unsigned int curr_peak=0; peak_it+curr_peak < peak_it_last-1;++curr_peak)
00341 {
00342
00343
00344 double curr_mz = (peak_it+curr_peak)->getMZ();
00345 double dist2nextpeak = (peak_it+curr_peak+1)->getMZ() - curr_mz;
00346
00347 if (dist2nextpeak <= max_peak_distance_)
00348 {
00349 #ifdef DEBUG_2D
00350 std::cout << "Isotopic pattern found ! " << std::endl;
00351 std::cout << "We are at: " << (peak_it+curr_peak)->getMZ() << " " << curr_mz << std::endl;
00352 #endif
00353 if (iso_last_scan.size() > 0)
00354 {
00355
00356 std::vector<double>::iterator it =
00357 searchInScan_(iso_last_scan.begin(),iso_last_scan.end(),curr_mz);
00358
00359 double delta_mz = fabs(*it - curr_mz);
00360 std::vector<double>::iterator itneu = iso_last_scan.begin();
00361
00362 if ( delta_mz > tolerance_mz)
00363 {
00364 mz_in_hash = curr_mz;
00365
00366
00367
00368
00369
00370 IsotopeCluster new_cluster;
00371 new_cluster.peaks_.charge_ = current_charge;
00372 new_cluster.scans_.push_back( curr_scan );
00373 cluster_iter = iso_map_.insert(std::pair<double,IsotopeCluster>(mz_in_hash,new_cluster));
00374
00375 }
00376 else
00377 {
00378
00379
00380
00381 cluster_iter = clusters_last_scan[distance(iso_last_scan.begin(),it)];
00382
00383
00384 if(find(cluster_iter->second.scans_.begin(),cluster_iter->second.scans_.end(),curr_scan)
00385 == cluster_iter->second.scans_.end())
00386 {
00387 cluster_iter->second.scans_.push_back( curr_scan );
00388 }
00389
00390
00391
00392
00393
00394 }
00395
00396 }
00397 else
00398 {
00399
00400
00401
00402
00403
00404 mz_in_hash = curr_mz;
00405
00406
00407 IsotopeCluster new_cluster;
00408 new_cluster.peaks_.charge_ = current_charge;
00409 new_cluster.scans_.push_back( curr_scan );
00410 cluster_iter = iso_map_.insert(std::pair<double,IsotopeCluster>(mz_in_hash,new_cluster));
00411
00412 }
00413
00414
00415
00416
00417
00418
00419
00420 cluster_iter->second.peaks_.insert(std::pair<UInt,UInt>(curr_scan,curr_peak));
00421
00422 iso_curr_scan.push_back( mz_in_hash );
00423 clusters_curr_scan.push_back(cluster_iter);
00424 ++curr_peak;
00425
00426 cluster_iter->second.peaks_.insert(std::pair<UInt,UInt>(curr_scan,curr_peak));
00427 iso_curr_scan.push_back((peak_it+curr_peak)->getMZ());
00428 clusters_curr_scan.push_back(cluster_iter);
00429
00430
00431 if ( (curr_peak+1) >= nr_peaks_in_scan ) break;
00432 dist2nextpeak = (peak_it+curr_peak+1)->getMZ() - (peak_it+curr_peak)->getMZ();
00433
00434
00435
00436 while (dist2nextpeak <= max_peak_distance_
00437 && curr_peak < (nr_peaks_in_scan-1) )
00438 {
00439 cluster_iter->second.peaks_.insert(std::pair<UInt,UInt>(curr_scan,curr_peak+1));
00440 iso_curr_scan.push_back((peak_it+curr_peak+1)->getMZ());
00441 clusters_curr_scan.push_back(cluster_iter);
00442
00443 ++curr_peak;
00444 if(curr_peak >= nr_peaks_in_scan-1) break;
00445 dist2nextpeak = (peak_it+curr_peak+1)->getMZ() - (peak_it+curr_peak)->getMZ();
00446
00447
00448 }
00449
00450
00451
00452 }
00453 else
00454 {
00455 if (iso_last_scan.size() > 0)
00456 {
00457
00458 std::vector<double>::iterator it =
00459 searchInScan_(iso_last_scan.begin(),iso_last_scan.end(),curr_mz);
00460
00461 double delta_mz = fabs(*it - curr_mz);
00462 std::vector<double>::iterator itneu = iso_last_scan.begin();
00463
00464 if ( delta_mz > tolerance_mz)
00465 {
00466 mz_in_hash = curr_mz;
00467
00468
00469
00470
00471
00472 IsotopeCluster new_cluster;
00473 new_cluster.peaks_.charge_ = current_charge;
00474 new_cluster.scans_.push_back( curr_scan );
00475 cluster_iter = iso_map_.insert(std::pair<double,IsotopeCluster>(mz_in_hash,new_cluster));
00476
00477 }
00478 else
00479 {
00480
00481
00482
00483 cluster_iter = clusters_last_scan[distance(iso_last_scan.begin(),it)];
00484
00485
00486 if(find(cluster_iter->second.scans_.begin(),cluster_iter->second.scans_.end(),curr_scan)
00487 == cluster_iter->second.scans_.end())
00488 {
00489 cluster_iter->second.scans_.push_back( curr_scan );
00490 }
00491
00492
00493
00494
00495
00496 }
00497
00498 }
00499 else
00500 {
00501
00502
00503
00504
00505
00506 mz_in_hash = curr_mz;
00507
00508
00509 IsotopeCluster new_cluster;
00510 new_cluster.peaks_.charge_ = current_charge;
00511 new_cluster.scans_.push_back( curr_scan );
00512 cluster_iter = iso_map_.insert(std::pair<double,IsotopeCluster>(mz_in_hash,new_cluster));
00513
00514 }
00515
00516
00517
00518
00519
00520
00521
00522 cluster_iter->second.peaks_.insert(std::pair<UInt,UInt>(curr_scan,curr_peak));
00523
00524 iso_curr_scan.push_back( mz_in_hash );
00525 clusters_curr_scan.push_back(cluster_iter);
00526
00527
00528 }
00529
00530 current_charge = 0;
00531 }
00532 }
00533 last_rt = current_rt;
00534 }
00535 curr_region_ = iso_map_.begin();
00536 #ifdef DEBUG_2D
00537 std::cout << iso_map_.size() << " isotopic clusters were found ! " << std::endl;
00538 #endif
00539
00540 if(real_2D_) optimizeRegions_(first,last,ms_exp);
00541 else optimizeRegionsScanwise_(first,last,ms_exp);
00542 }
00543
00544
00545
00546
00547 template < >
00548 void TwoDOptimization::optimizeRegions_<MSExperiment<RawDataPoint1D >::const_iterator, PickedPeak1D >
00549 (MSExperiment<RawDataPoint1D >::const_iterator& first,
00550 MSExperiment<RawDataPoint1D >::const_iterator& last,
00551 MSExperiment<PickedPeak1D >& ms_exp);
00552
00553 template < >
00554 void TwoDOptimization::optimizeRegionsScanwise_<MSExperiment<RawDataPoint1D >::const_iterator, PickedPeak1D >
00555 (MSExperiment<RawDataPoint1D >::const_iterator& first,
00556 MSExperiment<RawDataPoint1D >::const_iterator& last,
00557 MSExperiment<PickedPeak1D >& ms_exp);
00558
00559 template <typename InputSpectrumIterator,typename OutputPeakType>
00560 void TwoDOptimization::optimizeRegions_(InputSpectrumIterator& ,
00561 InputSpectrumIterator& ,
00562 MSExperiment< OutputPeakType >& )
00563 {
00564 throw Exception::IllegalArgument(__FILE__,
00565 __LINE__,
00566 __PRETTY_FUNCTION__,
00567 "wrong input peak type, must be PickedPeak1D,");
00568
00569
00570 }
00571
00572 template <typename InputSpectrumIterator,typename OutputPeakType>
00573 void TwoDOptimization::optimizeRegionsScanwise_(InputSpectrumIterator& ,
00574 InputSpectrumIterator& ,
00575 MSExperiment< OutputPeakType >& )
00576 {
00577 throw Exception::IllegalArgument(__FILE__,
00578 __LINE__,
00579 __PRETTY_FUNCTION__,
00580 "wrong input peak type, must be PickedPeak1D,");
00581
00582
00583 }
00584
00585
00586 template<typename InputSpectrumIterator,typename OutputPeakType>
00587 void TwoDOptimization::getRegionEndpoints_(MSExperiment<OutputPeakType>& exp,
00588 InputSpectrumIterator& first,
00589 InputSpectrumIterator& last,
00590 unsigned int iso_map_idx,
00591 double noise_level)
00592 {
00593 OptimizationFunctions::signal2D.clear();
00594 typedef typename InputSpectrumIterator::value_type InputExperimentType;
00595 typedef typename InputExperimentType::value_type InputPeakType;
00596 typedef std::multimap<double,IsotopeCluster> MapType;
00597
00598 double rt, first_peak_mz, last_peak_mz;
00599
00600
00601
00602 typename MSExperiment<InputPeakType>::SpectrumType spec;
00603 InputPeakType peak;
00604
00605 MapType::iterator iso_map_iter = iso_map_.begin();
00606 for(unsigned int i=0;i<iso_map_idx;++i) ++iso_map_iter;
00607
00608 #ifdef DEBUG2D
00609 std::cout << "rt begin: "<<exp[iso_map_iter->second.scans_[0]].getRT()
00610 << "\trt end: "<<exp[iso_map_iter->second.scans_[iso_map_iter->second.scans_.size()-1]].getRT()
00611 << " \t"<<iso_map_iter->second.scans_.size()<<" scans"
00612 << std::endl;
00613 #endif
00614
00615
00616 for(unsigned int i=0; i< iso_map_iter->second.scans_.size();++i)
00617 {
00618 typename MSExperiment<OutputPeakType>::iterator exp_it;
00619
00620
00621 rt = exp[iso_map_iter->second.scans_[i]].getRT();
00622 spec.setRT(rt);
00623 InputSpectrumIterator iter = lower_bound(first, last, spec, typename MSSpectrum<InputPeakType>::RTLess());
00624
00625 exp_it = exp.RTBegin(rt);
00626 #ifdef DEBUG2D
00627 std::cout << exp_it->getRT() << " vs "<< iter->getRT()<<std::endl;
00628 #endif
00629
00630 IndexSet::const_iterator j=(iso_map_iter->second.peaks_.begin());
00631
00632 Idx pair;
00633 pair.first = iso_map_iter->second.peaks_.begin()->first + i;
00634 IndexSet::const_iterator set_iter = lower_bound(iso_map_iter->second.peaks_.begin(),
00635 iso_map_iter->second.peaks_.end(),
00636 pair,IndexLess());
00637
00638
00639 first_peak_mz = (exp_it->begin() + set_iter->second)->getMZ() - 1;
00640
00641
00642 ++pair.first;
00643 IndexSet::const_iterator set_iter2 = lower_bound(iso_map_iter->second.peaks_.begin(),
00644 iso_map_iter->second.peaks_.end(),
00645 pair,IndexLess());
00646 --set_iter2;
00647 last_peak_mz = (exp_it->begin() + set_iter2->second)->getMZ() + 1;
00648
00649
00650 peak.setPosition(first_peak_mz);
00651 typename MSExperiment<InputPeakType>::SpectrumType::const_iterator raw_data_iter
00652 = lower_bound(iter->begin(), iter->end(), peak, typename InputPeakType::PositionLess());
00653 if(raw_data_iter != iter->begin())
00654 {
00655 --raw_data_iter;
00656 }
00657 double intensity = raw_data_iter->getIntensity();
00658
00659 while(raw_data_iter != iter->begin() && (raw_data_iter-1)->getIntensity() < intensity &&
00660 (raw_data_iter-1)->getIntensity() > noise_level)
00661 {
00662 --raw_data_iter;
00663 intensity = raw_data_iter->getIntensity();
00664 }
00665 ++raw_data_iter;
00666 Idx left,right;
00667 left.first = distance(first,iter);
00668 left.second = distance(iter->begin(),raw_data_iter);
00669 #ifdef DEBUG2D
00670 std::cout << "left: "<<iter->getRT()<<"\t"<<raw_data_iter->getMZ()<<std::endl;
00671 #endif
00672
00673 peak.setPosition(last_peak_mz + 1);
00674 raw_data_iter
00675 = upper_bound(iter->begin(), iter->end(), peak, typename InputPeakType::PositionLess());
00676 if(raw_data_iter == iter->end()) --raw_data_iter;
00677 intensity = raw_data_iter->getIntensity();
00678
00679 while(raw_data_iter+1 != iter->end() && (raw_data_iter+1)->getIntensity() < intensity)
00680 {
00681 ++raw_data_iter;
00682 intensity = raw_data_iter->getIntensity();
00683 if( (raw_data_iter+1!= iter->end()) && (raw_data_iter+1)->getIntensity() > noise_level ) break;
00684 }
00685 right.first = left.first;
00686 right.second = distance(iter->begin(),raw_data_iter);
00687 #ifdef DEBUG2D
00688 std::cout << "rightt: "<<iter->getRT()<<"\t"<<raw_data_iter->getMZ()<<std::endl;
00689 #endif
00690
00691 OptimizationFunctions::signal2D.push_back(left);
00692 OptimizationFunctions::signal2D.push_back(right);
00693 }
00694 #ifdef DEBUG2D
00695
00696 std::cout << first_peak_mz <<"\t"<<last_peak_mz<<std::endl;
00697 #endif
00698 }
00699
00700
00701
00702
00703 }
00704
00705 #endif //OPENMS_TRANSFORMATIONS_RAW2PEAK_TWODOPTIMIZATION_H