CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch12/src/RecoParticleFlow/PFTracking/src/ConvBremPFTrackFinder.cc

Go to the documentation of this file.
00001 #include "RecoParticleFlow/PFTracking/interface/ConvBremPFTrackFinder.h"
00002 #include "FWCore/Framework/interface/ESHandle.h"
00003 
00004 
00005 #include "DataFormats/ParticleFlowReco/interface/PFRecHit.h"
00006 #include "RecoParticleFlow/PFProducer/interface/Utils.h"
00007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00008 #include "DataFormats/EgammaReco/interface/ElectronSeed.h"
00009 #include "DataFormats/EgammaReco/interface/ElectronSeedFwd.h"
00010 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h" 
00011 #include "DataFormats/GsfTrackReco/interface/GsfTrackFwd.h" 
00012 #include "TrackingTools/IPTools/interface/IPTools.h"
00013 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
00014 #include "RecoParticleFlow/PFClusterTools/interface/PFEnergyCalibration.h"
00015 #include "TMath.h"
00016 #include "RecoParticleFlow/PFClusterTools/interface/LinkByRecHit.h"
00017 
00018 using namespace edm;
00019 using namespace std;
00020 using namespace reco;
00021 
00022 ConvBremPFTrackFinder::ConvBremPFTrackFinder(const TransientTrackBuilder& builder,
00023                                              double mvaBremConvCut,
00024                                              string mvaWeightFileConvBrem):
00025   builder_(builder),
00026   mvaBremConvCut_(mvaBremConvCut),
00027   mvaWeightFileConvBrem_(mvaWeightFileConvBrem)
00028 {
00029   tmvaReader_ = new TMVA::Reader("!Color:Silent");
00030   tmvaReader_->AddVariable("secR",&secR);
00031   tmvaReader_->AddVariable("sTIP",&sTIP);
00032   tmvaReader_->AddVariable("nHITS1",&nHITS1);
00033   tmvaReader_->AddVariable("secPin",&secPin);
00034   tmvaReader_->AddVariable("Epout",&Epout);
00035   tmvaReader_->AddVariable("detaBremKF",&detaBremKF);
00036   tmvaReader_->AddVariable("ptRatioGsfKF",&ptRatioGsfKF);
00037   tmvaReader_->BookMVA("BDT",mvaWeightFileConvBrem.c_str());
00038 }
00039 ConvBremPFTrackFinder::~ConvBremPFTrackFinder(){delete tmvaReader_;}
00040 
00041 void
00042 ConvBremPFTrackFinder::runConvBremFinder(const Handle<PFRecTrackCollection>& thePfRecTrackCol,
00043                                          const Handle<VertexCollection>& primaryVertex,
00044                                          const edm::Handle<reco::PFDisplacedTrackerVertexCollection>& pfNuclears,
00045                                          const edm::Handle<reco::PFConversionCollection >& pfConversions,
00046                                          const edm::Handle<reco::PFV0Collection >& pfV0,
00047                                          bool useNuclear,
00048                                          bool useConversions,
00049                                          bool useV0,
00050                                          const reco::PFClusterCollection & theEClus,
00051                                          reco::GsfPFRecTrack gsfpfrectk)
00052 {
00053   
00054 
00055   found_ = false;
00056   bool debug = false;
00057   bool debugRef = false;
00058   
00059   if(debug)
00060     cout << "runConvBremFinder:: Entering " << endl;
00061   
00062   
00063   
00064   reco::GsfTrackRef refGsf =  gsfpfrectk.gsfTrackRef();
00065   reco::PFRecTrackRef pfTrackRef = gsfpfrectk.kfPFRecTrackRef();
00066   vector<PFBrem> primPFBrem = gsfpfrectk.PFRecBrem();
00067  
00068   
00069   const PFRecTrackCollection& PfRTkColl = *(thePfRecTrackCol.product());
00070   reco::PFRecTrackCollection::const_iterator pft=PfRTkColl.begin();
00071   reco::PFRecTrackCollection::const_iterator pftend=PfRTkColl.end();
00072   PFEnergyCalibration pfcalib_;
00073 
00074 
00075 
00076 
00077   vector<PFRecTrackRef> AllPFRecTracks;
00078   AllPFRecTracks.clear();
00079   unsigned int ipft = 0;
00080   
00081 
00082   for(;pft!=pftend;++pft,ipft++){
00083     // do not consider the kf track already associated to the seed
00084     if(pfTrackRef.isNonnull())
00085       if(pfTrackRef->trackRef() == pft->trackRef()) continue;
00086     
00087     PFRecTrackRef pfRecTrRef(thePfRecTrackCol,ipft);  
00088     TrackRef trackRef = pfRecTrRef->trackRef();
00089     reco::TrackBaseRef selTrackBaseRef(trackRef);
00090     
00091     if(debug)
00092       cout << "runConvBremFinder:: pushing_back High Purity " << pft->trackRef()->pt()  
00093            << " eta,phi " << pft->trackRef()->eta() << ", " <<   pft->trackRef()->phi() 
00094            <<  " Memory Address Ref  " << &*trackRef << " Memory Address BaseRef  " << &*selTrackBaseRef << endl;  
00095     AllPFRecTracks.push_back(pfRecTrRef);
00096   }
00097   
00098 
00099   if(useConversions) {
00100     const PFConversionCollection& PfConvColl = *(pfConversions.product());
00101     for(unsigned i=0;i<PfConvColl.size(); i++) {
00102       reco::PFConversionRef convRef(pfConversions,i);
00103       
00104       unsigned int trackSize=(convRef->pfTracks()).size();
00105       if ( convRef->pfTracks().size() < 2) continue;
00106       for(unsigned iTk=0;iTk<trackSize; iTk++) {
00107         PFRecTrackRef compPFTkRef = convRef->pfTracks()[iTk];
00108         reco::TrackBaseRef newTrackBaseRef(compPFTkRef->trackRef());
00109         // do not consider the kf track already associated to the seed
00110         if(pfTrackRef.isNonnull()) {
00111           reco::TrackBaseRef primaryTrackBaseRef(pfTrackRef->trackRef());
00112           if(primaryTrackBaseRef == newTrackBaseRef) continue;
00113         }
00114         bool notFound = true;
00115         for(unsigned iPF = 0; iPF < AllPFRecTracks.size(); iPF++) {
00116           reco::TrackBaseRef selTrackBaseRef(AllPFRecTracks[iPF]->trackRef());
00117 
00118           if(debugRef) 
00119             cout << "## Track 1 HP pt " << AllPFRecTracks[iPF]->trackRef()->pt() << " eta, phi " << AllPFRecTracks[iPF]->trackRef()->eta() << ", " << AllPFRecTracks[iPF]->trackRef()->phi() 
00120                  << " Memory Address Ref  " << &(*AllPFRecTracks[iPF]->trackRef()) << " Memory Address BaseRef  " << &*selTrackBaseRef << endl;
00121           if(debugRef) 
00122             cout << "** Track 2 CONV pt " << compPFTkRef->trackRef()->pt() << " eta, phi " <<  compPFTkRef->trackRef()->eta() << ", " << compPFTkRef->trackRef()->phi() 
00123                  << " Memory Address Ref " << &*compPFTkRef->trackRef() <<  " Memory Address BaseRef " << &*newTrackBaseRef << endl;
00124           //if(selTrackBaseRef == newTrackBaseRef ||  AllPFRecTracks[iPF]->trackRef()== compPFTkRef->trackRef()) {
00125           if(AllPFRecTracks[iPF]->trackRef()== compPFTkRef->trackRef()) {
00126             if(debugRef) 
00127               cout << "  SAME BREM REF " << endl;
00128             notFound = false;
00129           }
00130         }
00131         if(notFound) {
00132           if(debug)
00133             cout << "runConvBremFinder:: pushing_back Conversions " << compPFTkRef->trackRef()->pt() 
00134                  << " eta,phi " << compPFTkRef->trackRef()->eta() << " phi " << compPFTkRef->trackRef()->phi() <<endl; 
00135           AllPFRecTracks.push_back(compPFTkRef);
00136         }
00137       }
00138     }
00139   }
00140 
00141   if(useNuclear) {
00142     const PFDisplacedTrackerVertexCollection& PfNuclColl = *(pfNuclears.product());
00143     for(unsigned i=0;i<PfNuclColl.size(); i++) {
00144       const reco::PFDisplacedTrackerVertexRef dispacedVertexRef(pfNuclears, i );
00145       unsigned int trackSize= dispacedVertexRef->pfRecTracks().size();
00146       for(unsigned iTk=0;iTk < trackSize; iTk++) {
00147         reco::PFRecTrackRef newPFRecTrackRef = dispacedVertexRef->pfRecTracks()[iTk]; 
00148         reco::TrackBaseRef newTrackBaseRef(newPFRecTrackRef->trackRef());
00149         // do not consider the kf track already associated to the seed
00150         if(pfTrackRef.isNonnull()) {
00151           reco::TrackBaseRef primaryTrackBaseRef(pfTrackRef->trackRef());
00152           if(primaryTrackBaseRef == newTrackBaseRef) continue;
00153         }
00154         bool notFound = true;
00155         for(unsigned iPF = 0; iPF < AllPFRecTracks.size(); iPF++) {
00156           reco::TrackBaseRef selTrackBaseRef(AllPFRecTracks[iPF]->trackRef());
00157           if(selTrackBaseRef == newTrackBaseRef) notFound = false;
00158         }
00159         if(notFound) {
00160           if(debug)
00161             cout << "runConvBremFinder:: pushing_back displaced Vertex pt " << newPFRecTrackRef->trackRef()->pt()  
00162                  << " eta,phi " << newPFRecTrackRef->trackRef()->eta() << ", " <<   newPFRecTrackRef->trackRef()->phi() <<  endl; 
00163           AllPFRecTracks.push_back(newPFRecTrackRef);
00164         }
00165       }
00166     }
00167   }
00168 
00169   if(useV0) {
00170     const PFV0Collection& PfV0Coll = *(pfV0.product());
00171     for(unsigned i=0;i<PfV0Coll.size(); i++) {
00172       reco::PFV0Ref v0Ref( pfV0, i );
00173       unsigned int trackSize=(v0Ref->pfTracks()).size();
00174       for(unsigned iTk=0;iTk<trackSize; iTk++) {
00175         reco::PFRecTrackRef newPFRecTrackRef = (v0Ref->pfTracks())[iTk]; 
00176         reco::TrackBaseRef newTrackBaseRef(newPFRecTrackRef->trackRef());
00177         // do not consider the kf track already associated to the seed
00178         if(pfTrackRef.isNonnull()) {
00179           reco::TrackBaseRef primaryTrackBaseRef(pfTrackRef->trackRef());
00180           if(primaryTrackBaseRef == newTrackBaseRef) continue;
00181         }
00182         bool notFound = true;
00183         for(unsigned iPF = 0; iPF < AllPFRecTracks.size(); iPF++) {
00184           reco::TrackBaseRef selTrackBaseRef(AllPFRecTracks[iPF]->trackRef());
00185           if(selTrackBaseRef == newTrackBaseRef) notFound = false;
00186         }
00187         if(notFound) {
00188           if(debug)
00189             cout << "runConvBremFinder:: pushing_back V0 " << newPFRecTrackRef->trackRef()->pt()  
00190                  << " eta,phi " << newPFRecTrackRef->trackRef()->eta() << ", " <<   newPFRecTrackRef->trackRef()->phi() << endl; 
00191           AllPFRecTracks.push_back(newPFRecTrackRef);
00192         }
00193       }
00194     }
00195   }
00196 
00197 
00198 
00199   pfRecTrRef_vec_.clear();
00200 
00201 
00202   for(unsigned iPF = 0; iPF < AllPFRecTracks.size(); iPF++) {
00203   
00204  
00205     double dphi= fabs(AllPFRecTracks[iPF]->trackRef()->phi()-refGsf->phi()); 
00206     if (dphi>TMath::Pi()) dphi-= TMath::TwoPi();
00207     double deta=fabs(AllPFRecTracks[iPF]->trackRef()->eta()-refGsf->eta());
00208      
00209     // limiting the phase space (just for saving cpu-time)
00210     if( fabs(dphi)> 1.0  || fabs(deta) > 0.4) continue;
00211     
00212   
00213     double minDEtaBremKF = 1000.;
00214     double minDPhiBremKF = 1000.;
00215     double minDRBremKF = 1000.;
00216     double minDEtaBremKFPos = 1000.;
00217     double minDPhiBremKFPos = 1000.;
00218     double minDRBremKFPos = 1000.;
00219     reco:: TrackRef  trkRef = AllPFRecTracks[iPF]->trackRef();
00220  
00221     double secEta = trkRef->innerMomentum().eta();
00222     double secPhi = trkRef->innerMomentum().phi();
00223     
00224     for(unsigned ipbrem = 0; ipbrem < primPFBrem.size(); ipbrem++) {
00225       if(primPFBrem[ipbrem].indTrajPoint() == 99) continue;
00226       const reco::PFTrajectoryPoint& atPrimECAL 
00227         = primPFBrem[ipbrem].extrapolatedPoint( reco::PFTrajectoryPoint::ECALEntrance );
00228       if( ! atPrimECAL.isValid() ) continue;
00229       double bremEta = atPrimECAL.momentum().Eta();
00230       double bremPhi = atPrimECAL.momentum().Phi();
00231 
00232       
00233       double deta = fabs(bremEta - secEta);
00234       double dphi = fabs(bremPhi - secPhi);
00235       if (dphi>TMath::Pi()) dphi-= TMath::TwoPi();
00236       double DR = sqrt(deta*deta + dphi*dphi);
00237       
00238       
00239       double detaPos = fabs(bremEta - trkRef->innerPosition().eta());
00240       double dphiPos = fabs(bremPhi - trkRef->innerPosition().phi());
00241       if (dphiPos>TMath::Pi()) dphiPos-= TMath::TwoPi();
00242       double DRPos = sqrt(detaPos*detaPos + dphiPos*dphiPos);
00243       
00244 
00245 
00246       // find the closest track tangent
00247       if(DR < minDRBremKF) {
00248         
00249         minDRBremKF = DR;
00250         minDEtaBremKF = deta;
00251         minDPhiBremKF = fabs(dphi);
00252       }
00253       
00254       if(DRPos < minDRBremKFPos) {
00255         minDRBremKFPos = DR;
00256         minDEtaBremKFPos = detaPos;
00257         minDPhiBremKFPos = fabs(dphiPos);
00258       }
00259   
00260     }
00261 
00262     //gsfR
00263     float gsfR = sqrt(refGsf->innerPosition().x()*refGsf->innerPosition().x() + 
00264                       refGsf->innerPosition().y()*refGsf->innerPosition().y() );        
00265     
00266     
00267     // secR
00268     secR = sqrt(trkRef->innerPosition().x()*trkRef->innerPosition().x() + 
00269                 trkRef->innerPosition().y()*trkRef->innerPosition().y() );   
00270     
00271   
00272     // apply loose selection (to be parallel) between the secondary track and brem-tangents.
00273     // Moreover if the secR is internal with respect to the GSF track by two pixel layers discard it.
00274     if( (minDPhiBremKF < 0.1 || minDPhiBremKFPos < 0.1) &&
00275         (minDEtaBremKF < 0.02 ||  minDEtaBremKFPos < 0.02)&&
00276         secR > (gsfR-8)) {
00277       
00278 
00279       if(debug)
00280         cout << "runConvBremFinder:: OK Find track and BREM close " 
00281              << " MinDphi " << minDPhiBremKF << " MinDeta " << minDEtaBremKF  << endl;
00282       
00283 
00284       float MinDist = 100000.;
00285       float EE_calib = 0.; 
00286       PFRecTrack pfrectrack = *AllPFRecTracks[iPF];
00287       pfrectrack.calculatePositionREP();
00288       // Find and ECAL associated cluster
00289       for (PFClusterCollection::const_iterator clus = theEClus.begin();
00290            clus != theEClus.end();
00291            clus++ ) {
00292         const math::XYZPoint gp_Clus = clus->position();
00293         double dist = -1.;
00294         PFCluster clust = *clus;
00295         clust.calculatePositionREP();
00296         dist = LinkByRecHit::testTrackAndClusterByRecHit(pfrectrack , clust );
00297         
00298         if(dist > 0.) {
00299           bool applyCrackCorrections = false;
00300           vector<double> ps1Ene(0);
00301           vector<double> ps2Ene(0);
00302           double ps1,ps2;
00303           ps1=ps2=0.;
00304           if(dist < MinDist) {
00305             MinDist = dist;
00306             EE_calib = pfcalib_.energyEm(*clus,ps1Ene,ps2Ene,ps1,ps2,applyCrackCorrections);
00307           }
00308         }
00309       }
00310       if(MinDist > 0. && MinDist < 100000.) {
00311 
00312         // compute all the input variables for conv brem selection
00313         
00314         secPout = sqrt(trkRef->outerMomentum().x()*trkRef->outerMomentum().x() +
00315                        trkRef->outerMomentum().y()*trkRef->outerMomentum().y() +
00316                        trkRef->outerMomentum().z()*trkRef->outerMomentum().z());
00317         
00318         secPin = sqrt(trkRef->innerMomentum().x()*trkRef->innerMomentum().x() +
00319                       trkRef->innerMomentum().y()*trkRef->innerMomentum().y() +
00320                       trkRef->innerMomentum().z()*trkRef->innerMomentum().z());
00321         
00322 
00323         // maybe put innter momentum pt? 
00324         ptRatioGsfKF = trkRef->pt()/(refGsf->ptMode());
00325         
00326         Vertex dummy;
00327         const Vertex *pv = &dummy;
00328         edm::Ref<VertexCollection> pvRef;
00329         if (primaryVertex->size() != 0) {
00330           pv = &*primaryVertex->begin();
00331           // we always use the first vertex (at the moment)
00332           pvRef = edm::Ref<VertexCollection>(primaryVertex, 0);
00333         } else { // create a dummy PV
00334           Vertex::Error e;
00335           e(0, 0) = 0.0015 * 0.0015;
00336           e(1, 1) = 0.0015 * 0.0015;
00337           e(2, 2) = 15. * 15.;
00338           Vertex::Point p(0, 0, 0);
00339           dummy = Vertex(p, e, 0, 0, 0);
00340         }
00341         
00342         
00343         // direction of the Gsf track
00344         GlobalVector direction(refGsf->innerMomentum().x(), 
00345                                refGsf->innerMomentum().y(), 
00346                                refGsf->innerMomentum().z());
00347         
00348         TransientTrack transientTrack = builder_.build(*trkRef);         
00349         sTIP = IPTools::signedTransverseImpactParameter(transientTrack, direction, *pv).second.significance();
00350         
00351 
00352         Epout = EE_calib/secPout;
00353         
00354         // eta distance brem-secondary kf track
00355         detaBremKF = minDEtaBremKF;
00356         
00357         // Number of commont hits
00358         trackingRecHit_iterator  nhit=refGsf->recHitsBegin();
00359         trackingRecHit_iterator  nhit_end=refGsf->recHitsEnd();
00360         unsigned int tmp_sh = 0;
00361         //uint ish=0;
00362         
00363         for (;nhit!=nhit_end;++nhit){
00364           if ((*nhit)->isValid()){
00365             trackingRecHit_iterator  ihit=trkRef->recHitsBegin();
00366             trackingRecHit_iterator  ihit_end=trkRef->recHitsEnd();
00367             for (;ihit!=ihit_end;++ihit){
00368               if ((*ihit)->isValid()) {
00369                 // method 1
00370                 if((*nhit)->sharesInput(&*(*ihit),TrackingRecHit::all))  tmp_sh++;
00371                 
00372                 // method 2 to switch in case of problem with rechit collections
00373                 //  if(((*ihit)->geographicalId()==(*nhit)->geographicalId())&&
00374                 //  (((*nhit)->localPosition()-(*ihit)->localPosition()).mag()<0.01)) ish++;
00375                 
00376                 
00377               }
00378             }
00379           }
00380         }
00381           
00382         nHITS1 = tmp_sh;
00383         
00384         double mvaValue = tmvaReader_->EvaluateMVA("BDT");
00385         
00386         if(debug) 
00387           cout << " The imput variables for conv brem tracks identification " << endl
00388                << " secR          " << secR << " gsfR " << gsfR  << endl
00389                << " N shared hits " << nHITS1 << endl
00390                << " sTIP          " << sTIP << endl
00391                << " detaBremKF    " << detaBremKF << endl
00392                << " E/pout        " << Epout << endl
00393                << " pin           " << secPin << endl
00394                << " ptRatioKFGsf  " << ptRatioGsfKF << endl
00395                << " ***** MVA ***** " << mvaValue << endl;
00396         
00397         if(mvaValue > mvaBremConvCut_) {
00398           found_ = true;
00399           pfRecTrRef_vec_.push_back(AllPFRecTracks[iPF]);
00400           
00401         }
00402       } // end MinDIST
00403     } // end selection kf - brem tangents
00404   } // loop on the kf tracks
00405 
00406 
00407 
00408 
00409 
00410     
00411 }
00412 
00413