CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch1/src/RecoParticleFlow/PFProducer/src/KDTreeLinkerTrackEcal.cc

Go to the documentation of this file.
00001 #include "RecoParticleFlow/PFProducer/interface/KDTreeLinkerTrackEcal.h"
00002 
00003 #include "DataFormats/ParticleFlowReco/interface/PFCluster.h"
00004 #include "TMath.h"
00005 
00006 KDTreeLinkerTrackEcal::KDTreeLinkerTrackEcal()
00007   : KDTreeLinkerBase()
00008 {}
00009 
00010 KDTreeLinkerTrackEcal::~KDTreeLinkerTrackEcal()
00011 {
00012   clear();
00013 }
00014 
00015 void
00016 KDTreeLinkerTrackEcal::insertTargetElt(reco::PFBlockElement     *track)
00017 {
00018   targetSet_.insert(track);
00019 }
00020 
00021 
00022 void
00023 KDTreeLinkerTrackEcal::insertFieldClusterElt(reco::PFBlockElement       *ecalCluster)
00024 {
00025   reco::PFClusterRef clusterref = ecalCluster->clusterRef();
00026 
00027   // This test is more or less done in PFBlockAlgo.h. In others cases, it should be switch on.
00028   //   if (!((clusterref->layer() == PFLayer::ECAL_ENDCAP) ||
00029   //    (clusterref->layer() == PFLayer::ECAL_BARREL)))
00030   //     return;
00031 
00032   const std::vector<reco::PFRecHitFraction> &fraction = clusterref->recHitFractions();
00033 
00034   // We create a list of ecalCluster
00035   fieldClusterSet_.insert(ecalCluster);
00036   for(size_t rhit = 0; rhit < fraction.size(); ++rhit) {
00037     const reco::PFRecHitRef& rh = fraction[rhit].recHitRef();
00038     double fract = fraction[rhit].fraction();
00039 
00040     if ((rh.isNull()) || (fract < 1E-4))
00041       continue;
00042       
00043     const reco::PFRecHit& rechit = *rh;
00044       
00045     // We save the links rechit to EcalClusters
00046     rechit2ClusterLinks_[&rechit].insert(ecalCluster);
00047     
00048     // We create a liste of rechits
00049     rechitsSet_.insert(&rechit);
00050   }
00051 }
00052 
00053 void 
00054 KDTreeLinkerTrackEcal::buildTree()
00055 {
00056   // List of pseudo-rechits that will be used to create the KDTree
00057   std::vector<KDTreeNodeInfo> eltList;
00058 
00059   // Filling of this list
00060   for(RecHitSet::const_iterator it = rechitsSet_.begin(); 
00061       it != rechitsSet_.end(); it++) {
00062     
00063     const reco::PFRecHit::REPPoint &posrep = (*it)->positionREP();
00064     
00065     KDTreeNodeInfo rh1 (*it, posrep.Eta(), posrep.Phi());
00066     eltList.push_back(rh1);
00067     
00068     // Here we solve the problem of phi circular set by duplicating some rechits
00069     // too close to -Pi (or to Pi) and adding (substracting) to them 2 * Pi.
00070     if (rh1.dim2 > (M_PI - getPhiOffset())) {
00071       double phi = rh1.dim2 - 2 * M_PI;
00072       KDTreeNodeInfo rh2(*it, posrep.Eta(), phi); 
00073       eltList.push_back(rh2);
00074     }
00075 
00076     if (rh1.dim2 < (M_PI * -1.0 + getPhiOffset())) {
00077       double phi = rh1.dim2 + 2 * M_PI;
00078       KDTreeNodeInfo rh3(*it, posrep.Eta(), phi); 
00079       eltList.push_back(rh3);
00080     }
00081   }
00082 
00083   // Here we define the upper/lower bounds of the 2D space (eta/phi).
00084   double phimin = -1.0 * M_PI - getPhiOffset();
00085   double phimax = M_PI + getPhiOffset();
00086 
00087   // etamin-etamax, phimin-phimax
00088   KDTreeBox region(-3.0, 3.0, phimin, phimax);
00089 
00090   // We may now build the KDTree
00091   tree_.build(eltList, region);
00092 }
00093 
00094 void
00095 KDTreeLinkerTrackEcal::searchLinks()
00096 {
00097   // Must of the code has been taken from LinkByRecHit.cc
00098 
00099   // We iterate over the tracks.
00100   for(BlockEltSet::iterator it = targetSet_.begin(); 
00101       it != targetSet_.end(); it++) {
00102         
00103     reco::PFRecTrackRef trackref = (*it)->trackRefPF();
00104 
00105     // We set the multilinks flag of the track to true. It will allow us to 
00106     // use in an optimized way our algo results in the recursive linking algo.
00107     (*it)->setIsValidMultilinks(true);
00108 
00109     const reco::PFTrajectoryPoint& atECAL = 
00110       trackref->extrapolatedPoint(reco::PFTrajectoryPoint::ECALShowerMax);
00111 
00112     // The track didn't reach ecal
00113     if( ! atECAL.isValid() ) continue;
00114     
00115     const reco::PFTrajectoryPoint& atVertex = 
00116       trackref->extrapolatedPoint( reco::PFTrajectoryPoint::ClosestApproach );
00117     
00118     double trackPt = sqrt(atVertex.momentum().Vect().Perp2());
00119     double tracketa = atECAL.positionREP().Eta();
00120     double trackphi = atECAL.positionREP().Phi();
00121     double trackx = atECAL.position().X();
00122     double tracky = atECAL.position().Y();
00123     double trackz = atECAL.position().Z();
00124     
00125     // Estimate the maximal envelope in phi/eta that will be used to find rechit candidates.
00126     // Same envelope for cap et barrel rechits.
00127     double range = getCristalPhiEtaMaxSize() * (2.0 + 1.0 / std::min(1., trackPt / 2.)); 
00128 
00129     // We search for all candidate recHits, ie all recHits contained in the maximal size envelope.
00130     std::vector<KDTreeNodeInfo> recHits;
00131     KDTreeBox trackBox(tracketa-range, tracketa+range, trackphi-range, trackphi+range);
00132     tree_.search(trackBox, recHits);
00133     
00134     // Here we check all rechit candidates using the non-approximated method.
00135     for(std::vector<KDTreeNodeInfo>::const_iterator rhit = recHits.begin(); 
00136         rhit != recHits.end(); ++rhit) {
00137            
00138       const std::vector< math::XYZPoint >& cornersxyz      = rhit->ptr->getCornersXYZ();
00139       const math::XYZPoint& posxyz                         = rhit->ptr->position();
00140       const reco::PFRecHit::REPPoint &rhrep                = rhit->ptr->positionREP();
00141       const std::vector<reco::PFRecHit::REPPoint>& corners = rhit->ptr->getCornersREP();
00142       if(corners.size() != 4) continue;
00143       
00144       double rhsizeEta = fabs(corners[0].Eta() - corners[2].Eta());
00145       double rhsizePhi = fabs(corners[0].Phi() - corners[2].Phi());
00146       if ( rhsizePhi > M_PI ) rhsizePhi = 2.*M_PI - rhsizePhi;
00147       
00148       double deta = fabs(rhrep.Eta() - tracketa);
00149       double dphi = fabs(rhrep.Phi() - trackphi);
00150       if ( dphi > M_PI ) dphi = 2.*M_PI - dphi;
00151       
00152       // Find all clusters associated to given rechit
00153       RecHit2BlockEltMap::iterator ret = rechit2ClusterLinks_.find(rhit->ptr);
00154       
00155       for(BlockEltSet::const_iterator clusterIt = ret->second.begin(); 
00156           clusterIt != ret->second.end(); clusterIt++) {
00157         
00158         reco::PFClusterRef clusterref = (*clusterIt)->clusterRef();
00159         double clusterz = clusterref->position().Z();
00160         int fracsNbr = clusterref->recHitFractions().size();
00161 
00162         if (clusterref->layer() == PFLayer::ECAL_BARREL){ // BARREL
00163           // Check if the track is in the barrel
00164           if (fabs(trackz) > 300.) continue;
00165 
00166           double _rhsizeEta = rhsizeEta * (2.00 + 1.0 / (fracsNbr * std::min(1.,trackPt/2.)));
00167           double _rhsizePhi = rhsizePhi * (2.00 + 1.0 / (fracsNbr * std::min(1.,trackPt/2.)));
00168           
00169           // Check if the track and the cluster are linked
00170           if(deta < (_rhsizeEta / 2.) && dphi < (_rhsizePhi / 2.))
00171             target2ClusterLinks_[*it].insert(*clusterIt);
00172 
00173           
00174         } else { // ENDCAP
00175 
00176           // Check if the track is in the cap
00177           if (fabs(trackz) < 300.) continue;
00178           if (trackz * clusterz < 0.) continue;
00179           
00180           double x[5];
00181           double y[5];
00182           for ( unsigned jc=0; jc<4; ++jc ) {
00183             math::XYZPoint cornerposxyz = cornersxyz[jc];
00184             x[jc] = cornerposxyz.X() + (cornerposxyz.X()-posxyz.X())
00185               * (1.00+0.50/fracsNbr /std::min(1.,trackPt/2.));
00186             y[jc] = cornerposxyz.Y() + (cornerposxyz.Y()-posxyz.Y())
00187               * (1.00+0.50/fracsNbr /std::min(1.,trackPt/2.));
00188           }
00189           
00190           x[4] = x[0];
00191           y[4] = y[0];
00192           
00193           bool isinside = TMath::IsInside(trackx,
00194                                           tracky,
00195                                           5,x,y);
00196           
00197           // Check if the track and the cluster are linked
00198           if( isinside )
00199             target2ClusterLinks_[*it].insert(*clusterIt);
00200         }
00201       }
00202     }
00203   }
00204 }
00205 
00206 void
00207 KDTreeLinkerTrackEcal::updatePFBlockEltWithLinks()
00208 {
00209   //TODO YG : Check if cluster positionREP() is valid ?
00210 
00211   // Here we save in each track the list of phi/eta values of linked clusters.
00212   for (BlockElt2BlockEltMap::iterator it = target2ClusterLinks_.begin();
00213        it != target2ClusterLinks_.end(); ++it) {
00214     reco::PFMultiLinksTC multitracks(true);
00215 
00216     for (BlockEltSet::iterator jt = it->second.begin();
00217          jt != it->second.end(); ++jt) {
00218 
00219       double clusterPhi = (*jt)->clusterRef()->positionREP().Phi();
00220       double clusterEta = (*jt)->clusterRef()->positionREP().Eta();
00221 
00222       multitracks.linkedClusters.push_back(std::make_pair(clusterPhi, clusterEta));
00223     }
00224 
00225     it->first->setMultilinks(multitracks);
00226   }
00227 }
00228 
00229 void
00230 KDTreeLinkerTrackEcal::clear()
00231 {
00232   targetSet_.clear();
00233   fieldClusterSet_.clear();
00234 
00235   rechitsSet_.clear();
00236 
00237   rechit2ClusterLinks_.clear();
00238   target2ClusterLinks_.clear();
00239 
00240   tree_.clear();
00241 }