CMS 3D CMS Logo

KDTreeLinkerTrackHcal.cc
Go to the documentation of this file.
2 
4 #include "TMath.h"
5 
6 // the text name is different so that we can easily
7 // construct it when calling the factory
10  "KDTreeTrackAndHCALLinker");
11 
12 
15 {
17  phiOffset_ = 0.32;
18 }
19 
21 {
22  clear();
23 }
24 
25 void
27 {
28  if( track->trackRefPF()->extrapolatedPoint( reco::PFTrajectoryPoint::HCALEntrance ).isValid() ) {
29  targetSet_.insert(track);
30  }
31 }
32 
33 
34 void
36 {
37  const reco::PFClusterRef& clusterref = hcalCluster->clusterRef();
38 
39  // This test is more or less done in PFBlockAlgo.h. In others cases, it should be switch on.
40  // if (!((clusterref->layer() == PFLayer::HCAL_ENDCAP) ||
41  // (clusterref->layer() == PFLayer::HCAL_BARREL1)))
42  // return;
43 
44  const std::vector<reco::PFRecHitFraction> &fraction = clusterref->recHitFractions();
45 
46  // We create a list of hcalCluster
47  fieldClusterSet_.insert(hcalCluster);
48  for(size_t rhit = 0; rhit < fraction.size(); ++rhit) {
49  const reco::PFRecHitRef& rh = fraction[rhit].recHitRef();
50  double fract = fraction[rhit].fraction();
51 
52  if ((rh.isNull()) || (fract < 1E-4))
53  continue;
54 
55  const reco::PFRecHit& rechit = *rh;
56 
57  // We save the links rechit to HcalClusters
58  rechit2ClusterLinks_[&rechit].insert(hcalCluster);
59 
60  // We create a liste of rechits
61  rechitsSet_.insert(&rechit);
62  }
63 }
64 
65 void
67 {
68  // List of pseudo-rechits that will be used to create the KDTree
69  std::vector<KDTreeNodeInfo<reco::PFRecHit const*>> eltList;
70 
71  // Filling of this list
72  for(RecHitSet::const_iterator it = rechitsSet_.begin();
73  it != rechitsSet_.end(); it++) {
74 
75  const reco::PFRecHit::REPPoint &posrep = (*it)->positionREP();
76 
77  KDTreeNodeInfo<reco::PFRecHit const*> rh1 (*it, posrep.eta(), posrep.phi());
78  eltList.push_back(rh1);
79 
80  // Here we solve the problem of phi circular set by duplicating some rechits
81  // too close to -Pi (or to Pi) and adding (substracting) to them 2 * Pi.
82  if (rh1.dim[1] > (M_PI - phiOffset_)) {
83  double phi = rh1.dim[1] - 2 * M_PI;
84  KDTreeNodeInfo rh2(*it, posrep.eta(), phi);
85  eltList.push_back(rh2);
86  }
87 
88  if (rh1.dim[1] < (M_PI * -1.0 + phiOffset_)) {
89  double phi = rh1.dim[1] + 2 * M_PI;
90  KDTreeNodeInfo rh3(*it, posrep.eta(), phi);
91  eltList.push_back(rh3);
92  }
93  }
94 
95  // Here we define the upper/lower bounds of the 2D space (eta/phi).
96  double phimin = -1.0 * M_PI - phiOffset_;
97  double phimax = M_PI + phiOffset_;
98 
99  // etamin-etamax, phimin-phimax
100  KDTreeBox region(-3.0, 3.0, phimin, phimax);
101 
102  // We may now build the KDTree
103  tree_.build(eltList, region);
104 }
105 
106 void
108 {
109  // Must of the code has been taken from LinkByRecHit.cc
110 
111  // We iterate over the tracks.
112  for(BlockEltSet::iterator it = targetSet_.begin();
113  it != targetSet_.end(); it++) {
114 
115  reco::PFRecTrackRef trackref = (*it)->trackRefPF();
116 
117  const reco::PFTrajectoryPoint& atHCAL =
118  trackref->extrapolatedPoint(reco::PFTrajectoryPoint::HCALEntrance);
119  const reco::PFTrajectoryPoint& atHCALExit =
120  trackref->extrapolatedPoint(reco::PFTrajectoryPoint::HCALExit);
121 
122  // The track didn't reach hcal
123  if( ! atHCAL.isValid()) continue;
124 
125  double dHeta = atHCALExit.positionREP().eta() - atHCAL.positionREP().eta();
126  double dHphi = atHCALExit.positionREP().phi() - atHCAL.positionREP().phi();
127  if ( dHphi > M_PI ) dHphi = dHphi - 2. * M_PI;
128  else if ( dHphi < -M_PI ) dHphi = dHphi + 2. * M_PI;
129 
130  double tracketa = atHCAL.positionREP().eta() + 0.1 * dHeta;
131  double trackphi = atHCAL.positionREP().phi() + 0.1 * dHphi;
132 
133  if (trackphi > M_PI) trackphi -= 2 * M_PI;
134  else if (trackphi < -M_PI) trackphi += 2 * M_PI;
135 
136  // Estimate the maximal envelope in phi/eta that will be used to find rechit candidates.
137  // Same envelope for cap et barrel rechits.
138  double inflation = 1.;
139  double rangeeta = (cristalPhiEtaMaxSize_ * (1.5 + 0.5) + 0.2 * fabs(dHeta)) * inflation;
140  double rangephi = (cristalPhiEtaMaxSize_ * (1.5 + 0.5) + 0.2 * fabs(dHphi)) * inflation;
141 
142  // We search for all candidate recHits, ie all recHits contained in the maximal size envelope.
143  std::vector<reco::PFRecHit const*> recHits;
144  KDTreeBox trackBox(tracketa - rangeeta, tracketa + rangeeta,
145  trackphi - rangephi, trackphi + rangephi);
146  tree_.search(trackBox, recHits);
147 
148  // Here we check all rechit candidates using the non-approximated method.
149  for(auto const& recHit : recHits) {
150 
151  const auto &rhrep = recHit->positionREP();
152  const auto & corners = recHit->getCornersREP();
153 
154  double rhsizeeta = fabs(corners[3].eta() - corners[1].eta());
155  double rhsizephi = fabs(corners[3].phi() - corners[1].phi());
156  if ( rhsizephi > M_PI ) rhsizephi = 2.*M_PI - rhsizephi;
157 
158  double deta = fabs(rhrep.eta() - tracketa);
159  double dphi = fabs(rhrep.phi() - trackphi);
160  if ( dphi > M_PI ) dphi = 2.*M_PI - dphi;
161 
162  // Find all clusters associated to given rechit
163  RecHit2BlockEltMap::iterator ret = rechit2ClusterLinks_.find(recHit);
164 
165  for(BlockEltSet::iterator clusterIt = ret->second.begin();
166  clusterIt != ret->second.end(); clusterIt++) {
167 
168  const reco::PFClusterRef clusterref = (*clusterIt)->clusterRef();
169  int fracsNbr = clusterref->recHitFractions().size();
170 
171  double _rhsizeeta = rhsizeeta * (1.5 + 0.5 / fracsNbr) + 0.2 * fabs(dHeta);
172  double _rhsizephi = rhsizephi * (1.5 + 0.5 / fracsNbr) + 0.2 * fabs(dHphi);
173 
174  // Check if the track and the cluster are linked
175  if(deta < (_rhsizeeta / 2.) && dphi < (_rhsizephi / 2.))
176  cluster2TargetLinks_[*clusterIt].insert(*it);
177  }
178  }
179  }
180 }
181 
182 void
184 {
185  //TODO YG : Check if cluster positionREP() is valid ?
186 
187  // Here we save in each HCAL cluster the list of phi/eta values of linked clusters.
188  for (BlockElt2BlockEltMap::iterator it = cluster2TargetLinks_.begin();
189  it != cluster2TargetLinks_.end(); ++it) {
190  reco::PFMultiLinksTC multitracks(true);
191 
192  for (BlockEltSet::iterator jt = it->second.begin();
193  jt != it->second.end(); ++jt) {
194 
195  reco::PFRecTrackRef trackref = (*jt)->trackRefPF();
196  const reco::PFTrajectoryPoint& atHCAL =
197  trackref->extrapolatedPoint(reco::PFTrajectoryPoint::HCALEntrance);
198  double tracketa = atHCAL.positionREP().eta();
199  double trackphi = atHCAL.positionREP().phi();
200 
201  multitracks.linkedClusters.push_back(std::make_pair(trackphi, tracketa));
202  }
203 
204  it->first->setMultilinks(multitracks);
205  }
206 
207  // We set the multilinks flag of the track to true. It will allow us to
208  // use in an optimized way our algo results in the recursive linking algo.
209  for (BlockEltSet::iterator it = fieldClusterSet_.begin();
210  it != fieldClusterSet_.end(); ++it)
211  (*it)->setIsValidMultilinks(true);
212 
213 }
214 
215 void
217 {
218  targetSet_.clear();
219  fieldClusterSet_.clear();
220 
221  rechitsSet_.clear();
222 
223  rechit2ClusterLinks_.clear();
224  cluster2TargetLinks_.clear();
225 
226  tree_.clear();
227 }
float phi() const
momentum azimuthal angle
Definition: PtEtaPhiMass.h:53
const REPPoint & positionREP() const
trajectory position in (rho, eta, phi) base
Abstract base class for a PFBlock element (track, cluster...)
virtual const PFClusterRef & clusterRef() const
KDTreeLinkerAlgo< reco::PFRecHit const * > tree_
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:32
RecHit2BlockEltMap rechit2ClusterLinks_
PFMultilinksType linkedClusters
void updatePFBlockEltWithLinks() override
void search(const KDTreeBox &searchBox, std::vector< DATA > &resRecHitList)
BlockElt2BlockEltMap cluster2TargetLinks_
const Fraction< n, m >::type & fract()
Definition: Fraction.h:37
virtual const PFRecTrackRef & trackRefPF() const
bool isNull() const
Checks for null.
Definition: Ref.h:248
#define M_PI
void build(std::vector< KDTreeNodeInfo< DATA > > &eltList, const KDTreeBox &region)
bool isValid() const
is this point valid ?
void insertTargetElt(reco::PFBlockElement *track) override
void insertFieldClusterElt(reco::PFBlockElement *hcalCluster) override
#define DEFINE_EDM_PLUGIN(factory, type, name)
A PFTrack holds several trajectory points, which basically contain the position and momentum of a tra...
float eta() const
momentum pseudorapidity
Definition: PtEtaPhiMass.h:51