CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TrackDetectorAssociator.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: TrackAssociator
4 // Class: TrackDetectorAssociator
5 //
6 /*
7 
8  Description: <one line class summary>
9 
10  Implementation:
11  <Notes on implementation>
12 */
13 //
14 // Original Author: Dmytro Kovalskyi
15 // Created: Fri Apr 21 10:59:41 PDT 2006
16 //
17 //
18 
23 
27 
42 
43 // calorimeter and muon infos
51 
53 #include <stack>
54 #include <set>
55 
57 #include "Math/VectorUtil.h"
58 #include <algorithm>
59 
66 
67 using namespace reco;
68 
70  ivProp_ = nullptr;
71  useDefaultPropagator_ = false;
72 }
73 
75 
77  ivProp_ = ptr;
78  cachedTrajectory_.setPropagator(ivProp_);
79 }
80 
81 void TrackDetectorAssociator::useDefaultPropagator() { useDefaultPropagator_ = true; }
82 
84  // access the calorimeter geometry
85  theCaloGeometry_ = &iSetup.getData(parameters.theCaloGeometryToken);
86 
87  // get the tracking Geometry
88  theTrackingGeometry_ = &iSetup.getData(parameters.theTrackingGeometryToken);
89 
90  if (useDefaultPropagator_ && (!defProp_ || theMagneticFieldWatcher_.check(iSetup))) {
91  // setup propagator
92  const MagneticField* bField = &iSetup.getData(parameters.bFieldToken);
93 
94  auto prop = std::make_unique<SteppingHelixPropagator>(bField, anyDirection);
95  prop->setMaterialMode(false);
96  prop->applyRadX0Correction(true);
97  // prop->setDebug(true); // tmp
98  defProp_ = std::move(prop);
99  setPropagator(defProp_.get());
100  }
101 
102  ecalDetIdAssociator_ = &iSetup.getData(parameters.ecalDetIdAssociatorToken);
103  hcalDetIdAssociator_ = &iSetup.getData(parameters.hcalDetIdAssociatorToken);
104  hoDetIdAssociator_ = &iSetup.getData(parameters.hoDetIdAssociatorToken);
105  caloDetIdAssociator_ = &iSetup.getData(parameters.caloDetIdAssociatorToken);
106  muonDetIdAssociator_ = &iSetup.getData(parameters.muonDetIdAssociatorToken);
107  preshowerDetIdAssociator_ = &iSetup.getData(parameters.preshowerDetIdAssociatorToken);
108 }
109 
111  const edm::EventSetup& iSetup,
112  const FreeTrajectoryState& fts,
114  return associate(iEvent, iSetup, parameters, &fts);
115 }
116 
118  const edm::EventSetup& iSetup,
120  const FreeTrajectoryState* innerState,
121  const FreeTrajectoryState* outerState) {
123  if (!parameters.useEcal && !parameters.useCalo && !parameters.useHcal && !parameters.useHO && !parameters.useMuon &&
124  !parameters.usePreshower)
125  throw cms::Exception("ConfigurationError")
126  << "Configuration error! No subdetector was selected for the track association.";
127 
128  SteppingHelixStateInfo trackOrigin(*innerState);
129  info.stateAtIP = *innerState;
130  cachedTrajectory_.setStateAtIP(trackOrigin);
131 
132  init(iSetup, parameters);
133  // get track trajectory
134  // ECAL points (EB+EE)
135  // If the phi angle between a track entrance and exit points is more
136  // than 2 crystals, it is possible that the track will cross 3 crystals
137  // and therefore one has to check at least 3 points along the track
138  // trajectory inside ECAL. In order to have a chance to cross 4 crystalls
139  // in the barrel, a track should have P_t as low as 3 GeV or smaller
140  // If it's necessary, number of points along trajectory can be increased
141 
142  info.setCaloGeometry(theCaloGeometry_);
143 
144  cachedTrajectory_.reset_trajectory();
145  // estimate propagation outer boundaries based on
146  // requested sub-detector information. For now limit
147  // propagation region only if muon matching is not
148  // requested.
149  double HOmaxR = hoDetIdAssociator_->volume().maxR();
150  double HOmaxZ = hoDetIdAssociator_->volume().maxZ();
151  double minR = ecalDetIdAssociator_->volume().minR();
152  double minZ = preshowerDetIdAssociator_->volume().minZ();
153  cachedTrajectory_.setMaxHORadius(HOmaxR);
154  cachedTrajectory_.setMaxHOLength(HOmaxZ * 2.);
155  cachedTrajectory_.setMinDetectorRadius(minR);
156  cachedTrajectory_.setMinDetectorLength(minZ * 2.);
157 
158  if (parameters.useMuon) {
159  double maxR = muonDetIdAssociator_->volume().maxR();
160  double maxZ = muonDetIdAssociator_->volume().maxZ();
161  cachedTrajectory_.setMaxDetectorRadius(maxR);
162  cachedTrajectory_.setMaxDetectorLength(maxZ * 2.);
163  } else {
164  cachedTrajectory_.setMaxDetectorRadius(HOmaxR);
165  cachedTrajectory_.setMaxDetectorLength(HOmaxZ * 2.);
166  }
167 
168  // If track extras exist and outerState is before HO maximum, then use outerState
169  if (outerState) {
170  if (outerState->position().perp() < HOmaxR && std::abs(outerState->position().z()) < HOmaxZ) {
171  LogTrace("TrackAssociator") << "Using outerState as trackOrigin at Rho=" << outerState->position().perp()
172  << " Z=" << outerState->position().z() << "\n";
173  trackOrigin = SteppingHelixStateInfo(*outerState);
174  } else if (innerState) {
175  LogTrace("TrackAssociator") << "Using innerState as trackOrigin at Rho=" << innerState->position().perp()
176  << " Z=" << innerState->position().z() << "\n";
177  trackOrigin = SteppingHelixStateInfo(*innerState);
178  }
179  }
180 
181  if (trackOrigin.momentum().mag() == 0)
182  return info;
183  if (edm::isNotFinite(trackOrigin.momentum().x()) or edm::isNotFinite(trackOrigin.momentum().y()) or
184  edm::isNotFinite(trackOrigin.momentum().z()))
185  return info;
186  if (!cachedTrajectory_.propagateAll(trackOrigin))
187  return info;
188 
189  // get trajectory in calorimeters
190  cachedTrajectory_.findEcalTrajectory(ecalDetIdAssociator_->volume());
191  cachedTrajectory_.findHcalTrajectory(hcalDetIdAssociator_->volume());
192  cachedTrajectory_.findHOTrajectory(hoDetIdAssociator_->volume());
193  cachedTrajectory_.findPreshowerTrajectory(preshowerDetIdAssociator_->volume());
194 
195  info.trkGlobPosAtEcal = getPoint(cachedTrajectory_.getStateAtEcal().position());
196  info.trkGlobPosAtHcal = getPoint(cachedTrajectory_.getStateAtHcal().position());
197  info.trkGlobPosAtHO = getPoint(cachedTrajectory_.getStateAtHO().position());
198 
199  info.trkMomAtEcal = cachedTrajectory_.getStateAtEcal().momentum();
200  info.trkMomAtHcal = cachedTrajectory_.getStateAtHcal().momentum();
201  info.trkMomAtHO = cachedTrajectory_.getStateAtHO().momentum();
202 
203  if (parameters.useEcal)
204  fillEcal(iEvent, info, parameters);
205  if (parameters.useCalo)
206  fillCaloTowers(iEvent, info, parameters);
207  if (parameters.useHcal)
208  fillHcal(iEvent, info, parameters);
209  if (parameters.useHO)
210  fillHO(iEvent, info, parameters);
211  if (parameters.usePreshower)
212  fillPreshower(iEvent, info, parameters);
213  if (parameters.useMuon)
214  fillMuon(iEvent, info, parameters);
215  if (parameters.truthMatch)
216  fillCaloTruth(iEvent, info, parameters);
217 
218  return info;
219 }
220 
224  const std::vector<SteppingHelixStateInfo>& trajectoryStates = cachedTrajectory_.getEcalTrajectory();
225 
226  for (std::vector<SteppingHelixStateInfo>::const_iterator itr = trajectoryStates.begin();
227  itr != trajectoryStates.end();
228  itr++)
229  LogTrace("TrackAssociator") << "ECAL trajectory point (rho, z, phi): " << itr->position().perp() << ", "
230  << itr->position().z() << ", " << itr->position().phi();
231 
232  std::vector<GlobalPoint> coreTrajectory;
233  for (std::vector<SteppingHelixStateInfo>::const_iterator itr = trajectoryStates.begin();
234  itr != trajectoryStates.end();
235  itr++)
236  coreTrajectory.push_back(itr->position());
237 
238  if (coreTrajectory.empty()) {
239  LogTrace("TrackAssociator") << "ECAL track trajectory is empty; moving on\n";
240  info.isGoodEcal = false;
241  return;
242  }
243  info.isGoodEcal = true;
244 
245  // Find ECAL crystals
247  iEvent.getByToken(parameters.EBRecHitsToken, EBRecHits);
248  if (!EBRecHits.isValid())
249  throw cms::Exception("FatalError") << "Unable to find EBRecHitCollection in the event!\n";
250 
252  iEvent.getByToken(parameters.EERecHitsToken, EERecHits);
253  if (!EERecHits.isValid())
254  throw cms::Exception("FatalError") << "Unable to find EERecHitCollection in event!\n";
255 
256  std::set<DetId> ecalIdsInRegion;
257  if (parameters.accountForTrajectoryChangeCalo) {
258  // get trajectory change with respect to initial state
259  DetIdAssociator::MapRange mapRange =
260  getMapRange(cachedTrajectory_.trajectoryDelta(CachedTrajectory::IpToEcal), parameters.dREcalPreselection);
261  ecalIdsInRegion = ecalDetIdAssociator_->getDetIdsCloseToAPoint(coreTrajectory[0], mapRange);
262  } else
263  ecalIdsInRegion = ecalDetIdAssociator_->getDetIdsCloseToAPoint(coreTrajectory[0], parameters.dREcalPreselection);
264  LogTrace("TrackAssociator") << "ECAL hits in the region: " << ecalIdsInRegion.size();
265  if (parameters.dREcalPreselection > parameters.dREcal)
266  ecalIdsInRegion = ecalDetIdAssociator_->getDetIdsInACone(ecalIdsInRegion, coreTrajectory, parameters.dREcal);
267  LogTrace("TrackAssociator") << "ECAL hits in the cone: " << ecalIdsInRegion.size();
268  info.crossedEcalIds = ecalDetIdAssociator_->getCrossedDetIds(ecalIdsInRegion, coreTrajectory);
269  const std::vector<DetId>& crossedEcalIds = info.crossedEcalIds;
270  LogTrace("TrackAssociator") << "ECAL crossed hits " << crossedEcalIds.size();
271 
272  // add EcalRecHits
273  for (std::vector<DetId>::const_iterator itr = crossedEcalIds.begin(); itr != crossedEcalIds.end(); itr++) {
274  std::vector<EcalRecHit>::const_iterator ebHit = (*EBRecHits).find(*itr);
275  std::vector<EcalRecHit>::const_iterator eeHit = (*EERecHits).find(*itr);
276  if (ebHit != (*EBRecHits).end())
277  info.crossedEcalRecHits.push_back(&*ebHit);
278  else if (eeHit != (*EERecHits).end())
279  info.crossedEcalRecHits.push_back(&*eeHit);
280  else
281  LogTrace("TrackAssociator") << "Crossed EcalRecHit is not found for DetId: " << itr->rawId();
282  }
283  for (std::set<DetId>::const_iterator itr = ecalIdsInRegion.begin(); itr != ecalIdsInRegion.end(); itr++) {
284  std::vector<EcalRecHit>::const_iterator ebHit = (*EBRecHits).find(*itr);
285  std::vector<EcalRecHit>::const_iterator eeHit = (*EERecHits).find(*itr);
286  if (ebHit != (*EBRecHits).end())
287  info.ecalRecHits.push_back(&*ebHit);
288  else if (eeHit != (*EERecHits).end())
289  info.ecalRecHits.push_back(&*eeHit);
290  else
291  LogTrace("TrackAssociator") << "EcalRecHit from the cone is not found for DetId: " << itr->rawId();
292  }
293 }
294 
298  // use ECAL and HCAL trajectories to match a tower. (HO isn't used for matching).
299  std::vector<GlobalPoint> trajectory;
300  const std::vector<SteppingHelixStateInfo>& ecalTrajectoryStates = cachedTrajectory_.getEcalTrajectory();
301  const std::vector<SteppingHelixStateInfo>& hcalTrajectoryStates = cachedTrajectory_.getHcalTrajectory();
302  for (std::vector<SteppingHelixStateInfo>::const_iterator itr = ecalTrajectoryStates.begin();
303  itr != ecalTrajectoryStates.end();
304  itr++)
305  trajectory.push_back(itr->position());
306  for (std::vector<SteppingHelixStateInfo>::const_iterator itr = hcalTrajectoryStates.begin();
307  itr != hcalTrajectoryStates.end();
308  itr++)
309  trajectory.push_back(itr->position());
310 
311  if (trajectory.empty()) {
312  LogTrace("TrackAssociator") << "HCAL trajectory is empty; moving on\n";
313  info.isGoodCalo = false;
314  return;
315  }
316  info.isGoodCalo = true;
317 
318  // find crossed CaloTowers
320  iEvent.getByToken(parameters.caloTowersToken, caloTowers);
321  if (!caloTowers.isValid())
322  throw cms::Exception("FatalError") << "Unable to find CaloTowers in event!\n";
323 
324  std::set<DetId> caloTowerIdsInRegion;
325  if (parameters.accountForTrajectoryChangeCalo) {
326  // get trajectory change with respect to initial state
327  DetIdAssociator::MapRange mapRange =
328  getMapRange(cachedTrajectory_.trajectoryDelta(CachedTrajectory::IpToHcal), parameters.dRHcalPreselection);
329  caloTowerIdsInRegion = caloDetIdAssociator_->getDetIdsCloseToAPoint(trajectory[0], mapRange);
330  } else
331  caloTowerIdsInRegion = caloDetIdAssociator_->getDetIdsCloseToAPoint(trajectory[0], parameters.dRHcalPreselection);
332 
333  LogTrace("TrackAssociator") << "Towers in the region: " << caloTowerIdsInRegion.size();
334 
335  auto caloTowerIdsInAConeBegin = caloTowerIdsInRegion.begin();
336  auto caloTowerIdsInAConeEnd = caloTowerIdsInRegion.end();
337  std::set<DetId> caloTowerIdsInAConeTmp;
338  if (!caloDetIdAssociator_->selectAllInACone(parameters.dRHcal)) {
339  caloTowerIdsInAConeTmp =
340  caloDetIdAssociator_->getDetIdsInACone(caloTowerIdsInRegion, trajectory, parameters.dRHcal);
341  caloTowerIdsInAConeBegin = caloTowerIdsInAConeTmp.begin();
342  caloTowerIdsInAConeEnd = caloTowerIdsInAConeTmp.end();
343  }
344  LogTrace("TrackAssociator") << "Towers in the cone: "
345  << std::distance(caloTowerIdsInAConeBegin, caloTowerIdsInAConeEnd);
346 
347  info.crossedTowerIds = caloDetIdAssociator_->getCrossedDetIds(caloTowerIdsInRegion, trajectory);
348  const std::vector<DetId>& crossedCaloTowerIds = info.crossedTowerIds;
349  LogTrace("TrackAssociator") << "Towers crossed: " << crossedCaloTowerIds.size();
350 
351  // add CaloTowers
352  for (std::vector<DetId>::const_iterator itr = crossedCaloTowerIds.begin(); itr != crossedCaloTowerIds.end(); itr++) {
353  CaloTowerCollection::const_iterator tower = (*caloTowers).find(*itr);
354  if (tower != (*caloTowers).end())
355  info.crossedTowers.push_back(&*tower);
356  else
357  LogTrace("TrackAssociator") << "Crossed CaloTower is not found for DetId: " << (*itr).rawId();
358  }
359 
360  for (std::set<DetId>::const_iterator itr = caloTowerIdsInAConeBegin; itr != caloTowerIdsInAConeEnd; itr++) {
361  CaloTowerCollection::const_iterator tower = (*caloTowers).find(*itr);
362  if (tower != (*caloTowers).end())
363  info.towers.push_back(&*tower);
364  else
365  LogTrace("TrackAssociator") << "CaloTower from the cone is not found for DetId: " << (*itr).rawId();
366  }
367 }
368 
372  std::vector<GlobalPoint> trajectory;
373  const std::vector<SteppingHelixStateInfo>& trajectoryStates = cachedTrajectory_.getPreshowerTrajectory();
374  for (std::vector<SteppingHelixStateInfo>::const_iterator itr = trajectoryStates.begin();
375  itr != trajectoryStates.end();
376  itr++)
377  trajectory.push_back(itr->position());
378 
379  if (trajectory.empty()) {
380  LogTrace("TrackAssociator") << "Preshower trajectory is empty; moving on\n";
381  return;
382  }
383 
384  std::set<DetId> idsInRegion =
385  preshowerDetIdAssociator_->getDetIdsCloseToAPoint(trajectory[0], parameters.dRPreshowerPreselection);
386 
387  LogTrace("TrackAssociator") << "Number of Preshower Ids in the region: " << idsInRegion.size();
388  info.crossedPreshowerIds = preshowerDetIdAssociator_->getCrossedDetIds(idsInRegion, trajectory);
389  LogTrace("TrackAssociator") << "Number of Preshower Ids in crossed: " << info.crossedPreshowerIds.size();
390 }
391 
395  const std::vector<SteppingHelixStateInfo>& trajectoryStates = cachedTrajectory_.getHcalTrajectory();
396 
397  std::vector<GlobalPoint> coreTrajectory;
398  for (std::vector<SteppingHelixStateInfo>::const_iterator itr = trajectoryStates.begin();
399  itr != trajectoryStates.end();
400  itr++)
401  coreTrajectory.push_back(itr->position());
402 
403  if (coreTrajectory.empty()) {
404  LogTrace("TrackAssociator") << "HCAL trajectory is empty; moving on\n";
405  info.isGoodHcal = false;
406  return;
407  }
408  info.isGoodHcal = true;
409 
410  // find crossed Hcals
412  iEvent.getByToken(parameters.HBHEcollToken, collection);
413  if (!collection.isValid())
414  throw cms::Exception("FatalError") << "Unable to find HBHERecHits in event!\n";
415 
416  std::set<DetId> idsInRegion;
417  if (parameters.accountForTrajectoryChangeCalo) {
418  // get trajectory change with respect to initial state
419  DetIdAssociator::MapRange mapRange =
420  getMapRange(cachedTrajectory_.trajectoryDelta(CachedTrajectory::IpToHcal), parameters.dRHcalPreselection);
421  idsInRegion = hcalDetIdAssociator_->getDetIdsCloseToAPoint(coreTrajectory[0], mapRange);
422  } else
423  idsInRegion = hcalDetIdAssociator_->getDetIdsCloseToAPoint(coreTrajectory[0], parameters.dRHcalPreselection);
424 
425  LogTrace("TrackAssociator") << "HCAL hits in the region: " << idsInRegion.size() << "\n"
426  << DetIdInfo::info(idsInRegion, nullptr);
427 
428  auto idsInAConeBegin = idsInRegion.begin();
429  auto idsInAConeEnd = idsInRegion.end();
430  std::set<DetId> idsInAConeTmp;
431  if (!hcalDetIdAssociator_->selectAllInACone(parameters.dRHcal)) {
432  idsInAConeTmp = hcalDetIdAssociator_->getDetIdsInACone(idsInRegion, coreTrajectory, parameters.dRHcal);
433  idsInAConeBegin = idsInAConeTmp.begin();
434  idsInAConeEnd = idsInAConeTmp.end();
435  }
436  LogTrace("TrackAssociator") << "HCAL hits in the cone: " << std::distance(idsInAConeBegin, idsInAConeEnd) << "\n"
437  << DetIdInfo::info(std::set<DetId>(idsInAConeBegin, idsInAConeEnd), nullptr);
438  info.crossedHcalIds = hcalDetIdAssociator_->getCrossedDetIds(idsInRegion, coreTrajectory);
439  const std::vector<DetId>& crossedIds = info.crossedHcalIds;
440  LogTrace("TrackAssociator") << "HCAL hits crossed: " << crossedIds.size() << "\n"
441  << DetIdInfo::info(crossedIds, nullptr);
442 
443  // add Hcal
444  for (std::vector<DetId>::const_iterator itr = crossedIds.begin(); itr != crossedIds.end(); itr++) {
445  HBHERecHitCollection::const_iterator hit = (*collection).find(*itr);
446  if (hit != (*collection).end())
447  info.crossedHcalRecHits.push_back(&*hit);
448  else
449  LogTrace("TrackAssociator") << "Crossed HBHERecHit is not found for DetId: " << itr->rawId();
450  }
451  for (std::set<DetId>::const_iterator itr = idsInAConeBegin; itr != idsInAConeEnd; itr++) {
452  HBHERecHitCollection::const_iterator hit = (*collection).find(*itr);
453  if (hit != (*collection).end())
454  info.hcalRecHits.push_back(&*hit);
455  else
456  LogTrace("TrackAssociator") << "HBHERecHit from the cone is not found for DetId: " << itr->rawId();
457  }
458 }
459 
463  const std::vector<SteppingHelixStateInfo>& trajectoryStates = cachedTrajectory_.getHOTrajectory();
464 
465  std::vector<GlobalPoint> coreTrajectory;
466  for (std::vector<SteppingHelixStateInfo>::const_iterator itr = trajectoryStates.begin();
467  itr != trajectoryStates.end();
468  itr++)
469  coreTrajectory.push_back(itr->position());
470 
471  if (coreTrajectory.empty()) {
472  LogTrace("TrackAssociator") << "HO trajectory is empty; moving on\n";
473  info.isGoodHO = false;
474  return;
475  }
476  info.isGoodHO = true;
477 
478  // find crossed HOs
480  iEvent.getByToken(parameters.HOcollToken, collection);
481  if (!collection.isValid())
482  throw cms::Exception("FatalError") << "Unable to find HORecHits in event!\n";
483 
484  std::set<DetId> idsInRegion;
485  if (parameters.accountForTrajectoryChangeCalo) {
486  // get trajectory change with respect to initial state
487  DetIdAssociator::MapRange mapRange =
488  getMapRange(cachedTrajectory_.trajectoryDelta(CachedTrajectory::IpToHO), parameters.dRHcalPreselection);
489  idsInRegion = hoDetIdAssociator_->getDetIdsCloseToAPoint(coreTrajectory[0], mapRange);
490  } else
491  idsInRegion = hoDetIdAssociator_->getDetIdsCloseToAPoint(coreTrajectory[0], parameters.dRHcalPreselection);
492 
493  LogTrace("TrackAssociator") << "idsInRegion.size(): " << idsInRegion.size();
494 
495  auto idsInAConeBegin = idsInRegion.begin();
496  auto idsInAConeEnd = idsInRegion.end();
497  std::set<DetId> idsInAConeTmp;
498  if (!hoDetIdAssociator_->selectAllInACone(parameters.dRHcal)) {
499  idsInAConeTmp = hoDetIdAssociator_->getDetIdsInACone(idsInRegion, coreTrajectory, parameters.dRHcal);
500  idsInAConeBegin = idsInAConeTmp.begin();
501  idsInAConeEnd = idsInAConeTmp.end();
502  }
503  LogTrace("TrackAssociator") << "idsInACone.size(): " << std::distance(idsInAConeBegin, idsInAConeEnd);
504  info.crossedHOIds = hoDetIdAssociator_->getCrossedDetIds(idsInRegion, coreTrajectory);
505  const std::vector<DetId>& crossedIds = info.crossedHOIds;
506 
507  // add HO
508  for (std::vector<DetId>::const_iterator itr = crossedIds.begin(); itr != crossedIds.end(); itr++) {
509  HORecHitCollection::const_iterator hit = (*collection).find(*itr);
510  if (hit != (*collection).end())
511  info.crossedHORecHits.push_back(&*hit);
512  else
513  LogTrace("TrackAssociator") << "Crossed HORecHit is not found for DetId: " << itr->rawId();
514  }
515 
516  for (std::set<DetId>::const_iterator itr = idsInAConeBegin; itr != idsInAConeEnd; itr++) {
517  HORecHitCollection::const_iterator hit = (*collection).find(*itr);
518  if (hit != (*collection).end())
519  info.hoRecHits.push_back(&*hit);
520  else
521  LogTrace("TrackAssociator") << "HORecHit from the cone is not found for DetId: " << itr->rawId();
522  }
523 }
524 
526  const SimTrack& track,
527  const SimVertex& vertex) {
528  GlobalVector vector(track.momentum().x(), track.momentum().y(), track.momentum().z());
529  GlobalPoint point(vertex.position().x(), vertex.position().y(), vertex.position().z());
530 
531  int charge = track.type() > 0 ? -1 : 1; // lepton convention
532  if (abs(track.type()) == 211 || // pion
533  abs(track.type()) == 321 || // kaon
534  abs(track.type()) == 2212)
535  charge = track.type() < 0 ? -1 : 1;
536  return getFreeTrajectoryState(bField, vector, point, charge);
537 }
538 
540  const GlobalVector& momentum,
541  const GlobalPoint& vertex,
542  const int charge) {
543  GlobalTrajectoryParameters tPars(vertex, momentum, charge, bField);
544 
545  ROOT::Math::SMatrixIdentity id;
546  AlgebraicSymMatrix66 covT(id);
547  covT *= 1e-6; // initialize to sigma=1e-3
548  CartesianTrajectoryError tCov(covT);
549 
550  return FreeTrajectoryState(tPars, tCov);
551 }
552 
554  const reco::Track& track) {
555  GlobalVector vector(track.momentum().x(), track.momentum().y(), track.momentum().z());
556 
557  GlobalPoint point(track.vertex().x(), track.vertex().y(), track.vertex().z());
558 
560 
561  // FIX THIS !!!
562  // need to convert from perigee to global or helix (curvilinear) frame
563  // for now just an arbitrary matrix.
564  ROOT::Math::SMatrixIdentity id;
565  AlgebraicSymMatrix66 covT(id);
566  covT *= 1e-6; // initialize to sigma=1e-3
567  CartesianTrajectoryError tCov(covT);
568 
569  return FreeTrajectoryState(tPars, tCov);
570 }
571 
572 DetIdAssociator::MapRange TrackDetectorAssociator::getMapRange(const std::pair<float, float>& delta, const float dR) {
573  DetIdAssociator::MapRange mapRange;
574  mapRange.dThetaPlus = dR;
575  mapRange.dThetaMinus = dR;
576  mapRange.dPhiPlus = dR;
577  mapRange.dPhiMinus = dR;
578  if (delta.first > 0)
579  mapRange.dThetaPlus += delta.first;
580  else
581  mapRange.dThetaMinus += std::abs(delta.first);
582  if (delta.second > 0)
583  mapRange.dPhiPlus += delta.second;
584  else
585  mapRange.dPhiMinus += std::abs(delta.second);
586  LogTrace("TrackAssociator") << "Selection range: (dThetaPlus, dThetaMinus, dPhiPlus, dPhiMinus, dRPreselection): "
587  << mapRange.dThetaPlus << ", " << mapRange.dThetaMinus << ", " << mapRange.dPhiPlus
588  << ", " << mapRange.dPhiMinus << ", " << dR;
589  return mapRange;
590 }
591 
592 void TrackDetectorAssociator::getTAMuonChamberMatches(std::vector<TAMuonChamberMatch>& matches,
594  // Strategy:
595  // Propagate through the whole detector, estimate change in eta and phi
596  // along the trajectory, add this to dRMuon and find DetIds around this
597  // direction using the map. Then propagate fast to each surface and apply
598  // final matching criteria.
599 
600  // get the direction first
601  SteppingHelixStateInfo trajectoryPoint = cachedTrajectory_.getStateAtHcal();
602  // If trajectory point at HCAL is not valid, try to use the outer most state of the
603  // trajectory instead.
604  if (!trajectoryPoint.isValid())
605  trajectoryPoint = cachedTrajectory_.getOuterState();
606  if (!trajectoryPoint.isValid()) {
607  LogTrace("TrackAssociator")
608  << "trajectory position at HCAL is not valid. Assume the track cannot reach muon detectors and skip it";
609  return;
610  }
611 
612  GlobalVector direction = trajectoryPoint.momentum().unit();
613  LogTrace("TrackAssociator") << "muon direction: " << direction
614  << "\n\t and corresponding point: " << trajectoryPoint.position() << "\n";
615 
616  DetIdAssociator::MapRange mapRange =
617  getMapRange(cachedTrajectory_.trajectoryDelta(CachedTrajectory::FullTrajectory), parameters.dRMuonPreselection);
618 
619  // and find chamber DetIds
620 
621  std::set<DetId> muonIdsInRegion = muonDetIdAssociator_->getDetIdsCloseToAPoint(trajectoryPoint.position(), mapRange);
622  LogTrace("TrackAssociator") << "Number of chambers to check: " << muonIdsInRegion.size();
623  for (std::set<DetId>::const_iterator detId = muonIdsInRegion.begin(); detId != muonIdsInRegion.end(); detId++) {
624  const GeomDet* geomDet = muonDetIdAssociator_->getGeomDet(*detId);
625  TrajectoryStateOnSurface stateOnSurface = cachedTrajectory_.propagate(&geomDet->surface());
626  if (!stateOnSurface.isValid()) {
627  LogTrace("TrackAssociator") << "Failed to propagate the track; moving on\n\t"
628  << "Element is not crosssed: " << DetIdInfo::info(*detId, nullptr) << "\n";
629  continue;
630  }
631  LocalPoint localPoint = geomDet->surface().toLocal(stateOnSurface.freeState()->position());
632  LocalError localError = stateOnSurface.localError().positionError();
633  float distanceX = 0.f;
634  float distanceY = 0.f;
635  if (const CSCChamber* cscChamber = dynamic_cast<const CSCChamber*>(geomDet)) {
636  const CSCChamberSpecs* chamberSpecs = cscChamber->specs();
637  if (!chamberSpecs) {
638  LogTrace("TrackAssociator") << "Failed to get CSCChamberSpecs from CSCChamber; moving on\n";
639  continue;
640  }
641  const CSCLayerGeometry* layerGeometry = chamberSpecs->oddLayerGeometry(1);
642  if (!layerGeometry) {
643  LogTrace("TrackAssociator") << "Failed to get CSCLayerGeometry from CSCChamberSpecs; moving on\n";
644  continue;
645  }
646  const CSCWireTopology* wireTopology = layerGeometry->wireTopology();
647  if (!wireTopology) {
648  LogTrace("TrackAssociator") << "Failed to get CSCWireTopology from CSCLayerGeometry; moving on\n";
649  continue;
650  }
651 
652  float wideWidth = wireTopology->wideWidthOfPlane();
653  float narrowWidth = wireTopology->narrowWidthOfPlane();
654  float length = wireTopology->lengthOfPlane();
655  // If slanted, there is no y offset between local origin and symmetry center of wire plane
656  float yOfFirstWire = std::abs(wireTopology->wireAngle()) > 1.E-06f ? -0.5 * length : wireTopology->yOfWire(1);
657  // y offset between local origin and symmetry center of wire plane
658  float yCOWPOffset = yOfFirstWire + 0.5f * length;
659 
660  // tangent of the incline angle from inside the trapezoid
661  float tangent = (wideWidth - narrowWidth) / (2.f * length);
662  // y position wrt bottom of trapezoid
663  float yPrime = localPoint.y() + std::abs(yOfFirstWire);
664  // half trapezoid width at y' is 0.5 * narrowWidth + x side of triangle with the above tangent and side y'
665  float halfWidthAtYPrime = 0.5f * narrowWidth + yPrime * tangent;
666  distanceX = std::abs(localPoint.x()) - halfWidthAtYPrime;
667  distanceY = std::abs(localPoint.y() - yCOWPOffset) - 0.5f * length;
668  } else {
669  distanceX = std::abs(localPoint.x()) - 0.5f * geomDet->surface().bounds().width();
670  distanceY = std::abs(localPoint.y()) - 0.5f * geomDet->surface().bounds().length();
671  }
672  if ((distanceX < parameters.muonMaxDistanceX && distanceY < parameters.muonMaxDistanceY) ||
673  (distanceX * distanceX <
674  localError.xx() * parameters.muonMaxDistanceSigmaX * parameters.muonMaxDistanceSigmaX &&
675  distanceY * distanceY <
676  localError.yy() * parameters.muonMaxDistanceSigmaY * parameters.muonMaxDistanceSigmaY)) {
677  LogTrace("TrackAssociator") << "found a match: " << DetIdInfo::info(*detId, nullptr) << "\n";
679  match.tState = stateOnSurface;
680  match.localDistanceX = distanceX;
681  match.localDistanceY = distanceY;
682  match.id = *detId;
683  matches.push_back(match);
684  } else {
685  LogTrace("TrackAssociator") << "chamber is too far: " << DetIdInfo::info(*detId, nullptr)
686  << "\n\tdistanceX: " << distanceX << "\t distanceY: " << distanceY
687  << "\t sigmaX: " << distanceX / sqrt(localError.xx())
688  << "\t sigmaY: " << distanceY / sqrt(localError.yy()) << "\n";
689  }
690  }
691 }
692 
696  // Get the segments from the event
698  iEvent.getByToken(parameters.dtSegmentsToken, dtSegments);
699  if (!dtSegments.isValid())
700  throw cms::Exception("FatalError") << "Unable to find DTRecSegment4DCollection in event!\n";
701 
703  iEvent.getByToken(parameters.cscSegmentsToken, cscSegments);
704  if (!cscSegments.isValid())
705  throw cms::Exception("FatalError") << "Unable to find CSCSegmentCollection in event!\n";
706 
708  if (parameters.useGEM)
709  iEvent.getByToken(parameters.gemSegmentsToken, gemSegments);
711  if (parameters.useME0)
712  iEvent.getByToken(parameters.me0SegmentsToken, me0Segments);
713 
715 
716  // check the map of available segments
717  // if there is no segments in a given direction at all,
718  // then there is no point to fly there.
719  //
720  // MISSING
721  // Possible solution: quick search for presence of segments
722  // for the set of DetIds
723 
724  // get a set of matches corresponding to muon chambers
725  std::vector<TAMuonChamberMatch> matchedChambers;
726  LogTrace("TrackAssociator") << "Trying to Get ChamberMatches" << std::endl;
727  getTAMuonChamberMatches(matchedChambers, parameters);
728  LogTrace("TrackAssociator") << "Chambers matched: " << matchedChambers.size() << "\n";
729 
730  // Iterate over all chamber matches and fill segment matching
731  // info if it's available
732  for (std::vector<TAMuonChamberMatch>::iterator matchedChamber = matchedChambers.begin();
733  matchedChamber != matchedChambers.end();
734  matchedChamber++) {
735  const GeomDet* geomDet = muonDetIdAssociator_->getGeomDet((*matchedChamber).id);
736  // DT chamber
737  if (const DTChamber* chamber = dynamic_cast<const DTChamber*>(geomDet)) {
738  // Get the range for the corresponding segments
739  DTRecSegment4DCollection::range range = dtSegments->get(chamber->id());
740  // Loop over the segments of this chamber
741  for (DTRecSegment4DCollection::const_iterator segment = range.first; segment != range.second; segment++) {
742  if (addTAMuonSegmentMatch(*matchedChamber, &(*segment), parameters)) {
743  matchedChamber->segments.back().dtSegmentRef = DTRecSegment4DRef(dtSegments, segment - dtSegments->begin());
744  }
745  }
746  }
747  // CSC Chamber
748  else if (const CSCChamber* chamber = dynamic_cast<const CSCChamber*>(geomDet)) {
749  // Get the range for the corresponding segments
750  CSCSegmentCollection::range range = cscSegments->get(chamber->id());
751  // Loop over the segments
752  for (CSCSegmentCollection::const_iterator segment = range.first; segment != range.second; segment++) {
753  if (addTAMuonSegmentMatch(*matchedChamber, &(*segment), parameters)) {
754  matchedChamber->segments.back().cscSegmentRef = CSCSegmentRef(cscSegments, segment - cscSegments->begin());
755  }
756  }
757  } else {
758  // GEM Chamber
759  if (parameters.useGEM) {
760  if (const GEMSuperChamber* chamber = dynamic_cast<const GEMSuperChamber*>(geomDet)) {
761  // Get the range for the corresponding segments
762  GEMSegmentCollection::range range = gemSegments->get(chamber->id());
763  // Loop over the segments
764  for (GEMSegmentCollection::const_iterator segment = range.first; segment != range.second; segment++) {
765  if (addTAMuonSegmentMatch(*matchedChamber, &(*segment), parameters)) {
766  matchedChamber->segments.back().gemSegmentRef =
767  GEMSegmentRef(gemSegments, segment - gemSegments->begin());
768  }
769  }
770  }
771  }
772  // ME0 Chamber
773  if (parameters.useME0) {
774  if (const ME0Chamber* chamber = dynamic_cast<const ME0Chamber*>(geomDet)) {
775  // Get the range for the corresponding segments
776  ME0SegmentCollection::range range = me0Segments->get(chamber->id());
777  // Loop over the segments
778  for (ME0SegmentCollection::const_iterator segment = range.first; segment != range.second; segment++) {
779  if (addTAMuonSegmentMatch(*matchedChamber, &(*segment), parameters)) {
780  matchedChamber->segments.back().me0SegmentRef =
781  ME0SegmentRef(me0Segments, segment - me0Segments->begin());
782  }
783  }
784  }
785  }
786  }
787  info.chambers.push_back(*matchedChamber);
788  }
789 }
790 
792  const RecSegment* segment,
794  LogTrace("TrackAssociator") << "Segment local position: " << segment->localPosition() << "\n"
795  << std::hex << segment->geographicalId().rawId() << "\n";
796 
797  const GeomDet* chamber = muonDetIdAssociator_->getGeomDet(matchedChamber.id);
798  TrajectoryStateOnSurface trajectoryStateOnSurface = matchedChamber.tState;
799  GlobalPoint segmentGlobalPosition = chamber->toGlobal(segment->localPosition());
800 
801  LogTrace("TrackAssociator") << "Segment global position: " << segmentGlobalPosition
802  << " \t (R_xy,eta,phi): " << segmentGlobalPosition.perp() << ","
803  << segmentGlobalPosition.eta() << "," << segmentGlobalPosition.phi() << "\n";
804 
805  LogTrace("TrackAssociator") << "\teta hit: " << segmentGlobalPosition.eta()
806  << " \tpropagator: " << trajectoryStateOnSurface.freeState()->position().eta() << "\n"
807  << "\tphi hit: " << segmentGlobalPosition.phi()
808  << " \tpropagator: " << trajectoryStateOnSurface.freeState()->position().phi()
809  << std::endl;
810 
811  bool isGood = false;
812  bool isDTWithoutY = false;
813  const DTRecSegment4D* dtseg = dynamic_cast<const DTRecSegment4D*>(segment);
814  if (dtseg && (!dtseg->hasZed()))
815  isDTWithoutY = true;
816 
817  float deltaPhi(std::abs(segmentGlobalPosition.phi() - trajectoryStateOnSurface.freeState()->position().phi()));
818  if (deltaPhi > float(M_PI))
819  deltaPhi = std::abs(deltaPhi - float(M_PI) * 2.f);
820  float deltaEta = std::abs(segmentGlobalPosition.eta() - trajectoryStateOnSurface.freeState()->position().eta());
821 
822  if (isDTWithoutY) {
823  isGood = deltaPhi < parameters.dRMuon;
824  // Be in chamber
825  isGood &= deltaEta < .3f;
826  } else
827  isGood = deltaEta * deltaEta + deltaPhi * deltaPhi < parameters.dRMuon * parameters.dRMuon;
828 
829  if (isGood) {
830  TAMuonSegmentMatch muonSegment;
831  muonSegment.segmentGlobalPosition = getPoint(segmentGlobalPosition);
832  muonSegment.segmentLocalPosition = getPoint(segment->localPosition());
833  muonSegment.segmentLocalDirection = getVector(segment->localDirection());
834  muonSegment.segmentLocalErrorXX = segment->localPositionError().xx();
835  muonSegment.segmentLocalErrorYY = segment->localPositionError().yy();
836  muonSegment.segmentLocalErrorXY = segment->localPositionError().xy();
837  muonSegment.segmentLocalErrorDxDz = segment->localDirectionError().xx();
838  muonSegment.segmentLocalErrorDyDz = segment->localDirectionError().yy();
839 
840  // DANGEROUS - compiler cannot guaranty parameters ordering
841  // AlgebraicSymMatrix segmentCovMatrix = segment->parametersError();
842  // muonSegment.segmentLocalErrorXDxDz = segmentCovMatrix[2][0];
843  // muonSegment.segmentLocalErrorYDyDz = segmentCovMatrix[3][1];
844  muonSegment.segmentLocalErrorXDxDz = -999.f;
845  muonSegment.segmentLocalErrorYDyDz = -999.f;
846  muonSegment.hasZed = true;
847  muonSegment.hasPhi = true;
848 
849  // timing information
850  muonSegment.t0 = 0.f;
851  if (dtseg) {
852  if ((dtseg->hasPhi()) && (!isDTWithoutY)) {
853  int phiHits = dtseg->phiSegment()->specificRecHits().size();
854  // int zHits = dtseg->zSegment()->specificRecHits().size();
855  int hits = 0;
856  double t0 = 0.;
857  // TODO: cuts on hit numbers not optimized in any way yet...
858  if (phiHits > 5 && dtseg->phiSegment()->ist0Valid()) {
859  t0 += dtseg->phiSegment()->t0() * phiHits;
860  hits += phiHits;
861  LogTrace("TrackAssociator") << " Phi t0: " << dtseg->phiSegment()->t0() << " hits: " << phiHits;
862  }
863  // the z segments seem to contain little useful information...
864  // if (zHits>3) {
865  // t0+=s->zSegment()->t0()*zHits;
866  // hits+=zHits;
867  // LogTrace("TrackAssociator") << " Z t0: " << s->zSegment()->t0() << " hits: " << zHits << std::endl;
868  // }
869  if (hits)
870  muonSegment.t0 = t0 / hits;
871  // LogTrace("TrackAssociator") << " --- t0: " << muonSegment.t0 << std::endl;
872  } else {
873  // check and set dimensionality
874  if (isDTWithoutY)
875  muonSegment.hasZed = false;
876  if (!dtseg->hasPhi())
877  muonSegment.hasPhi = false;
878  }
879  }
880  matchedChamber.segments.push_back(muonSegment);
881  }
882 
883  return isGood;
884 }
885 
886 //********************** NON-CORE CODE ******************************//
887 
891  // get list of simulated tracks and their vertices
892  using namespace edm;
893  Handle<SimTrackContainer> simTracks;
894  iEvent.getByToken(parameters.simTracksToken, simTracks);
895  if (!simTracks.isValid())
896  throw cms::Exception("FatalError") << "No simulated tracks found\n";
897 
898  Handle<SimVertexContainer> simVertices;
899  iEvent.getByToken(parameters.simVerticesToken, simVertices);
900  if (!simVertices.isValid())
901  throw cms::Exception("FatalError") << "No simulated vertices found\n";
902 
903  // get sim calo hits
904  Handle<PCaloHitContainer> simEcalHitsEB;
905  iEvent.getByToken(parameters.simEcalHitsEBToken, simEcalHitsEB);
906  if (!simEcalHitsEB.isValid())
907  throw cms::Exception("FatalError") << "No simulated ECAL EB hits found\n";
908 
909  Handle<PCaloHitContainer> simEcalHitsEE;
910  iEvent.getByToken(parameters.simEcalHitsEEToken, simEcalHitsEE);
911  if (!simEcalHitsEE.isValid())
912  throw cms::Exception("FatalError") << "No simulated ECAL EE hits found\n";
913 
914  Handle<PCaloHitContainer> simHcalHits;
915  iEvent.getByToken(parameters.simHcalHitsToken, simHcalHits);
916  if (!simHcalHits.isValid())
917  throw cms::Exception("FatalError") << "No simulated HCAL hits found\n";
918 
919  // find truth partner
920  SimTrackContainer::const_iterator simTrack = simTracks->begin();
921  for (; simTrack != simTracks->end(); ++simTrack) {
922  math::XYZVector simP3(simTrack->momentum().x(), simTrack->momentum().y(), simTrack->momentum().z());
923  math::XYZVector recoP3(info.stateAtIP.momentum().x(), info.stateAtIP.momentum().y(), info.stateAtIP.momentum().z());
924  if (ROOT::Math::VectorUtil::DeltaR(recoP3, simP3) < 0.1)
925  break;
926  }
927  if (simTrack != simTracks->end()) {
928  info.simTrack = &(*simTrack);
929  float ecalTrueEnergy(0);
930  float hcalTrueEnergy(0);
931 
932  // loop over calo hits
933  for (PCaloHitContainer::const_iterator hit = simEcalHitsEB->begin(); hit != simEcalHitsEB->end(); ++hit)
934  if (hit->geantTrackId() == info.simTrack->genpartIndex())
935  ecalTrueEnergy += hit->energy();
936 
937  for (PCaloHitContainer::const_iterator hit = simEcalHitsEE->begin(); hit != simEcalHitsEE->end(); ++hit)
938  if (hit->geantTrackId() == info.simTrack->genpartIndex())
939  ecalTrueEnergy += hit->energy();
940 
941  for (PCaloHitContainer::const_iterator hit = simHcalHits->begin(); hit != simHcalHits->end(); ++hit)
942  if (hit->geantTrackId() == info.simTrack->genpartIndex())
943  hcalTrueEnergy += hit->energy();
944 
945  info.ecalTrueEnergy = ecalTrueEnergy;
946  info.hcalTrueEnergy = hcalTrueEnergy;
947  info.hcalTrueEnergyCorrected = hcalTrueEnergy;
948  if (std::abs(info.trkGlobPosAtHcal.eta()) < 1.3f)
949  info.hcalTrueEnergyCorrected = hcalTrueEnergy * 113.2f;
950  else if (std::abs(info.trkGlobPosAtHcal.eta()) < 3.0f)
951  info.hcalTrueEnergyCorrected = hcalTrueEnergy * 167.2f;
952  }
953 }
954 
956  const edm::EventSetup& iSetup,
957  const reco::Track& track,
959  Direction direction /*= Any*/) {
960  double currentStepSize = cachedTrajectory_.getPropagationStep();
961 
962  const MagneticField* bField = &iSetup.getData(parameters.bFieldToken);
963 
964  if (track.extra().isNull()) {
965  if (direction != InsideOut)
966  throw cms::Exception("FatalError") << "No TrackExtra information is available and association is done with "
967  "something else than InsideOut track.\n"
968  << "Either change the parameter or provide needed data!\n";
969  LogTrace("TrackAssociator") << "Track Extras not found\n";
971  return associate(iEvent, iSetup, parameters, &initialState); // 5th argument is null pointer
972  }
973 
974  LogTrace("TrackAssociator") << "Track Extras found\n";
977  FreeTrajectoryState referenceState = trajectoryStateTransform::initialFreeState(track, bField);
978 
979  LogTrace("TrackAssociator") << "inner track state (rho, z, phi):" << track.innerPosition().Rho() << ", "
980  << track.innerPosition().z() << ", " << track.innerPosition().phi() << "\n";
981  LogTrace("TrackAssociator") << "innerFreeState (rho, z, phi):" << innerState.position().perp() << ", "
982  << innerState.position().z() << ", " << innerState.position().phi() << "\n";
983 
984  LogTrace("TrackAssociator") << "outer track state (rho, z, phi):" << track.outerPosition().Rho() << ", "
985  << track.outerPosition().z() << ", " << track.outerPosition().phi() << "\n";
986  LogTrace("TrackAssociator") << "outerFreeState (rho, z, phi):" << outerState.position().perp() << ", "
987  << outerState.position().z() << ", " << outerState.position().phi() << "\n";
988 
989  // InsideOut first
990  if (crossedIP(track)) {
991  switch (direction) {
992  case InsideOut:
993  case Any:
994  return associate(iEvent, iSetup, parameters, &referenceState, &outerState);
995  break;
996  case OutsideIn: {
997  cachedTrajectory_.setPropagationStep(-std::abs(currentStepSize));
998  TrackDetMatchInfo result = associate(iEvent, iSetup, parameters, &innerState, &referenceState);
999  cachedTrajectory_.setPropagationStep(currentStepSize);
1000  return result;
1001  break;
1002  }
1003  }
1004  } else {
1005  switch (direction) {
1006  case InsideOut:
1007  return associate(iEvent, iSetup, parameters, &innerState, &outerState);
1008  break;
1009  case OutsideIn: {
1010  cachedTrajectory_.setPropagationStep(-std::abs(currentStepSize));
1011  TrackDetMatchInfo result = associate(iEvent, iSetup, parameters, &outerState, &innerState);
1012  cachedTrajectory_.setPropagationStep(currentStepSize);
1013  return result;
1014  break;
1015  }
1016  case Any: {
1017  // check if we deal with clear outside-in case
1018  if (track.innerPosition().Dot(track.innerMomentum()) < 0 &&
1019  track.outerPosition().Dot(track.outerMomentum()) < 0) {
1020  cachedTrajectory_.setPropagationStep(-std::abs(currentStepSize));
1022  if (track.innerPosition().Mag2() < track.outerPosition().Mag2())
1023  result = associate(iEvent, iSetup, parameters, &innerState, &outerState);
1024  else
1025  result = associate(iEvent, iSetup, parameters, &outerState, &innerState);
1026  cachedTrajectory_.setPropagationStep(currentStepSize);
1027  return result;
1028  }
1029  }
1030  }
1031  }
1032 
1033  // all other cases
1034  return associate(iEvent, iSetup, parameters, &innerState, &outerState);
1035 }
1036 
1038  const edm::EventSetup& iSetup,
1039  const SimTrack& track,
1040  const SimVertex& vertex,
1042  auto const* bField = &iSetup.getData(parameters.bFieldToken);
1043  return associate(iEvent, iSetup, getFreeTrajectoryState(bField, track, vertex), parameters);
1044 }
1045 
1047  const edm::EventSetup& iSetup,
1048  const GlobalVector& momentum,
1049  const GlobalPoint& vertex,
1050  const int charge,
1052  auto const* bField = &iSetup.getData(parameters.bFieldToken);
1053  return associate(iEvent, iSetup, getFreeTrajectoryState(bField, momentum, vertex, charge), parameters);
1054 }
1055 
1057  bool crossed = true;
1058  crossed &= (track.innerPosition().rho() > 3); // something close to active volume
1059  crossed &= (track.outerPosition().rho() > 3); // something close to active volume
1060  crossed &=
1061  ((track.innerPosition().x() * track.innerMomentum().x() + track.innerPosition().y() * track.innerMomentum().y() <
1062  0) !=
1063  (track.outerPosition().x() * track.outerMomentum().x() + track.outerPosition().y() * track.outerMomentum().y() <
1064  0));
1065  return crossed;
1066 }
edm::ESGetToken< DetIdAssociator, DetIdAssociatorRecord > hoDetIdAssociatorToken
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
edm::EDGetTokenT< CaloTowerCollection > caloTowersToken
std::vector< DetId > crossedPreshowerIds
math::XYZPoint trkGlobPosAtHO
edm::ESGetToken< DetIdAssociator, DetIdAssociatorRecord > caloDetIdAssociatorToken
float xx() const
Definition: LocalError.h:22
static const TGPicture * info(bool iBackgroundIsBlack)
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:50
const CSCWireTopology * wireTopology() const
virtual float length() const =0
T perp() const
Definition: PV3DBase.h:69
static bool crossedIP(const reco::Track &track)
std::vector< const CaloTower * > crossedTowers
virtual LocalError localDirectionError() const =0
Error on the local direction.
std::vector< const CaloTower * > towers
uint16_t *__restrict__ id
TrajectoryStateOnSurface tState
const TrackExtraRef & extra() const
reference to &quot;extra&quot; object
Definition: Track.h:139
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventIDconst &, edm::Timestampconst & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
std::vector< const HBHERecHit * > crossedHcalRecHits
std::vector< DetId > crossedTowerIds
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
math::XYZPoint segmentGlobalPosition
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
edm::Ref< CSCSegmentCollection > CSCSegmentRef
edm::ESGetToken< DetIdAssociator, DetIdAssociatorRecord > preshowerDetIdAssociatorToken
const DTChamberRecSegment2D * phiSegment() const
The superPhi segment: 0 if no phi projection available.
std::vector< const EcalRecHit * > ecalRecHits
hits in the cone
std::vector< DetId > crossedEcalIds
edm::Ref< ME0SegmentCollection > ME0SegmentRef
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49
math::XYZPoint segmentLocalPosition
virtual LocalVector localDirection() const =0
Local direction.
edm::EDGetTokenT< GEMSegmentCollection > gemSegmentsToken
int init
Definition: HydjetWrapper.h:64
void useDefaultPropagator()
use the default propagator
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
std::vector< CaloTower >::const_iterator const_iterator
T y() const
Definition: PV3DBase.h:60
void fillHcal(const edm::Event &, TrackDetMatchInfo &, const AssociatorParameters &)
const Bounds & bounds() const
Definition: Surface.h:87
bool addTAMuonSegmentMatch(TAMuonChamberMatch &, const RecSegment *, const AssociatorParameters &)
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:664
const math::XYZPoint & outerPosition() const
position of the outermost hit
Definition: Track.h:62
std::vector< DetId > crossedHcalIds
GlobalVector momentum() const
std::vector< const EcalRecHit * > crossedEcalRecHits
hits in detector elements crossed by a track
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
math::XYZVector segmentLocalDirection
LocalError positionError() const
edm::EDGetTokenT< EBRecHitCollection > EBRecHitsToken
void setCaloGeometry(const CaloGeometry *geometry)
static const double deltaEta
Definition: CaloConstants.h:8
#define LogTrace(id)
math::XYZPoint trkGlobPosAtHcal
const uint16_t range(const Frame &aFrame)
tuple result
Definition: mps_fire.py:311
FreeTrajectoryState stateAtIP
track info
void fillCaloTowers(const edm::Event &, TrackDetMatchInfo &, const AssociatorParameters &)
edm::ESGetToken< DetIdAssociator, DetIdAssociatorRecord > ecalDetIdAssociatorToken
const math::XYZPoint & innerPosition() const
position of the innermost hit
Definition: Track.h:56
void fillCaloTruth(const edm::Event &, TrackDetMatchInfo &, const AssociatorParameters &)
bool getData(T &iHolder) const
Definition: EventSetup.h:122
edm::EDGetTokenT< edm::PCaloHitContainer > simEcalHitsEBToken
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
edm::EDGetTokenT< HORecHitCollection > HOcollToken
float xy() const
Definition: LocalError.h:23
const Point & vertex() const
reference point on the track. This method is DEPRECATED, please use referencePoint() instead ...
Definition: TrackBase.h:676
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > theCaloGeometryToken
void init(const edm::EventSetup &, const AssociatorParameters &)
float wireAngle() const override
int iEvent
Definition: GenABIO.cc:224
T mag() const
Definition: PV3DBase.h:64
GlobalVector trkMomAtHO
void setPropagator(const Propagator *)
use a user configured propagator
void getTAMuonChamberMatches(std::vector< TAMuonChamberMatch > &matches, const AssociatorParameters &parameters)
float yy() const
Definition: LocalError.h:24
GlobalPoint position() const
T sqrt(T t)
Definition: SSEVec.h:19
LocalPoint toLocal(const GlobalPoint &gp) const
edm::Ref< GEMSegmentCollection > GEMSegmentRef
FreeTrajectoryState const * freeState(bool withErrors=true) const
T z() const
Definition: PV3DBase.h:61
std::vector< TAMuonSegmentMatch > segments
distance sign convention: negative - crossed chamber, positive - missed chamber
int genpartIndex() const
index of the corresponding Generator particle in the Event container (-1 if no Genpart) ...
Definition: SimTrack.h:37
edm::EDGetTokenT< ME0SegmentCollection > me0SegmentsToken
std::vector< TAMuonChamberMatch > chambers
def move
Definition: eostools.py:511
std::vector< const HBHERecHit * > hcalRecHits
edm::ESGetToken< GlobalTrackingGeometry, GlobalTrackingGeometryRecord > theTrackingGeometryToken
edm::EDGetTokenT< DTRecSegment4DCollection > dtSegmentsToken
std::vector< DetId > crossedHOIds
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > bFieldToken
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const math::XYZTLorentzVectorD & position() const
Definition: CoreSimVertex.h:21
edm::ESGetToken< DetIdAssociator, DetIdAssociatorRecord > muonDetIdAssociatorToken
bool hasPhi() const
Does it have the Phi projection?
double lengthOfPlane() const
const LocalTrajectoryError & localError() const
edm::EDGetTokenT< edm::SimTrackContainer > simTracksToken
float yOfWire(float wire, float x=0.) const
edm::EDGetTokenT< EERecHitCollection > EERecHitsToken
bool isValid() const
Definition: HandleBase.h:70
edm::EDGetTokenT< HBHERecHitCollection > HBHEcollToken
void fillPreshower(const edm::Event &iEvent, TrackDetMatchInfo &info, const AssociatorParameters &)
bool isNull() const
Checks for null.
Definition: Ref.h:235
GlobalVector momentum() const
std::vector< DTRecHit1D > specificRecHits() const
Access to specific components.
bool ist0Valid() const
double narrowWidthOfPlane() const
#define M_PI
edm::Ref< DTRecSegment4DCollection > DTRecSegment4DRef
ROOT::Math::SMatrix< double, 6, 6, ROOT::Math::MatRepSym< double, 6 > > AlgebraicSymMatrix66
Vector3DBase unit() const
Definition: Vector3DBase.h:54
bool hasZed() const
Does it have the Z projection?
GlobalPoint position() const
const math::XYZVector & outerMomentum() const
momentum vector at the outermost hit position
Definition: Track.h:65
double wideWidthOfPlane() const
GlobalVector trkMomAtEcal
edm::EDGetTokenT< edm::SimVertexContainer > simVerticesToken
static std::string info(const DetId &, const TrackerTopology *tTopo)
Definition: DetIdInfo.cc:25
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
virtual LocalError localPositionError() const =0
void fillMuon(const edm::Event &, TrackDetMatchInfo &, const AssociatorParameters &)
static FreeTrajectoryState getFreeTrajectoryState(const MagneticField *, const reco::Track &)
get FreeTrajectoryState from different track representations
const CSCLayerGeometry * oddLayerGeometry(int iendcap) const
Accessors for LayerGeometry&#39;s.
DetIdAssociator::MapRange getMapRange(const std::pair< float, float > &delta, const float dR)
T eta() const
Definition: PV3DBase.h:73
GlobalVector trkMomAtHcal
int type() const
particle type (HEP PDT convension)
Definition: CoreSimTrack.h:22
void fillEcal(const edm::Event &, TrackDetMatchInfo &, const AssociatorParameters &)
const math::XYZTLorentzVectorD & momentum() const
Definition: CoreSimTrack.h:19
const math::XYZVector & innerMomentum() const
momentum vector at the innermost hit position
Definition: Track.h:59
std::vector< const HORecHit * > crossedHORecHits
TrackDetMatchInfo associate(const edm::Event &, const edm::EventSetup &, const FreeTrajectoryState &, const AssociatorParameters &)
math::XYZPoint trkGlobPosAtEcal
Track position at different parts of the calorimeter.
edm::EDGetTokenT< edm::PCaloHitContainer > simHcalHitsToken
int charge() const
track electric charge
Definition: TrackBase.h:596
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
DetId geographicalId() const
FreeTrajectoryState innerFreeState(const reco::Track &tk, const MagneticField *field, bool withErr=true)
double t0() const
Get the segment t0 (if recomputed, 0 is returned otherwise)
std::vector< const HORecHit * > hoRecHits
edm::ESGetToken< DetIdAssociator, DetIdAssociatorRecord > hcalDetIdAssociatorToken
T x() const
Definition: PV3DBase.h:59
virtual LocalPoint localPosition() const =0
virtual float width() const =0
FreeTrajectoryState outerFreeState(const reco::Track &tk, const MagneticField *field, bool withErr=true)
const SimTrack * simTrack
MC truth info.
FreeTrajectoryState initialFreeState(const reco::Track &tk, const MagneticField *field, bool withErr=true)
edm::EDGetTokenT< CSCSegmentCollection > cscSegmentsToken
void fillHO(const edm::Event &, TrackDetMatchInfo &, const AssociatorParameters &)
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
edm::EDGetTokenT< edm::PCaloHitContainer > simEcalHitsEEToken