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