CMS 3D CMS Logo

KDTreeLinkerTrackHcal.cc
Go to the documentation of this file.
4 
5 // This class is used to find all links between Tracks and HCAL clusters
6 // using a KDTree algorithm.
7 // It is used in PFBlockAlgo.cc in the function links().
9 public:
11  ~KDTreeLinkerTrackHcal() override;
12 
13  // With this method, we create the list of psCluster that we want to link.
15 
16  // Here, we create the list of hcalCluster that we want to link. From hcalCluster
17  // and fraction, we will create a second list of rechits that will be used to
18  // build the KDTree.
19  void insertFieldClusterElt(reco::PFBlockElement* hcalCluster) override;
20 
21  // The KDTree building from rechits list.
22  void buildTree() override;
23 
24  // Here we will iterate over all tracks. For each track intersection point with HCAL,
25  // we will search the closest rechits in the KDTree, from rechits we will find the
26  // hcalClusters and after that we will check the links between the track and
27  // all closest hcalClusters.
28  void searchLinks() override;
29 
30  // Here, we will store all PS/HCAL founded links in the PFBlockElement class
31  // of each psCluster in the PFmultilinks field.
32  void updatePFBlockEltWithLinks() override;
33 
34  // Here we free all allocated structures.
35  void clear() override;
36 
37 private:
38  // Data used by the KDTree algorithm : sets of Tracks and HCAL clusters.
41 
42  // Sets of rechits that compose the HCAL clusters.
44 
45  // Map of linked Track/HCAL clusters.
47 
48  // Map of the HCAL clusters associated to a rechit.
50 
51  // KD trees
53 
54  // TrajectoryPoints
59  bool checkExit_;
60 };
61 
62 // the text name is different so that we can easily
63 // construct it when calling the factory
64 DEFINE_EDM_PLUGIN(KDTreeLinkerFactory, KDTreeLinkerTrackHcal, "KDTreeTrackAndHCALLinker");
65 
67  : KDTreeLinkerBase(conf),
68  trajectoryLayerEntranceString_(conf.getParameter<std::string>("trajectoryLayerEntrance")),
69  trajectoryLayerExitString_(conf.getParameter<std::string>("trajectoryLayerExit")) {
70  // Initialization
72  phiOffset_ = 0.32;
73  // convert TrajectoryLayers info from string to enum
76  // make sure the requested setting is supported
83  // flag if exit layer should be checked or not
85 }
86 
88 
90  if (track->trackRefPF()->extrapolatedPoint(trajectoryLayerEntrance_).isValid()) {
91  targetSet_.insert(track);
92  }
93 }
94 
96  const reco::PFClusterRef& clusterref = hcalCluster->clusterRef();
97 
98  // This test is more or less done in PFBlockAlgo.h. In others cases, it should be switch on.
99  // if (!((clusterref->layer() == PFLayer::HCAL_ENDCAP) ||
100  // (clusterref->layer() == PFLayer::HCAL_BARREL1)))
101  // return;
102 
103  const std::vector<reco::PFRecHitFraction>& fraction = clusterref->recHitFractions();
104 
105  // We create a list of hcalCluster
106  fieldClusterSet_.insert(hcalCluster);
107  for (size_t rhit = 0; rhit < fraction.size(); ++rhit) {
108  const reco::PFRecHitRef& rh = fraction[rhit].recHitRef();
109  double fract = fraction[rhit].fraction();
110 
111  if ((rh.isNull()) || (fract < cutOffFrac))
112  continue;
113 
114  const reco::PFRecHit& rechit = *rh;
115 
116  // We save the links rechit to HcalClusters
117  rechit2ClusterLinks_[&rechit].insert(hcalCluster);
118 
119  // We create a liste of rechits
120  rechitsSet_.insert(&rechit);
121  }
122 }
123 
125  // List of pseudo-rechits that will be used to create the KDTree
126  std::vector<KDTreeNodeInfo<reco::PFRecHit const*, 2>> eltList;
127 
128  // Filling of this list
129  for (RecHitSet::const_iterator it = rechitsSet_.begin(); it != rechitsSet_.end(); it++) {
130  const reco::PFRecHit::REPPoint& posrep = (*it)->positionREP();
131 
132  KDTreeNodeInfo<reco::PFRecHit const*, 2> rh1(*it, posrep.eta(), posrep.phi());
133  eltList.push_back(rh1);
134 
135  // Here we solve the problem of phi circular set by duplicating some rechits
136  // too close to -Pi (or to Pi) and adding (substracting) to them 2 * Pi.
137  if (rh1.dims[1] > (M_PI - phiOffset_)) {
138  float phi = rh1.dims[1] - 2 * M_PI;
139  KDTreeNodeInfo<reco::PFRecHit const*, 2> rh2(*it, float(posrep.eta()), phi);
140  eltList.push_back(rh2);
141  }
142 
143  if (rh1.dims[1] < (M_PI * -1.0 + phiOffset_)) {
144  float phi = rh1.dims[1] + 2 * M_PI;
145  KDTreeNodeInfo<reco::PFRecHit const*, 2> rh3(*it, float(posrep.eta()), phi);
146  eltList.push_back(rh3);
147  }
148  }
149 
150  // Here we define the upper/lower bounds of the 2D space (eta/phi).
151  float phimin = -1.0 * M_PI - phiOffset_;
152  float phimax = M_PI + phiOffset_;
153 
154  // etamin-etamax, phimin-phimax
155  KDTreeBox region(-3.0f, 3.0f, phimin, phimax);
156 
157  // We may now build the KDTree
158  tree_.build(eltList, region);
159 }
160 
162  // Must of the code has been taken from LinkByRecHit.cc
163 
164  // We iterate over the tracks.
165  for (BlockEltSet::iterator it = targetSet_.begin(); it != targetSet_.end(); it++) {
166  reco::PFRecTrackRef trackref = (*it)->trackRefPF();
167 
168  const reco::PFTrajectoryPoint& atHCAL = trackref->extrapolatedPoint(trajectoryLayerEntrance_);
169 
170  // The track didn't reach hcal
171  if (!atHCAL.isValid())
172  continue;
173 
174  // In case the exit point check is requested, check eta and phi differences between entrance and exit
175  double dHeta = 0.0;
176  float dHphi = 0.0;
177  if (checkExit_) {
178  const reco::PFTrajectoryPoint& atHCALExit = trackref->extrapolatedPoint(trajectoryLayerExit_);
179  dHeta = atHCALExit.positionREP().eta() - atHCAL.positionREP().eta();
180  dHphi = atHCALExit.positionREP().phi() - atHCAL.positionREP().phi();
181  if (dHphi > M_PI)
182  dHphi = dHphi - 2. * M_PI;
183  else if (dHphi < -M_PI)
184  dHphi = dHphi + 2. * M_PI;
185  } // checkExit_
186 
187  float tracketa = atHCAL.positionREP().eta() + 0.1 * dHeta;
188  float trackphi = atHCAL.positionREP().phi() + 0.1 * dHphi;
189 
190  if (trackphi > M_PI)
191  trackphi -= 2 * M_PI;
192  else if (trackphi < -M_PI)
193  trackphi += 2 * M_PI;
194 
195  // Estimate the maximal envelope in phi/eta that will be used to find rechit candidates.
196  // Same envelope for cap et barrel rechits.
197  double inflation = 1.;
198  float rangeeta = (cristalPhiEtaMaxSize_ * (1.5 + 0.5) + 0.2 * fabs(dHeta)) * inflation;
199  float rangephi = (cristalPhiEtaMaxSize_ * (1.5 + 0.5) + 0.2 * fabs(dHphi)) * inflation;
200 
201  // We search for all candidate recHits, ie all recHits contained in the maximal size envelope.
202  std::vector<reco::PFRecHit const*> recHits;
203  KDTreeBox trackBox(tracketa - rangeeta, tracketa + rangeeta, trackphi - rangephi, trackphi + rangephi);
204  tree_.search(trackBox, recHits);
205 
206  // Here we check all rechit candidates using the non-approximated method.
207  for (auto const& recHit : recHits) {
208  const auto& rhrep = recHit->positionREP();
209  const auto& corners = recHit->getCornersREP();
210 
211  double rhsizeeta = fabs(corners[3].eta() - corners[1].eta());
212  double rhsizephi = fabs(corners[3].phi() - corners[1].phi());
213  if (rhsizephi > M_PI)
214  rhsizephi = 2. * M_PI - rhsizephi;
215 
216  double deta = fabs(rhrep.eta() - tracketa);
217  double dphi = fabs(rhrep.phi() - trackphi);
218  if (dphi > M_PI)
219  dphi = 2. * M_PI - dphi;
220 
221  // Find all clusters associated to given rechit
222  RecHit2BlockEltMap::iterator ret = rechit2ClusterLinks_.find(recHit);
223 
224  for (BlockEltSet::iterator clusterIt = ret->second.begin(); clusterIt != ret->second.end(); clusterIt++) {
225  const reco::PFClusterRef clusterref = (*clusterIt)->clusterRef();
226  int fracsNbr = clusterref->recHitFractions().size();
227 
228  double _rhsizeeta = rhsizeeta * (1.5 + 0.5 / fracsNbr) + 0.2 * fabs(dHeta);
229  double _rhsizephi = rhsizephi * (1.5 + 0.5 / fracsNbr) + 0.2 * fabs(dHphi);
230 
231  // Check if the track and the cluster are linked
232  if (deta < (_rhsizeeta / 2.) && dphi < (_rhsizephi / 2.))
233  cluster2TargetLinks_[*clusterIt].insert(*it);
234  }
235  }
236  }
237 }
238 
240  //TODO YG : Check if cluster positionREP() is valid ?
241 
242  // Here we save in each HCAL cluster the list of phi/eta values of linked clusters.
243  for (BlockElt2BlockEltMap::iterator it = cluster2TargetLinks_.begin(); it != cluster2TargetLinks_.end(); ++it) {
244  reco::PFMultiLinksTC multitracks(true);
245 
246  for (BlockEltSet::iterator jt = it->second.begin(); jt != it->second.end(); ++jt) {
247  reco::PFRecTrackRef trackref = (*jt)->trackRefPF();
248  const reco::PFTrajectoryPoint& atHCAL = trackref->extrapolatedPoint(trajectoryLayerEntrance_);
249  double tracketa = atHCAL.positionREP().eta();
250  double trackphi = atHCAL.positionREP().phi();
251 
252  multitracks.linkedClusters.push_back(std::make_pair(trackphi, tracketa));
253  }
254 
255  it->first->setMultilinks(multitracks);
256  }
257 
258  // We set the multilinks flag of the track to true. It will allow us to
259  // use in an optimized way our algo results in the recursive linking algo.
260  for (BlockEltSet::iterator it = fieldClusterSet_.begin(); it != fieldClusterSet_.end(); ++it)
261  (*it)->setIsValidMultilinks(true);
262 }
263 
265  targetSet_.clear();
266  fieldClusterSet_.clear();
267 
268  rechitsSet_.clear();
269 
270  rechit2ClusterLinks_.clear();
271  cluster2TargetLinks_.clear();
272 
273  tree_.clear();
274 }
float phi() const
momentum azimuthal angle
Definition: PtEtaPhiMass.h:54
const REPPoint & positionREP() const
trajectory position in (rho, eta, phi) base
std::set< reco::PFBlockElement * > BlockEltSet
Abstract base class for a PFBlock element (track, cluster...)
reco::PFTrajectoryPoint::LayerType trajectoryLayerEntrance_
void search(const KDTreeBox< DIM > &searchBox, std::vector< DATA > &resRecHitList)
ret
prodAgent to be discontinued
virtual const PFClusterRef & clusterRef() const
std::map< const reco::PFRecHit *, BlockEltSet > RecHit2BlockEltMap
std::map< reco::PFBlockElement *, BlockEltSet > BlockElt2BlockEltMap
KDTreeLinkerTrackHcal(const edm::ParameterSet &conf)
KDTreeLinkerAlgo< reco::PFRecHit const * > tree_
std::string trajectoryLayerEntranceString_
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:31
RecHit2BlockEltMap rechit2ClusterLinks_
PFMultilinksType linkedClusters
void updatePFBlockEltWithLinks() override
double f[11][100]
BlockElt2BlockEltMap cluster2TargetLinks_
std::set< const reco::PFRecHit * > RecHitSet
const Fraction< n, m >::type & fract()
Definition: Fraction.h:36
static LayerType layerTypeByName(const std::string &name)
virtual const PFRecTrackRef & trackRefPF() const
bool isNull() const
Checks for null.
Definition: Ref.h:235
void build(std::vector< KDTreeNodeInfo< DATA, DIM > > &eltList, const KDTreeBox< DIM > &region)
#define M_PI
const float cutOffFrac
bool isValid() const
is this point valid ?
void insertTargetElt(reco::PFBlockElement *track) override
reco::PFTrajectoryPoint::LayerType trajectoryLayerExit_
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:52
LayerType
Define the different layers where the track can be propagated.