CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC4_patch1/src/RecoTracker/FinalTrackSelectors/src/DuplicateTrackMerger.cc

Go to the documentation of this file.
00001 #include "RecoTracker/FinalTrackSelectors/interface/DuplicateTrackMerger.h"
00002 
00003 #include "FWCore/Framework/interface/Event.h"
00004 #include "TrackingTools/TrajectoryState/interface/FreeTrajectoryState.h"
00005 #include "TrackingTools/PatternTools/interface/TSCPBuilderNoMaterial.h"
00006 #include "TrackingTools/PatternTools/interface/TwoTrackMinimumDistance.h"
00007 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00008 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
00009 
00010 using reco::modules::DuplicateTrackMerger;
00011 
00012 DuplicateTrackMerger::DuplicateTrackMerger(const edm::ParameterSet& iPara) : merger_(iPara)
00013 {
00014   weightFileName_ = "RecoTracker/FinalTrackSelectors/data/DuplicateWeights.xml";
00015   minDeltaR3d_ = -4.0;
00016   minBDTG_ = -0.96;
00017   minpT_ = 0.2;
00018   minP_ = 0.4;
00019   maxDCA_ = 30.0;
00020   maxDPhi_ = 0.30;
00021   maxDLambda_ = 0.30;
00022   maxDdsz_ = 10.0;
00023   maxDdxy_ = 10.0;
00024   maxDQoP_ = 0.25;
00025   if(iPara.exists("minpT"))minpT_ = iPara.getParameter<double>("minpT");
00026   if(iPara.exists("minP"))minP_ = iPara.getParameter<double>("minP");
00027   if(iPara.exists("maxDCA"))maxDCA_ = iPara.getParameter<double>("maxDCA");
00028   if(iPara.exists("maxDPhi"))maxDPhi_ = iPara.getParameter<double>("maxDPhi");
00029   if(iPara.exists("maxDLambda"))maxDLambda_ = iPara.getParameter<double>("maxDLambda");
00030   if(iPara.exists("maxDdsz"))maxDdsz_ = iPara.getParameter<double>("maxDdsz");
00031   if(iPara.exists("maxDdxy"))maxDdxy_ = iPara.getParameter<double>("maxDdxy");
00032   if(iPara.exists("maxDQoP"))maxDQoP_ = iPara.getParameter<double>("maxDQoP");
00033   if(iPara.exists("source"))trackSource_ = iPara.getParameter<edm::InputTag>("source");
00034   if(iPara.exists("minDeltaR3d"))minDeltaR3d_ = iPara.getParameter<double>("minDeltaR3d");
00035   if(iPara.exists("minBDTG"))minBDTG_ = iPara.getParameter<double>("minBDTG");
00036   if(iPara.exists("weightsFile"))weightFileName_ = iPara.getParameter<std::string>("weightsFile");
00037 
00038   produces<std::vector<TrackCandidate> >("candidates");
00039   produces<CandidateToDuplicate>("candidateMap");
00040 
00041   std::string mvaFilePath = edm::FileInPath ( weightFileName_.c_str() ).fullPath();
00042 
00043   tmvaReader_ = new TMVA::Reader("!Color:Silent");
00044   tmvaReader_->AddVariable("ddsz",&tmva_ddsz_);
00045   tmvaReader_->AddVariable("ddxy",&tmva_ddxy_);
00046   tmvaReader_->AddVariable("dphi",&tmva_dphi_);
00047   tmvaReader_->AddVariable("dlambda",&tmva_dlambda_);
00048   tmvaReader_->AddVariable("dqoverp",&tmva_dqoverp_);
00049   tmvaReader_->AddVariable("delta3d_r",&tmva_d3dr_);
00050   tmvaReader_->AddVariable("delta3d_z",&tmva_d3dz_);
00051   tmvaReader_->AddVariable("outer_nMissingInner",&tmva_outer_nMissingInner_);
00052   tmvaReader_->AddVariable("inner_nMissingOuter",&tmva_inner_nMissingOuter_);
00053   tmvaReader_->BookMVA("BDTG",mvaFilePath);
00054 
00055 
00056 }
00057 
00058 DuplicateTrackMerger::~DuplicateTrackMerger()
00059 {
00060 
00061   /* no op */
00062 
00063 }
00064 
00065 void DuplicateTrackMerger::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00066 {
00067 
00068   merger_.init(iSetup);
00069 
00070   edm::Handle<edm::View<reco::Track> >handle;
00071   iEvent.getByLabel(trackSource_,handle);
00072 
00073   iSetup.get<IdealMagneticFieldRecord>().get(magfield_);
00074 
00075   TwoTrackMinimumDistance ttmd;
00076   TSCPBuilderNoMaterial tscpBuilder;
00077   std::auto_ptr<std::vector<TrackCandidate> > out_duplicateCandidates(new std::vector<TrackCandidate>());
00078 
00079   std::auto_ptr<CandidateToDuplicate> out_candidateMap(new CandidateToDuplicate());
00080 
00081   for(int i = 0; i < (int)handle->size(); i++){
00082     const reco::Track *rt1 = &(handle->at(i));
00083     if(rt1->innerMomentum().Rho() < minpT_)continue;
00084     if(rt1->innerMomentum().R() < minP_)continue;
00085     for(int j = i+1; j < (int)handle->size();j++){
00086       const reco::Track *rt2 = &(handle->at(j));
00087       if(rt2->innerMomentum().Rho() < minpT_)continue;
00088       if(rt2->innerMomentum().R() < minP_)continue;
00089       if(rt1->charge() != rt2->charge())continue;
00090       const reco::Track* t1,*t2;
00091       if(rt1->outerPosition().Rho() < rt2->outerPosition().Rho()){
00092         t1 = rt1;
00093         t2 = rt2;
00094       }else{
00095         t1 = rt2;
00096         t2 = rt1;
00097       }
00098       double deltaR3d = sqrt(pow(t1->outerPosition().x()-t2->innerPosition().x(),2) + pow(t1->outerPosition().y()-t2->innerPosition().y(),2) + pow(t1->outerPosition().z()-t2->innerPosition().z(),2));
00099 
00100       if(t1->outerPosition().Rho() > t2->innerPosition().Rho())deltaR3d *= -1.0;
00101       if(deltaR3d < minDeltaR3d_)continue;
00102       
00103       FreeTrajectoryState fts1 = trajectoryStateTransform::outerFreeState(*t1, &*magfield_);
00104       FreeTrajectoryState fts2 = trajectoryStateTransform::innerFreeState(*t2, &*magfield_);
00105       GlobalPoint avgPoint((t1->outerPosition().x()+t2->innerPosition().x())*0.5,(t1->outerPosition().y()+t2->innerPosition().y())*0.5,(t1->outerPosition().z()+t2->innerPosition().z())*0.5);
00106       TrajectoryStateClosestToPoint TSCP1 = tscpBuilder(fts1, avgPoint);
00107       TrajectoryStateClosestToPoint TSCP2 = tscpBuilder(fts2, avgPoint);
00108       if(!TSCP1.isValid())continue;
00109       if(!TSCP2.isValid())continue;
00110 
00111       const FreeTrajectoryState ftsn1 = TSCP1.theState();
00112       const FreeTrajectoryState ftsn2 = TSCP2.theState();
00113  
00114       if ( (ftsn2.position()-ftsn1.position()).mag() > maxDCA_ ) continue;
00115 
00116       double qoverp1 = ftsn1.signedInverseMomentum();
00117       double qoverp2 = ftsn2.signedInverseMomentum();
00118       tmva_dqoverp_ = qoverp1-qoverp2;
00119       if ( fabs(tmva_dqoverp_) > maxDQoP_ ) continue;
00120 
00121       double lambda1 =  M_PI/2 - ftsn1.momentum().theta();
00122       double lambda2 =  M_PI/2 - ftsn2.momentum().theta();
00123       tmva_dlambda_ = lambda1-lambda2;
00124       if ( fabs(tmva_dlambda_) > maxDLambda_ ) continue;
00125 
00126       double phi1 = ftsn1.momentum().phi();
00127       double phi2 = ftsn2.momentum().phi();
00128       tmva_dphi_ = phi1-phi2;
00129       if(fabs(tmva_dphi_) > M_PI) tmva_dphi_ = 2*M_PI - fabs(tmva_dphi_);
00130       if ( fabs(tmva_dphi_) > maxDPhi_ ) continue;
00131 
00132       double dxy1 = (-ftsn1.position().x() * ftsn1.momentum().y() + ftsn1.position().y() * ftsn1.momentum().x())/TSCP1.pt();
00133       double dxy2 = (-ftsn2.position().x() * ftsn2.momentum().y() + ftsn2.position().y() * ftsn2.momentum().x())/TSCP2.pt();
00134       tmva_ddxy_ = dxy1-dxy2;
00135       if ( fabs(tmva_ddxy_) > maxDdxy_ ) continue;
00136 
00137       double dsz1 = ftsn1.position().z() * TSCP1.pt() / TSCP1.momentum().mag() - (ftsn1.position().x() * ftsn1.momentum().y() + ftsn1.position().y() * ftsn1.momentum().x())/TSCP1.pt() * ftsn1.momentum().z()/ftsn1.momentum().mag();
00138       double dsz2 = ftsn2.position().z() * TSCP2.pt() / TSCP2.momentum().mag() - (ftsn2.position().x() * ftsn2.momentum().y() + ftsn2.position().y() * ftsn2.momentum().x())/TSCP2.pt() * ftsn2.momentum().z()/ftsn2.momentum().mag();
00139       tmva_ddsz_ = dsz1-dsz2;
00140       if ( fabs(tmva_ddsz_) > maxDdsz_ ) continue;
00141 
00142       tmva_d3dr_ = avgPoint.perp();
00143       tmva_d3dz_ = avgPoint.z();
00144       tmva_outer_nMissingInner_ = t2->trackerExpectedHitsInner().numberOfLostHits();
00145       tmva_inner_nMissingOuter_ = t1->trackerExpectedHitsOuter().numberOfLostHits();
00146       
00147       double mvaBDTG = tmvaReader_->EvaluateMVA("BDTG");
00148       if(mvaBDTG < minBDTG_)continue;
00149       
00150       
00151       TrackCandidate mergedTrack = merger_.merge(*t1,*t2);
00152       out_duplicateCandidates->push_back(mergedTrack);
00153       std::pair<Track,Track> trackPair(*t1,*t2);
00154       std::pair<TrackCandidate, std::pair<Track,Track> > cp(mergedTrack,trackPair);
00155       out_candidateMap->push_back(cp);
00156     }
00157   }
00158   iEvent.put(out_duplicateCandidates,"candidates");
00159   iEvent.put(out_candidateMap,"candidateMap");
00160 
00161 }