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  if (parameters.useMuon) {
207  double maxR = muonDetIdAssociator_->volume().maxR();
208  double maxZ = muonDetIdAssociator_->volume().maxZ();
209  cachedTrajectory_.setMaxDetectorRadius(maxR);
210  cachedTrajectory_.setMaxDetectorLength(maxZ*2.);
211  }
212  else {
213  cachedTrajectory_.setMaxDetectorRadius(HOmaxR);
214  cachedTrajectory_.setMaxDetectorLength(HOmaxZ*2.);
215  }
216 
217  // If track extras exist and outerState is before HO maximum, then use outerState
218  if (outerState) {
219  if (outerState->position().perp()<HOmaxR && fabs(outerState->position().z())<HOmaxZ) {
220  LogTrace("TrackAssociator") << "Using outerState as trackOrigin at Rho=" << outerState->position().perp()
221  << " Z=" << outerState->position().z() << "\n";
222  trackOrigin = SteppingHelixStateInfo(*outerState);
223  }
224  else if(innerState) {
225  LogTrace("TrackAssociator") << "Using innerState as trackOrigin at Rho=" << innerState->position().perp()
226  << " Z=" << innerState->position().z() << "\n";
227  trackOrigin = SteppingHelixStateInfo(*innerState);
228  }
229  }
230 
231  if ( trackOrigin.momentum().mag() == 0 ) return info;
232  if ( std::isnan(trackOrigin.momentum().x()) or std::isnan(trackOrigin.momentum().y()) or std::isnan(trackOrigin.momentum().z()) ) return info;
233  if ( ! cachedTrajectory_.propagateAll(trackOrigin) ) return info;
234 
235  // get trajectory in calorimeters
236  cachedTrajectory_.findEcalTrajectory( ecalDetIdAssociator_->volume() );
237  cachedTrajectory_.findHcalTrajectory( hcalDetIdAssociator_->volume() );
238  cachedTrajectory_.findHOTrajectory( hoDetIdAssociator_->volume() );
239  cachedTrajectory_.findPreshowerTrajectory( preshowerDetIdAssociator_->volume() );
240 
241  info.trkGlobPosAtEcal = getPoint( cachedTrajectory_.getStateAtEcal().position() );
242  info.trkGlobPosAtHcal = getPoint( cachedTrajectory_.getStateAtHcal().position() );
243  info.trkGlobPosAtHO = getPoint( cachedTrajectory_.getStateAtHO().position() );
244 
245  info.trkMomAtEcal = cachedTrajectory_.getStateAtEcal().momentum();
246  info.trkMomAtHcal = cachedTrajectory_.getStateAtHcal().momentum();
247  info.trkMomAtHO = cachedTrajectory_.getStateAtHO().momentum();
248 
249  if (parameters.useEcal) fillEcal( iEvent, info, parameters);
250  if (parameters.useCalo) fillCaloTowers( iEvent, info, parameters);
251  if (parameters.useHcal) fillHcal( iEvent, info, parameters);
252  if (parameters.useHO) fillHO( iEvent, info, parameters);
253  if (parameters.usePreshower) fillPreshower( iEvent, info, parameters);
254  if (parameters.useMuon) fillMuon( iEvent, info, parameters);
255  if (parameters.truthMatch) fillCaloTruth( iEvent, info, parameters);
256 
257  return info;
258 }
259 
263 {
264  const std::vector<SteppingHelixStateInfo>& trajectoryStates = cachedTrajectory_.getEcalTrajectory();
265 
266  for(std::vector<SteppingHelixStateInfo>::const_iterator itr = trajectoryStates.begin();
267  itr != trajectoryStates.end(); itr++)
268  LogTrace("TrackAssociator") << "ECAL trajectory point (rho, z, phi): " << itr->position().perp() <<
269  ", " << itr->position().z() << ", " << itr->position().phi();
270 
271  std::vector<GlobalPoint> coreTrajectory;
272  for(std::vector<SteppingHelixStateInfo>::const_iterator itr = trajectoryStates.begin();
273  itr != trajectoryStates.end(); itr++) coreTrajectory.push_back(itr->position());
274 
275  if(coreTrajectory.empty()) {
276  LogTrace("TrackAssociator") << "ECAL track trajectory is empty; moving on\n";
277  info.isGoodEcal = 0;
278  return;
279  }
280  info.isGoodEcal = 1;
281 
282  // Find ECAL crystals
284  iEvent.getByToken(parameters.EBRecHitsToken, EBRecHits);
285  if (!EBRecHits.isValid()) throw cms::Exception("FatalError") << "Unable to find EBRecHitCollection in the event!\n";
286 
288  iEvent.getByToken(parameters.EERecHitsToken, EERecHits);
289  if (!EERecHits.isValid()) throw cms::Exception("FatalError") << "Unable to find EERecHitCollection in event!\n";
290 
291  std::set<DetId> ecalIdsInRegion;
292  if (parameters.accountForTrajectoryChangeCalo){
293  // get trajectory change with respect to initial state
294  DetIdAssociator::MapRange mapRange = getMapRange(cachedTrajectory_.trajectoryDelta(CachedTrajectory::IpToEcal),
295  parameters.dREcalPreselection);
296  ecalIdsInRegion = ecalDetIdAssociator_->getDetIdsCloseToAPoint(coreTrajectory[0],mapRange);
297  } else ecalIdsInRegion = ecalDetIdAssociator_->getDetIdsCloseToAPoint(coreTrajectory[0], parameters.dREcalPreselection);
298  LogTrace("TrackAssociator") << "ECAL hits in the region: " << ecalIdsInRegion.size();
299  if (parameters.dREcalPreselection > parameters.dREcal)
300  ecalIdsInRegion = ecalDetIdAssociator_->getDetIdsInACone(ecalIdsInRegion, coreTrajectory, parameters.dREcal);
301  LogTrace("TrackAssociator") << "ECAL hits in the cone: " << ecalIdsInRegion.size();
302  std::vector<DetId> crossedEcalIds =
303  ecalDetIdAssociator_->getCrossedDetIds(ecalIdsInRegion, coreTrajectory);
304  LogTrace("TrackAssociator") << "ECAL crossed hits " << crossedEcalIds.size();
305 
306  info.crossedEcalIds = crossedEcalIds;
307 
308  // add EcalRecHits
309  for(std::vector<DetId>::const_iterator itr=crossedEcalIds.begin(); itr!=crossedEcalIds.end();itr++)
310  {
311  std::vector<EcalRecHit>::const_iterator ebHit = (*EBRecHits).find(*itr);
312  std::vector<EcalRecHit>::const_iterator eeHit = (*EERecHits).find(*itr);
313  if(ebHit != (*EBRecHits).end())
314  info.crossedEcalRecHits.push_back(&*ebHit);
315  else if(eeHit != (*EERecHits).end())
316  info.crossedEcalRecHits.push_back(&*eeHit);
317  else
318  LogTrace("TrackAssociator") << "Crossed EcalRecHit is not found for DetId: " << itr->rawId();
319  }
320  for(std::set<DetId>::const_iterator itr=ecalIdsInRegion.begin(); itr!=ecalIdsInRegion.end();itr++)
321  {
322  std::vector<EcalRecHit>::const_iterator ebHit = (*EBRecHits).find(*itr);
323  std::vector<EcalRecHit>::const_iterator eeHit = (*EERecHits).find(*itr);
324  if(ebHit != (*EBRecHits).end())
325  info.ecalRecHits.push_back(&*ebHit);
326  else if(eeHit != (*EERecHits).end())
327  info.ecalRecHits.push_back(&*eeHit);
328  else
329  LogTrace("TrackAssociator") << "EcalRecHit from the cone is not found for DetId: " << itr->rawId();
330  }
331 }
332 
336 {
337  // use ECAL and HCAL trajectories to match a tower. (HO isn't used for matching).
338  std::vector<GlobalPoint> trajectory;
339  const std::vector<SteppingHelixStateInfo>& ecalTrajectoryStates = cachedTrajectory_.getEcalTrajectory();
340  const std::vector<SteppingHelixStateInfo>& hcalTrajectoryStates = cachedTrajectory_.getHcalTrajectory();
341  for(std::vector<SteppingHelixStateInfo>::const_iterator itr = ecalTrajectoryStates.begin();
342  itr != ecalTrajectoryStates.end(); itr++) trajectory.push_back(itr->position());
343  for(std::vector<SteppingHelixStateInfo>::const_iterator itr = hcalTrajectoryStates.begin();
344  itr != hcalTrajectoryStates.end(); itr++) trajectory.push_back(itr->position());
345 
346  if(trajectory.empty()) {
347  LogTrace("TrackAssociator") << "HCAL trajectory is empty; moving on\n";
348  info.isGoodCalo = 0;
349  return;
350  }
351  info.isGoodCalo = 1;
352 
353  // find crossed CaloTowers
355  iEvent.getByToken(parameters.caloTowersToken, caloTowers);
356  if (!caloTowers.isValid()) throw cms::Exception("FatalError") << "Unable to find CaloTowers in event!\n";
357 
358  std::set<DetId> caloTowerIdsInRegion;
359  if (parameters.accountForTrajectoryChangeCalo){
360  // get trajectory change with respect to initial state
361  DetIdAssociator::MapRange mapRange = getMapRange(cachedTrajectory_.trajectoryDelta(CachedTrajectory::IpToHcal),
362  parameters.dRHcalPreselection);
363  caloTowerIdsInRegion = caloDetIdAssociator_->getDetIdsCloseToAPoint(trajectory[0],mapRange);
364  } else caloTowerIdsInRegion = caloDetIdAssociator_->getDetIdsCloseToAPoint(trajectory[0], parameters.dRHcalPreselection);
365 
366  LogTrace("TrackAssociator") << "Towers in the region: " << caloTowerIdsInRegion.size();
367  std::set<DetId> caloTowerIdsInACone = caloDetIdAssociator_->getDetIdsInACone(caloTowerIdsInRegion, trajectory, parameters.dRHcal);
368  LogTrace("TrackAssociator") << "Towers in the cone: " << caloTowerIdsInACone.size();
369  std::vector<DetId> crossedCaloTowerIds = caloDetIdAssociator_->getCrossedDetIds(caloTowerIdsInRegion, trajectory);
370  LogTrace("TrackAssociator") << "Towers crossed: " << crossedCaloTowerIds.size();
371 
372  info.crossedTowerIds = crossedCaloTowerIds;
373 
374  // add CaloTowers
375  for(std::vector<DetId>::const_iterator itr=crossedCaloTowerIds.begin(); itr!=crossedCaloTowerIds.end();itr++)
376  {
377  CaloTowerCollection::const_iterator tower = (*caloTowers).find(*itr);
378  if(tower != (*caloTowers).end())
379  info.crossedTowers.push_back(&*tower);
380  else
381  LogTrace("TrackAssociator") << "Crossed CaloTower is not found for DetId: " << (*itr).rawId();
382  }
383 
384  for(std::set<DetId>::const_iterator itr=caloTowerIdsInACone.begin(); itr!=caloTowerIdsInACone.end();itr++)
385  {
386  CaloTowerCollection::const_iterator tower = (*caloTowers).find(*itr);
387  if(tower != (*caloTowers).end())
388  info.towers.push_back(&*tower);
389  else
390  LogTrace("TrackAssociator") << "CaloTower from the cone is not found for DetId: " << (*itr).rawId();
391  }
392 
393 }
394 
398 {
399  std::vector<GlobalPoint> trajectory;
400  const std::vector<SteppingHelixStateInfo>& trajectoryStates = cachedTrajectory_.getPreshowerTrajectory();
401  for(std::vector<SteppingHelixStateInfo>::const_iterator itr = trajectoryStates.begin();
402  itr != trajectoryStates.end(); itr++) trajectory.push_back(itr->position());
403 
404  if(trajectory.empty()) {
405  LogTrace("TrackAssociator") << "Preshower trajectory is empty; moving on\n";
406  return;
407  }
408 
409  std::set<DetId> idsInRegion =
410  preshowerDetIdAssociator_->getDetIdsCloseToAPoint(trajectory[0],
411  parameters.dRPreshowerPreselection);
412 
413  LogTrace("TrackAssociator") << "Number of Preshower Ids in the region: " << idsInRegion.size();
414  std::vector<DetId> crossedIds = preshowerDetIdAssociator_->getCrossedDetIds(idsInRegion, trajectory);
415  LogTrace("TrackAssociator") << "Number of Preshower Ids in crossed: " << crossedIds.size();
416  info.crossedPreshowerIds = crossedIds;
417 }
418 
419 
423 {
424  const std::vector<SteppingHelixStateInfo>& trajectoryStates = cachedTrajectory_.getHcalTrajectory();
425 
426  std::vector<GlobalPoint> coreTrajectory;
427  for(std::vector<SteppingHelixStateInfo>::const_iterator itr = trajectoryStates.begin();
428  itr != trajectoryStates.end(); itr++) coreTrajectory.push_back(itr->position());
429 
430  if(coreTrajectory.empty()) {
431  LogTrace("TrackAssociator") << "HCAL trajectory is empty; moving on\n";
432  info.isGoodHcal = 0;
433  return;
434  }
435  info.isGoodHcal = 1;
436 
437  // find crossed Hcals
439  iEvent.getByToken(parameters.HBHEcollToken, collection);
440  if ( ! collection.isValid() ) throw cms::Exception("FatalError") << "Unable to find HBHERecHits in event!\n";
441 
442  std::set<DetId> idsInRegion;
443  if (parameters.accountForTrajectoryChangeCalo){
444  // get trajectory change with respect to initial state
445  DetIdAssociator::MapRange mapRange = getMapRange(cachedTrajectory_.trajectoryDelta(CachedTrajectory::IpToHcal),
446  parameters.dRHcalPreselection);
447  idsInRegion = hcalDetIdAssociator_->getDetIdsCloseToAPoint(coreTrajectory[0], mapRange);
448  } else idsInRegion = hcalDetIdAssociator_->getDetIdsCloseToAPoint(coreTrajectory[0], parameters.dRHcalPreselection);
449 
450  LogTrace("TrackAssociator") << "HCAL hits in the region: " << idsInRegion.size() << "\n" << DetIdInfo::info(idsInRegion,0);
451  std::set<DetId> idsInACone = hcalDetIdAssociator_->getDetIdsInACone(idsInRegion, coreTrajectory, parameters.dRHcal);
452  LogTrace("TrackAssociator") << "HCAL hits in the cone: " << idsInACone.size() << "\n" << DetIdInfo::info(idsInACone,0);
453  std::vector<DetId> crossedIds =
454  hcalDetIdAssociator_->getCrossedDetIds(idsInRegion, coreTrajectory);
455  LogTrace("TrackAssociator") << "HCAL hits crossed: " << crossedIds.size() << "\n" << DetIdInfo::info(crossedIds,0);
456 
457  info.crossedHcalIds = crossedIds;
458  // add Hcal
459  for(std::vector<DetId>::const_iterator itr=crossedIds.begin(); itr!=crossedIds.end();itr++)
460  {
461  HBHERecHitCollection::const_iterator hit = (*collection).find(*itr);
462  if( hit != (*collection).end() )
463  info.crossedHcalRecHits.push_back(&*hit);
464  else
465  LogTrace("TrackAssociator") << "Crossed HBHERecHit is not found for DetId: " << itr->rawId();
466  }
467  for(std::set<DetId>::const_iterator itr=idsInACone.begin(); itr!=idsInACone.end();itr++)
468  {
469  HBHERecHitCollection::const_iterator hit = (*collection).find(*itr);
470  if( hit != (*collection).end() )
471  info.hcalRecHits.push_back(&*hit);
472  else
473  LogTrace("TrackAssociator") << "HBHERecHit from the cone is not found for DetId: " << itr->rawId();
474  }
475 }
476 
480 {
481  const std::vector<SteppingHelixStateInfo>& trajectoryStates = cachedTrajectory_.getHOTrajectory();
482 
483  std::vector<GlobalPoint> coreTrajectory;
484  for(std::vector<SteppingHelixStateInfo>::const_iterator itr = trajectoryStates.begin();
485  itr != trajectoryStates.end(); itr++) coreTrajectory.push_back(itr->position());
486 
487  if(coreTrajectory.empty()) {
488  LogTrace("TrackAssociator") << "HO trajectory is empty; moving on\n";
489  info.isGoodHO = 0;
490  return;
491  }
492  info.isGoodHO = 1;
493 
494  // find crossed HOs
496  iEvent.getByToken(parameters.HOcollToken, collection);
497  if ( ! collection.isValid() ) throw cms::Exception("FatalError") << "Unable to find HORecHits in event!\n";
498 
499  std::set<DetId> idsInRegion;
500  if (parameters.accountForTrajectoryChangeCalo){
501  // get trajectory change with respect to initial state
502  DetIdAssociator::MapRange mapRange = getMapRange(cachedTrajectory_.trajectoryDelta(CachedTrajectory::IpToHO),
503  parameters.dRHcalPreselection);
504  idsInRegion = hoDetIdAssociator_->getDetIdsCloseToAPoint(coreTrajectory[0], mapRange);
505  } else idsInRegion = hoDetIdAssociator_->getDetIdsCloseToAPoint(coreTrajectory[0], parameters.dRHcalPreselection);
506 
507  LogTrace("TrackAssociator") << "idsInRegion.size(): " << idsInRegion.size();
508  std::set<DetId> idsInACone = hoDetIdAssociator_->getDetIdsInACone(idsInRegion, coreTrajectory, parameters.dRHcal);
509  LogTrace("TrackAssociator") << "idsInACone.size(): " << idsInACone.size();
510  std::vector<DetId> crossedIds =
511  hoDetIdAssociator_->getCrossedDetIds(idsInRegion, coreTrajectory);
512 
513  info.crossedHOIds = crossedIds;
514 
515  // add HO
516  for(std::vector<DetId>::const_iterator itr=crossedIds.begin(); itr!=crossedIds.end();itr++)
517  {
518  HORecHitCollection::const_iterator hit = (*collection).find(*itr);
519  if( hit != (*collection).end() )
520  info.crossedHORecHits.push_back(&*hit);
521  else
522  LogTrace("TrackAssociator") << "Crossed HORecHit is not found for DetId: " << itr->rawId();
523  }
524 
525  for(std::set<DetId>::const_iterator itr=idsInACone.begin(); itr!=idsInACone.end();itr++)
526  {
527  HORecHitCollection::const_iterator hit = (*collection).find(*itr);
528  if( hit != (*collection).end() )
529  info.hoRecHits.push_back(&*hit);
530  else
531  LogTrace("TrackAssociator") << "HORecHit from the cone is not found for DetId: " << itr->rawId();
532  }
533 }
534 
536  const SimTrack& track,
537  const SimVertex& vertex )
538 {
539  GlobalVector vector( track.momentum().x(), track.momentum().y(), track.momentum().z() );
540  GlobalPoint point( vertex.position().x(), vertex.position().y(), vertex.position().z() );
541 
542  int charge = track.type( )> 0 ? -1 : 1; // lepton convention
543  if ( abs(track.type( )) == 211 || // pion
544  abs(track.type( )) == 321 || // kaon
545  abs(track.type( )) == 2212 )
546  charge = track.type( )< 0 ? -1 : 1;
547  return getFreeTrajectoryState(iSetup, vector, point, charge);
548 }
549 
551  const GlobalVector& momentum,
552  const GlobalPoint& vertex,
553  const int charge)
554 {
556  iSetup.get<IdealMagneticFieldRecord>().get(bField);
557 
558  GlobalTrajectoryParameters tPars(vertex, momentum, charge, &*bField);
559 
560  ROOT::Math::SMatrixIdentity id;
561  AlgebraicSymMatrix66 covT(id); covT *= 1e-6; // initialize to sigma=1e-3
562  CartesianTrajectoryError tCov(covT);
563 
564  return FreeTrajectoryState(tPars, tCov);
565 }
566 
567 
569  const reco::Track& track )
570 {
572  iSetup.get<IdealMagneticFieldRecord>().get(bField);
573 
574  GlobalVector vector( track.momentum().x(), track.momentum().y(), track.momentum().z() );
575 
576  GlobalPoint point( track.vertex().x(), track.vertex().y(), track.vertex().z() );
577 
578  GlobalTrajectoryParameters tPars(point, vector, track.charge(), &*bField);
579 
580  // FIX THIS !!!
581  // need to convert from perigee to global or helix (curvilinear) frame
582  // for now just an arbitrary matrix.
583  ROOT::Math::SMatrixIdentity id;
584  AlgebraicSymMatrix66 covT(id); covT *= 1e-6; // initialize to sigma=1e-3
585  CartesianTrajectoryError tCov(covT);
586 
587  return FreeTrajectoryState(tPars, tCov);
588 }
589 
591  const float dR )
592 {
593  DetIdAssociator::MapRange mapRange;
594  mapRange.dThetaPlus = dR;
595  mapRange.dThetaMinus = dR;
596  mapRange.dPhiPlus = dR;
597  mapRange.dPhiMinus = dR;
598  if ( delta.first > 0 )
599  mapRange.dThetaPlus += delta.first;
600  else
601  mapRange.dThetaMinus += fabs(delta.first);
602  if ( delta.second > 0 )
603  mapRange.dPhiPlus += delta.second;
604  else
605  mapRange.dPhiMinus += fabs(delta.second);
606  LogTrace("TrackAssociator") << "Selection range: (dThetaPlus, dThetaMinus, dPhiPlus, dPhiMinus, dRPreselection): " <<
607  mapRange.dThetaPlus << ", " << mapRange.dThetaMinus << ", " <<
608  mapRange.dPhiPlus << ", " << mapRange.dPhiMinus << ", " << dR;
609  return mapRange;
610 }
611 
612 void TrackDetectorAssociator::getTAMuonChamberMatches(std::vector<TAMuonChamberMatch>& matches,
614 {
615  // Strategy:
616  // Propagate through the whole detector, estimate change in eta and phi
617  // along the trajectory, add this to dRMuon and find DetIds around this
618  // direction using the map. Then propagate fast to each surface and apply
619  // final matching criteria.
620 
621  // get the direction first
622  SteppingHelixStateInfo trajectoryPoint = cachedTrajectory_.getStateAtHcal();
623  // If trajectory point at HCAL is not valid, try to use the outer most state of the
624  // trajectory instead.
625  if(! trajectoryPoint.isValid() ) trajectoryPoint = cachedTrajectory_.getOuterState();
626  if(! trajectoryPoint.isValid() ) {
627  LogTrace("TrackAssociator") <<
628  "trajectory position at HCAL is not valid. Assume the track cannot reach muon detectors and skip it";
629  return;
630  }
631 
632  GlobalVector direction = trajectoryPoint.momentum().unit();
633  LogTrace("TrackAssociator") << "muon direction: " << direction << "\n\t and corresponding point: " <<
634  trajectoryPoint.position() <<"\n";
635 
636  DetIdAssociator::MapRange mapRange = getMapRange(cachedTrajectory_.trajectoryDelta(CachedTrajectory::FullTrajectory),
637  parameters.dRMuonPreselection);
638 
639  // and find chamber DetIds
640 
641  std::set<DetId> muonIdsInRegion =
642  muonDetIdAssociator_->getDetIdsCloseToAPoint(trajectoryPoint.position(), mapRange);
643  LogTrace("TrackAssociator") << "Number of chambers to check: " << muonIdsInRegion.size();
644  for(std::set<DetId>::const_iterator detId = muonIdsInRegion.begin(); detId != muonIdsInRegion.end(); detId++)
645  {
646  const GeomDet* geomDet = muonDetIdAssociator_->getGeomDet(*detId);
647  TrajectoryStateOnSurface stateOnSurface = cachedTrajectory_.propagate( &geomDet->surface() );
648  if (! stateOnSurface.isValid()) {
649  LogTrace("TrackAssociator") << "Failed to propagate the track; moving on\n\t"<<
650  "Element is not crosssed: " << DetIdInfo::info(*detId,0) << "\n";
651  continue;
652  }
653  LocalPoint localPoint = geomDet->surface().toLocal(stateOnSurface.freeState()->position());
654  LocalError localError = stateOnSurface.localError().positionError();
655  float distanceX = 0;
656  float distanceY = 0;
657  float sigmaX = 0.0;
658  float sigmaY = 0.0;
659  if(const CSCChamber* cscChamber = dynamic_cast<const CSCChamber*>(geomDet) ) {
660  const CSCChamberSpecs* chamberSpecs = cscChamber->specs();
661  if(! chamberSpecs) {
662  LogTrace("TrackAssociator") << "Failed to get CSCChamberSpecs from CSCChamber; moving on\n";
663  continue;
664  }
665  const CSCLayerGeometry* layerGeometry = chamberSpecs->oddLayerGeometry(1);
666  if(! layerGeometry) {
667  LogTrace("TrackAssociator") << "Failed to get CSCLayerGeometry from CSCChamberSpecs; moving on\n";
668  continue;
669  }
670  const CSCWireTopology* wireTopology = layerGeometry->wireTopology();
671  if(! wireTopology) {
672  LogTrace("TrackAssociator") << "Failed to get CSCWireTopology from CSCLayerGeometry; moving on\n";
673  continue;
674  }
675 
676  float wideWidth = wireTopology->wideWidthOfPlane();
677  float narrowWidth = wireTopology->narrowWidthOfPlane();
678  float length = wireTopology->lengthOfPlane();
679  // If slanted, there is no y offset between local origin and symmetry center of wire plane
680  float yOfFirstWire = fabs(wireTopology->wireAngle())>1.E-06 ? -0.5*length : wireTopology->yOfWire(1);
681  // y offset between local origin and symmetry center of wire plane
682  float yCOWPOffset = yOfFirstWire+0.5*length;
683 
684  // tangent of the incline angle from inside the trapezoid
685  float tangent = (wideWidth-narrowWidth)/(2.*length);
686  // y position wrt bottom of trapezoid
687  float yPrime = localPoint.y()+fabs(yOfFirstWire);
688  // half trapezoid width at y' is 0.5 * narrowWidth + x side of triangle with the above tangent and side y'
689  float halfWidthAtYPrime = 0.5*narrowWidth+yPrime*tangent;
690  distanceX = fabs(localPoint.x()) - halfWidthAtYPrime;
691  distanceY = fabs(localPoint.y()-yCOWPOffset) - 0.5*length;
692  sigmaX = distanceX/sqrt(localError.xx());
693  sigmaY = distanceY/sqrt(localError.yy());
694  } else {
695  distanceX = fabs(localPoint.x()) - geomDet->surface().bounds().width()/2.;
696  distanceY = fabs(localPoint.y()) - geomDet->surface().bounds().length()/2.;
697  sigmaX = distanceX/sqrt(localError.xx());
698  sigmaY = distanceY/sqrt(localError.yy());
699  }
700  if ( (distanceX < parameters.muonMaxDistanceX && distanceY < parameters.muonMaxDistanceY) ||
701  (sigmaX < parameters.muonMaxDistanceSigmaX && sigmaY < parameters.muonMaxDistanceSigmaY) ) {
702  LogTrace("TrackAssociator") << "found a match: " << DetIdInfo::info(*detId,0) << "\n";
704  match.tState = stateOnSurface;
705  match.localDistanceX = distanceX;
706  match.localDistanceY = distanceY;
707  match.id = *detId;
708  matches.push_back(match);
709  } else {
710  LogTrace("TrackAssociator") << "chamber is too far: " <<
711  DetIdInfo::info(*detId,0) << "\n\tdistanceX: " << distanceX << "\t distanceY: " << distanceY <<
712  "\t sigmaX: " << sigmaX << "\t sigmaY: " << sigmaY << "\n";
713  }
714  }
715 
716 }
717 
721 {
722  // Get the segments from the event
724  iEvent.getByToken(parameters.dtSegmentsToken, dtSegments);
725  if (! dtSegments.isValid())
726  throw cms::Exception("FatalError") << "Unable to find DTRecSegment4DCollection in event!\n";
727 
729  iEvent.getByToken(parameters.cscSegmentsToken, cscSegments);
730  if (! cscSegments.isValid())
731  throw cms::Exception("FatalError") << "Unable to find CSCSegmentCollection in event!\n";
732 
734 
735  // check the map of available segments
736  // if there is no segments in a given direction at all,
737  // then there is no point to fly there.
738  //
739  // MISSING
740  // Possible solution: quick search for presence of segments
741  // for the set of DetIds
742 
743  // get a set of matches corresponding to muon chambers
744  std::vector<TAMuonChamberMatch> matchedChambers;
745  LogTrace("TrackAssociator") << "Trying to Get ChamberMatches" << std::endl;
746  getTAMuonChamberMatches(matchedChambers, parameters);
747  LogTrace("TrackAssociator") << "Chambers matched: " << matchedChambers.size() << "\n";
748 
749  // Iterate over all chamber matches and fill segment matching
750  // info if it's available
751  for(std::vector<TAMuonChamberMatch>::iterator matchedChamber = matchedChambers.begin();
752  matchedChamber != matchedChambers.end(); matchedChamber++)
753  {
754  const GeomDet* geomDet = muonDetIdAssociator_->getGeomDet((*matchedChamber).id);
755  // DT chamber
756  if(const DTChamber* chamber = dynamic_cast<const DTChamber*>(geomDet) ) {
757  // Get the range for the corresponding segments
758  DTRecSegment4DCollection::range range = dtSegments->get(chamber->id());
759  // Loop over the segments of this chamber
760  for (DTRecSegment4DCollection::const_iterator segment = range.first; segment!=range.second; segment++) {
761  if (addTAMuonSegmentMatch(*matchedChamber, &(*segment), parameters)) {
762  matchedChamber->segments.back().dtSegmentRef = DTRecSegment4DRef(dtSegments, segment - dtSegments->begin());
763  }
764  }
765  }else{
766  // CSC Chamber
767  if(const CSCChamber* chamber = dynamic_cast<const CSCChamber*>(geomDet) ) {
768  // Get the range for the corresponding segments
769  CSCSegmentCollection::range range = cscSegments->get(chamber->id());
770  // Loop over the segments
771  for (CSCSegmentCollection::const_iterator segment = range.first; segment!=range.second; segment++) {
772  if (addTAMuonSegmentMatch(*matchedChamber, &(*segment), parameters)) {
773  matchedChamber->segments.back().cscSegmentRef = CSCSegmentRef(cscSegments, segment - cscSegments->begin());
774  }
775  }
776  }else{
777  // throw cms::Exception("FatalError") << "Failed to cast GeomDet object to either DTChamber or CSCChamber. Who is this guy anyway?\n";
778  }
779  }
780  info.chambers.push_back(*matchedChamber);
781  }
782 }
783 
784 
786  const RecSegment* segment,
788 {
789  LogTrace("TrackAssociator")
790  << "Segment local position: " << segment->localPosition() << "\n"
791  << std::hex << segment->geographicalId().rawId() << "\n";
792 
793  const GeomDet* chamber = muonDetIdAssociator_->getGeomDet(matchedChamber.id);
794  TrajectoryStateOnSurface trajectoryStateOnSurface = matchedChamber.tState;
795  GlobalPoint segmentGlobalPosition = chamber->toGlobal(segment->localPosition());
796 
797  LogTrace("TrackAssociator")
798  << "Segment global position: " << segmentGlobalPosition << " \t (R_xy,eta,phi): "
799  << segmentGlobalPosition.perp() << "," << segmentGlobalPosition.eta() << "," << segmentGlobalPosition.phi() << "\n";
800 
801  LogTrace("TrackAssociator")
802  << "\teta hit: " << segmentGlobalPosition.eta() << " \tpropagator: " << trajectoryStateOnSurface.freeState()->position().eta() << "\n"
803  << "\tphi hit: " << segmentGlobalPosition.phi() << " \tpropagator: " << trajectoryStateOnSurface.freeState()->position().phi() << std::endl;
804 
805  bool isGood = false;
806  bool isDTWithoutY = false;
807  const DTRecSegment4D* dtseg = dynamic_cast<const DTRecSegment4D*>(segment);
808  if ( dtseg && (! dtseg->hasZed()) )
809  isDTWithoutY = true;
810 
811  double deltaPhi(fabs(segmentGlobalPosition.phi()-trajectoryStateOnSurface.freeState()->position().phi()));
812  if(deltaPhi>M_PI) deltaPhi = fabs(deltaPhi-M_PI*2.);
813 
814  if( isDTWithoutY )
815  {
816  isGood = deltaPhi < parameters.dRMuon;
817  // Be in chamber
818  isGood &= fabs(segmentGlobalPosition.eta()-trajectoryStateOnSurface.freeState()->position().eta()) < .3;
819  } else isGood = sqrt( pow(segmentGlobalPosition.eta()-trajectoryStateOnSurface.freeState()->position().eta(),2) +
820  deltaPhi*deltaPhi) < parameters.dRMuon;
821 
822  if(isGood) {
823  TAMuonSegmentMatch muonSegment;
824  muonSegment.segmentGlobalPosition = getPoint(segmentGlobalPosition);
825  muonSegment.segmentLocalPosition = getPoint( segment->localPosition() );
826  muonSegment.segmentLocalDirection = getVector( segment->localDirection() );
827  muonSegment.segmentLocalErrorXX = segment->localPositionError().xx();
828  muonSegment.segmentLocalErrorYY = segment->localPositionError().yy();
829  muonSegment.segmentLocalErrorXY = segment->localPositionError().xy();
830  muonSegment.segmentLocalErrorDxDz = segment->localDirectionError().xx();
831  muonSegment.segmentLocalErrorDyDz = segment->localDirectionError().yy();
832 
833  // DANGEROUS - compiler cannot guaranty parameters ordering
834  // AlgebraicSymMatrix segmentCovMatrix = segment->parametersError();
835  // muonSegment.segmentLocalErrorXDxDz = segmentCovMatrix[2][0];
836  // muonSegment.segmentLocalErrorYDyDz = segmentCovMatrix[3][1];
837  muonSegment.segmentLocalErrorXDxDz = -999;
838  muonSegment.segmentLocalErrorYDyDz = -999;
839  muonSegment.hasZed = true;
840  muonSegment.hasPhi = true;
841 
842  // timing information
843  muonSegment.t0 = 0;
844  if ( dtseg ) {
845  if ( (dtseg->hasPhi()) && (! isDTWithoutY) ) {
846  int phiHits = dtseg->phiSegment()->specificRecHits().size();
847  // int zHits = dtseg->zSegment()->specificRecHits().size();
848  int hits=0;
849  double t0=0.;
850  // TODO: cuts on hit numbers not optimized in any way yet...
851  if (phiHits>5 && dtseg->phiSegment()->ist0Valid()) {
852  t0+=dtseg->phiSegment()->t0()*phiHits;
853  hits+=phiHits;
854  LogTrace("TrackAssociator") << " Phi t0: " << dtseg->phiSegment()->t0() << " hits: " << phiHits;
855  }
856  // the z segments seem to contain little useful information...
857 // if (zHits>3) {
858 // t0+=s->zSegment()->t0()*zHits;
859 // hits+=zHits;
860 // std::cout << " Z t0: " << s->zSegment()->t0() << " hits: " << zHits << std::endl;
861 // }
862  if (hits) muonSegment.t0 = t0/hits;
863 // std::cout << " --- t0: " << muonSegment.t0 << std::endl;
864  } else {
865  // check and set dimensionality
866  if (isDTWithoutY) muonSegment.hasZed = false;
867  if (! dtseg->hasPhi()) muonSegment.hasPhi = false;
868  }
869  }
870  matchedChamber.segments.push_back(muonSegment);
871  }
872 
873  return isGood;
874 }
875 
876 //********************** NON-CORE CODE ******************************//
877 
881 {
882  // get list of simulated tracks and their vertices
883  using namespace edm;
884  Handle<SimTrackContainer> simTracks;
885  iEvent.getByToken(parameters.simTracksToken, simTracks);
886  if (! simTracks.isValid() ) throw cms::Exception("FatalError") << "No simulated tracks found\n";
887 
888  Handle<SimVertexContainer> simVertices;
889  iEvent.getByToken(parameters.simVerticesToken, simVertices);
890  if (! simVertices.isValid() ) throw cms::Exception("FatalError") << "No simulated vertices found\n";
891 
892  // get sim calo hits
893  Handle<PCaloHitContainer> simEcalHitsEB;
894  iEvent.getByToken(parameters.simEcalHitsEBToken, simEcalHitsEB);
895  if (! simEcalHitsEB.isValid() ) throw cms::Exception("FatalError") << "No simulated ECAL EB hits found\n";
896 
897  Handle<PCaloHitContainer> simEcalHitsEE;
898  iEvent.getByToken(parameters.simEcalHitsEEToken, simEcalHitsEE);
899  if (! simEcalHitsEE.isValid() ) throw cms::Exception("FatalError") << "No simulated ECAL EE hits found\n";
900 
901  Handle<PCaloHitContainer> simHcalHits;
902  iEvent.getByToken(parameters.simHcalHitsToken, simHcalHits);
903  if (! simHcalHits.isValid() ) throw cms::Exception("FatalError") << "No simulated HCAL hits found\n";
904 
905  // find truth partner
906  SimTrackContainer::const_iterator simTrack = simTracks->begin();
907  for( ; simTrack != simTracks->end(); ++simTrack){
908  math::XYZVector simP3( simTrack->momentum().x(), simTrack->momentum().y(), simTrack->momentum().z() );
909  math::XYZVector recoP3( info.stateAtIP.momentum().x(), info.stateAtIP.momentum().y(), info.stateAtIP.momentum().z() );
910  if ( ROOT::Math::VectorUtil::DeltaR(recoP3, simP3) < 0.1 ) break;
911  }
912  if ( simTrack != simTracks->end() ) {
913  info.simTrack = &(*simTrack);
914  double ecalTrueEnergy(0);
915  double hcalTrueEnergy(0);
916 
917  // loop over calo hits
918  for( PCaloHitContainer::const_iterator hit = simEcalHitsEB->begin(); hit != simEcalHitsEB->end(); ++hit )
919  if ( hit->geantTrackId() == info.simTrack->genpartIndex() ) ecalTrueEnergy += hit->energy();
920 
921  for( PCaloHitContainer::const_iterator hit = simEcalHitsEE->begin(); hit != simEcalHitsEE->end(); ++hit )
922  if ( hit->geantTrackId() == info.simTrack->genpartIndex() ) ecalTrueEnergy += hit->energy();
923 
924  for( PCaloHitContainer::const_iterator hit = simHcalHits->begin(); hit != simHcalHits->end(); ++hit )
925  if ( hit->geantTrackId() == info.simTrack->genpartIndex() ) hcalTrueEnergy += hit->energy();
926 
927  info.ecalTrueEnergy = ecalTrueEnergy;
928  info.hcalTrueEnergy = hcalTrueEnergy;
929  info.hcalTrueEnergyCorrected = hcalTrueEnergy;
930  if ( fabs(info.trkGlobPosAtHcal.eta()) < 1.3 )
931  info.hcalTrueEnergyCorrected = hcalTrueEnergy*113.2;
932  else
933  if ( fabs(info.trkGlobPosAtHcal.eta()) < 3.0 )
934  info.hcalTrueEnergyCorrected = hcalTrueEnergy*167.2;
935  }
936 }
937 
939  const edm::EventSetup& iSetup,
940  const reco::Track& track,
942  Direction direction /*= Any*/ )
943 {
944  double currentStepSize = cachedTrajectory_.getPropagationStep();
945 
947  iSetup.get<IdealMagneticFieldRecord>().get(bField);
948 
949  if(track.extra().isNull()) {
950  if ( direction != InsideOut )
951  throw cms::Exception("FatalError") <<
952  "No TrackExtra information is available and association is done with something else than InsideOut track.\n" <<
953  "Either change the parameter or provide needed data!\n";
954  LogTrace("TrackAssociator") << "Track Extras not found\n";
956  return associate(iEvent, iSetup, parameters, &initialState); // 5th argument is null pointer
957  }
958 
959  LogTrace("TrackAssociator") << "Track Extras found\n";
962  FreeTrajectoryState referenceState = trajectoryStateTransform::initialFreeState(track,&*bField);
963 
964  LogTrace("TrackAssociator") << "inner track state (rho, z, phi):" <<
965  track.innerPosition().Rho() << ", " << track.innerPosition().z() <<
966  ", " << track.innerPosition().phi() << "\n";
967  LogTrace("TrackAssociator") << "innerFreeState (rho, z, phi):" <<
968  innerState.position().perp() << ", " << innerState.position().z() <<
969  ", " << innerState.position().phi() << "\n";
970 
971  LogTrace("TrackAssociator") << "outer track state (rho, z, phi):" <<
972  track.outerPosition().Rho() << ", " << track.outerPosition().z() <<
973  ", " << track.outerPosition().phi() << "\n";
974  LogTrace("TrackAssociator") << "outerFreeState (rho, z, phi):" <<
975  outerState.position().perp() << ", " << outerState.position().z() <<
976  ", " << outerState.position().phi() << "\n";
977 
978  // InsideOut first
979  if ( crossedIP( track ) ) {
980  switch ( direction ) {
981  case InsideOut:
982  case Any:
983  return associate(iEvent, iSetup, parameters, &referenceState, &outerState);
984  break;
985  case OutsideIn:
986  {
987  cachedTrajectory_.setPropagationStep( -fabs(currentStepSize) );
988  TrackDetMatchInfo result = associate(iEvent, iSetup, parameters, &innerState, &referenceState);
989  cachedTrajectory_.setPropagationStep( currentStepSize );
990  return result;
991  break;
992  }
993  }
994  } else {
995  switch ( direction ) {
996  case InsideOut:
997  return associate(iEvent, iSetup, parameters, &innerState, &outerState);
998  break;
999  case OutsideIn:
1000  {
1001  cachedTrajectory_.setPropagationStep( -fabs(currentStepSize) );
1002  TrackDetMatchInfo result = associate(iEvent, iSetup, parameters, &outerState, &innerState);
1003  cachedTrajectory_.setPropagationStep( currentStepSize );
1004  return result;
1005  break;
1006  }
1007  case Any:
1008  {
1009  // check if we deal with clear outside-in case
1010  if ( track.innerPosition().Dot( track.innerMomentum() ) < 0 &&
1011  track.outerPosition().Dot( track.outerMomentum() ) < 0 )
1012  {
1013  cachedTrajectory_.setPropagationStep( -fabs(currentStepSize) );
1015  if ( track.innerPosition().R() < track.outerPosition().R() )
1016  result = associate(iEvent, iSetup, parameters, &innerState, &outerState);
1017  else
1018  result = associate(iEvent, iSetup, parameters, &outerState, &innerState);
1019  cachedTrajectory_.setPropagationStep( currentStepSize );
1020  return result;
1021  }
1022  }
1023  }
1024  }
1025 
1026  // all other cases
1027  return associate(iEvent, iSetup, parameters, &innerState, &outerState);
1028 }
1029 
1031  const edm::EventSetup& iSetup,
1032  const SimTrack& track,
1033  const SimVertex& vertex,
1035 {
1036  return associate(iEvent, iSetup, getFreeTrajectoryState(iSetup, track, vertex), parameters);
1037 }
1038 
1040  const edm::EventSetup& iSetup,
1041  const GlobalVector& momentum,
1042  const GlobalPoint& vertex,
1043  const int charge,
1045 {
1046  return associate(iEvent, iSetup, getFreeTrajectoryState(iSetup, momentum, vertex, charge), parameters);
1047 }
1048 
1050 {
1051  bool crossed = true;
1052  crossed &= (track.innerPosition().rho() > 3 ); // something close to active volume
1053  crossed &= (track.outerPosition().rho() > 3 ); // something close to active volume
1054  crossed &= ( ( track.innerPosition().x()*track.innerMomentum().x() +
1055  track.innerPosition().y()*track.innerMomentum().y() < 0 ) !=
1056  ( track.outerPosition().x()*track.outerMomentum().x() +
1057  track.outerPosition().y()*track.outerMomentum().y() < 0 ) );
1058  return crossed;
1059 }
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
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:464
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:67
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:662
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:674
float wireAngle() const
int iEvent
Definition: GenABIO.cc:230
T mag() const
Definition: PV3DBase.h:67
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)
bool isnan(float x)
Definition: math.h:13
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
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:75
edm::EDGetTokenT< HBHERecHitCollection > HBHEcollToken
void fillPreshower(const edm::Event &iEvent, TrackDetMatchInfo &info, const AssociatorParameters &)
bool isNull() const
Checks for null.
Definition: Ref.h:249
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
static std::string info(const DetId &, const TrackerTopology *tTopo)
Definition: DetIdInfo.cc:29
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
virtual LocalError localPositionError() const =0
const T & get() const
Definition: EventSetup.h:56
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:554
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
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