CMS 3D CMS Logo

CMSSW_4_4_3_patch1/src/RecoParticleFlow/PFProducer/src/KDTreeLinkerTrackHcal.cc

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