CMS 3D CMS Logo

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