CMS 3D CMS Logo

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