CMS 3D CMS Logo

KDTreeLinkerTrackEcal.cc
Go to the documentation of this file.
5 
6 #include "TMath.h"
7 
8 // This class is used to find all links between Tracks and ECAL clusters
9 // using a KDTree algorithm.
10 // It is used in PFBlockAlgo.cc in the function links().
12 public:
14  ~KDTreeLinkerTrackEcal() override;
15 
16  // With this method, we create the list of track that we want to link.
18 
19  // Here, we create the list of ecalCluster that we want to link. From ecalCluster
20  // and fraction, we will create a second list of rechits that will be used to
21  // build the KDTree.
22  void insertFieldClusterElt(reco::PFBlockElement *ecalCluster) override;
23 
24  // The KDTree building from rechits list.
25  void buildTree() override;
26 
27  // Here we will iterate over all tracks. For each track intersection point with ECAL,
28  // we will search the closest rechits in the KDTree, from rechits we will find the
29  // ecalClusters and after that we will check the links between the track and
30  // all closest ecalClusters.
31  void searchLinks() override;
32 
33  // Here, we will store all Track/ECAL founded links in the PFBlockElement class
34  // of each psCluster in the PFmultilinks field.
35  void updatePFBlockEltWithLinks() override;
36 
37  // Here we free all allocated structures.
38  void clear() override;
39 
40 private:
41  // Data used by the KDTree algorithm : sets of Tracks and ECAL clusters.
44 
45  // Sets of rechits that compose the ECAL clusters.
47 
48  // Map of linked Track/ECAL clusters.
50 
51  // Map of the ECAL clusters associated to a rechit.
53 
54  // KD trees
56 };
57 
58 // the text name is different so that we can easily
59 // construct it when calling the factory
60 DEFINE_EDM_PLUGIN(KDTreeLinkerFactory, KDTreeLinkerTrackEcal, "KDTreeTrackAndECALLinker");
61 
63 
65 
67  if (track->trackRefPF()->extrapolatedPoint(reco::PFTrajectoryPoint::ECALShowerMax).isValid()) {
68  targetSet_.insert(track);
69  }
70 }
71 
73  const reco::PFClusterRef &clusterref = ecalCluster->clusterRef();
74 
75  // This test is more or less done in PFBlockAlgo.h. In others cases, it should be switch on.
76  // if (!((clusterref->layer() == PFLayer::ECAL_ENDCAP) ||
77  // (clusterref->layer() == PFLayer::ECAL_BARREL)))
78  // return;
79 
80  const std::vector<reco::PFRecHitFraction> &fraction = clusterref->recHitFractions();
81 
82  // We create a list of ecalCluster
83  fieldClusterSet_.insert(ecalCluster);
84  for (size_t rhit = 0; rhit < fraction.size(); ++rhit) {
85  const reco::PFRecHitRef &rh = fraction[rhit].recHitRef();
86  double fract = fraction[rhit].fraction();
87 
88  if ((rh.isNull()) || (fract < cutOffFrac))
89  continue;
90 
91  const reco::PFRecHit &rechit = *rh;
92 
93  // We save the links rechit to EcalClusters
94  rechit2ClusterLinks_[&rechit].insert(ecalCluster);
95 
96  // We create a liste of rechits
97  rechitsSet_.insert(&rechit);
98  }
99 }
100 
102  // List of pseudo-rechits that will be used to create the KDTree
103  std::vector<KDTreeNodeInfo<reco::PFRecHit const *, 2>> eltList;
104 
105  // Filling of this list
106  for (RecHitSet::const_iterator it = rechitsSet_.begin(); it != rechitsSet_.end(); it++) {
107  const reco::PFRecHit::REPPoint &posrep = (*it)->positionREP();
108 
109  KDTreeNodeInfo<reco::PFRecHit const *, 2> rh1(*it, posrep.eta(), posrep.phi());
110  eltList.push_back(rh1);
111 
112  // Here we solve the problem of phi circular set by duplicating some rechits
113  // too close to -Pi (or to Pi) and adding (substracting) to them 2 * Pi.
114  if (rh1.dims[1] > (M_PI - phiOffset_)) {
115  float phi = rh1.dims[1] - 2 * M_PI;
116  KDTreeNodeInfo<reco::PFRecHit const *, 2> rh2(*it, float(posrep.eta()), phi);
117  eltList.push_back(rh2);
118  }
119 
120  if (rh1.dims[1] < (M_PI * -1.0 + phiOffset_)) {
121  float phi = rh1.dims[1] + 2 * M_PI;
122  KDTreeNodeInfo<reco::PFRecHit const *, 2> rh3(*it, float(posrep.eta()), phi);
123  eltList.push_back(rh3);
124  }
125  }
126 
127  // Here we define the upper/lower bounds of the 2D space (eta/phi).
128  float phimin = -1.0 * M_PI - phiOffset_;
129  float phimax = M_PI + phiOffset_;
130 
131  // etamin-etamax, phimin-phimax
132  KDTreeBox region(-3.0f, 3.0f, phimin, phimax);
133 
134  // We may now build the KDTree
135  tree_.build(eltList, region);
136 }
137 
139  // Most of the code has been taken from LinkByRecHit.cc
140 
141  // We iterate over the tracks.
142  for (BlockEltSet::iterator it = targetSet_.begin(); it != targetSet_.end(); it++) {
143  reco::PFRecTrackRef trackref = (*it)->trackRefPF();
144 
145  const reco::PFTrajectoryPoint &atECAL = trackref->extrapolatedPoint(reco::PFTrajectoryPoint::ECALShowerMax);
146 
147  // The track didn't reach ecal
148  if (!atECAL.isValid())
149  continue;
150 
151  const reco::PFTrajectoryPoint &atVertex = trackref->extrapolatedPoint(reco::PFTrajectoryPoint::ClosestApproach);
152 
153  double trackPt = sqrt(atVertex.momentum().Vect().Perp2());
154  float tracketa = atECAL.positionREP().eta();
155  float trackphi = atECAL.positionREP().phi();
156  double trackx = atECAL.position().X();
157  double tracky = atECAL.position().Y();
158  double trackz = atECAL.position().Z();
159 
160  // Estimate the maximal envelope in phi/eta that will be used to find rechit candidates.
161  // Same envelope for cap et barrel rechits.
162  float range = cristalPhiEtaMaxSize_ * (2.0 + 1.0 / std::min(1., trackPt / 2.));
163 
164  // We search for all candidate recHits, ie all recHits contained in the maximal size envelope.
165  std::vector<reco::PFRecHit const *> recHits;
166  KDTreeBox trackBox(tracketa - range, tracketa + range, trackphi - range, trackphi + range);
167  tree_.search(trackBox, recHits);
168 
169  // Here we check all rechit candidates using the non-approximated method.
170  for (auto const &recHit : recHits) {
171  const auto &cornersxyz = recHit->getCornersXYZ();
172  const auto &posxyz = recHit->position();
173  const auto &rhrep = recHit->positionREP();
174  const auto &corners = recHit->getCornersREP();
175 
176  double rhsizeeta = fabs(corners[3].eta() - corners[1].eta());
177  double rhsizephi = fabs(corners[3].phi() - corners[1].phi());
178  if (rhsizephi > M_PI)
179  rhsizephi = 2. * M_PI - rhsizephi;
180 
181  double deta = fabs(rhrep.eta() - tracketa);
182  double dphi = fabs(rhrep.phi() - trackphi);
183  if (dphi > M_PI)
184  dphi = 2. * M_PI - dphi;
185 
186  // Find all clusters associated to given rechit
187  RecHit2BlockEltMap::iterator ret = rechit2ClusterLinks_.find(recHit);
188 
189  for (BlockEltSet::const_iterator clusterIt = ret->second.begin(); clusterIt != ret->second.end(); clusterIt++) {
190  reco::PFClusterRef clusterref = (*clusterIt)->clusterRef();
191  double clusterz = clusterref->position().z();
192  int fracsNbr = clusterref->recHitFractions().size();
193 
194  if (clusterref->layer() == PFLayer::ECAL_BARREL) { // BARREL
195  // Check if the track is in the barrel
196  if (fabs(trackz) > 300.)
197  continue;
198 
199  double _rhsizeeta = rhsizeeta * (2.00 + 1.0 / (fracsNbr * std::min(1., trackPt / 2.)));
200  double _rhsizephi = rhsizephi * (2.00 + 1.0 / (fracsNbr * std::min(1., trackPt / 2.)));
201 
202  // Check if the track and the cluster are linked
203  if (deta < (_rhsizeeta / 2.) && dphi < (_rhsizephi / 2.))
204  cluster2TargetLinks_[*clusterIt].insert(*it);
205 
206  } else { // ENDCAP
207 
208  // Check if the track is in the cap
209  if (fabs(trackz) < 300.)
210  continue;
211  if (trackz * clusterz < 0.)
212  continue;
213 
214  double x[5];
215  double y[5];
216  for (unsigned jc = 0; jc < 4; ++jc) {
217  auto cornerposxyz = cornersxyz[jc];
218  x[3 - jc] = cornerposxyz.x() +
219  (cornerposxyz.x() - posxyz.x()) * (1.00 + 0.50 / fracsNbr / std::min(1., trackPt / 2.));
220  y[3 - jc] = cornerposxyz.y() +
221  (cornerposxyz.y() - posxyz.y()) * (1.00 + 0.50 / fracsNbr / std::min(1., trackPt / 2.));
222  }
223 
224  x[4] = x[0];
225  y[4] = y[0];
226 
227  bool isinside = TMath::IsInside(trackx, tracky, 5, x, y);
228 
229  // Check if the track and the cluster are linked
230  if (isinside)
231  cluster2TargetLinks_[*clusterIt].insert(*it);
232  }
233  }
234  }
235  }
236 }
237 
239  //TODO YG : Check if cluster positionREP() is valid ?
240 
241  // Here we save in each ECAL cluster the list of phi/eta values of linked tracks.
242  for (BlockElt2BlockEltMap::iterator it = cluster2TargetLinks_.begin(); it != cluster2TargetLinks_.end(); ++it) {
243  const auto &ecalElt = it->first;
244  const auto &trackEltSet = it->second;
245  reco::PFMultiLinksTC multitracks(true);
246 
247  for (const auto &trackElt : trackEltSet) {
248  const reco::PFRecTrackRef &trackref = trackElt->trackRefPF();
249 
250  reco::PFMultilink multiLink(trackref);
251  multitracks.linkedPFObjects.push_back(multiLink);
252 
253  // We set the multilinks flag of the track (for links to ECAL) to true. It will allow us to
254  // use it in an optimized way in prefilter
255  trackElt->setIsValidMultilinks(true, _fieldType);
256  }
257 
258  // We set multilinks of the ECAL element (for links to tracks)
259  ecalElt->setMultilinks(multitracks, _targetType);
260  }
261 }
262 
264  targetSet_.clear();
265  fieldClusterSet_.clear();
266 
267  rechitsSet_.clear();
268 
269  rechit2ClusterLinks_.clear();
270  cluster2TargetLinks_.clear();
271 
272  tree_.clear();
273 }
void insertTargetElt(reco::PFBlockElement *track) override
std::set< reco::PFBlockElement * > BlockEltSet
Abstract base class for a PFBlock element (track, cluster...)
const REPPoint & positionREP() const
trajectory position in (rho, eta, phi) base
void search(const KDTreeBox< DIM > &searchBox, std::vector< DATA > &resRecHitList)
KDTreeLinkerAlgo< reco::PFRecHit const * > tree_
reco::PFBlockElement::Type _fieldType
ret
prodAgent to be discontinued
const math::XYZTLorentzVector & momentum() const
4-momenta quadrivector
std::map< const reco::PFRecHit *, BlockEltSet > RecHit2BlockEltMap
std::map< reco::PFBlockElement *, BlockEltSet > BlockElt2BlockEltMap
void insertFieldClusterElt(reco::PFBlockElement *ecalCluster) override
bool isValid() const
is this point valid ?
KDTreeLinkerTrackEcal(const edm::ParameterSet &conf)
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:31
Point of closest approach from beam axis (initial point in the case of PFSimParticle) ...
T sqrt(T t)
Definition: SSEVec.h:23
float phi() const
momentum azimuthal angle
Definition: PtEtaPhiMass.h:54
virtual const PFClusterRef & clusterRef() const
void updatePFBlockEltWithLinks() override
double f[11][100]
std::set< const reco::PFRecHit * > RecHitSet
const Fraction< n, m >::type & fract()
Definition: Fraction.h:36
bool isNull() const
Checks for null.
Definition: Ref.h:229
void build(std::vector< KDTreeNodeInfo< DATA, DIM > > &eltList, const KDTreeBox< DIM > &region)
#define M_PI
const float cutOffFrac
const math::XYZPoint & position() const
cartesian position (x, y, z)
PFMultilinksType linkedPFObjects
float eta() const
momentum pseudorapidity
Definition: PtEtaPhiMass.h:52
RecHit2BlockEltMap rechit2ClusterLinks_
reco::PFBlockElement::Type _targetType
#define DEFINE_EDM_PLUGIN(factory, type, name)
A PFTrack holds several trajectory points, which basically contain the position and momentum of a tra...
BlockElt2BlockEltMap cluster2TargetLinks_