CMS 3D CMS Logo

MuonIdTruthInfo.cc
Go to the documentation of this file.
10 
12  iC.mayConsume<edm::SimTrackContainer>(edm::InputTag("g4SimHits", ""));
13  iC.mayConsume<edm::PSimHitContainer>(edm::InputTag("g4SimHits", "MuonDTHits"));
14  iC.mayConsume<edm::PSimHitContainer>(edm::InputTag("g4SimHits", "MuonCSCHits"));
15 }
16 
18  // get a list of simulated track and find a track with the best match to
19  // the muon.track(). Use its id and chamber id to localize hits
20  // If a hit has non-zero local z coordinate, it's position wrt
21  // to the center of a chamber is extrapolated by a straight line
22 
24  iEvent.getByLabel<edm::SimTrackContainer>("g4SimHits", "", simTracks);
25  if (!simTracks.isValid()) {
26  LogTrace("MuonIdentification") << "No tracks found";
27  return;
28  }
29 
30  // get the tracking Geometry
32  iSetup.get<GlobalTrackingGeometryRecord>().get(geometry);
33  if (!geometry.isValid())
34  throw cms::Exception("FatalError") << "Unable to find GlobalTrackingGeometryRecord in event!\n";
35 
36  float bestMatchChi2 = 9999; //minimization creteria
37  unsigned int bestMatch = 0;
38  const unsigned int offset = 0; // kludge to fix a problem in trackId matching between tracks and hits.
39 
40  for (edm::SimTrackContainer::const_iterator simTrk = simTracks->begin(); simTrk != simTracks->end(); simTrk++) {
41  float chi2 = matchChi2(*aMuon.track().get(), *simTrk);
42  if (chi2 > bestMatchChi2)
43  continue;
44  bestMatchChi2 = chi2;
45  bestMatch = simTrk->trackId();
46  }
47 
48  bestMatch -= offset;
49 
50  std::vector<reco::MuonChamberMatch>& matches = aMuon.matches();
51  int numberOfTruthMatchedChambers = 0;
52 
53  // loop over chambers
54  for (std::vector<reco::MuonChamberMatch>::iterator chamberMatch = matches.begin(); chamberMatch != matches.end();
55  chamberMatch++) {
56  if (chamberMatch->id.det() != DetId::Muon) {
57  edm::LogWarning("MuonIdentification") << "Detector id of a muon chamber corresponds to not a muon detector";
58  continue;
59  }
60  reco::MuonSegmentMatch bestSegmentMatch;
61  double distance = 99999;
62 
63  if (chamberMatch->id.subdetId() == MuonSubdetId::DT) {
64  DTChamberId detId(chamberMatch->id.rawId());
65 
67  iEvent.getByLabel("g4SimHits", "MuonDTHits", simHits);
68  if (simHits.isValid()) {
69  for (edm::PSimHitContainer::const_iterator hit = simHits->begin(); hit != simHits->end(); hit++)
70  if (hit->trackId() == bestMatch)
71  checkSimHitForBestMatch(bestSegmentMatch, distance, *hit, detId, geometry);
72  } else
73  LogTrace("MuonIdentification") << "No DT simulated hits are found";
74  }
75 
76  if (chamberMatch->id.subdetId() == MuonSubdetId::CSC) {
77  CSCDetId detId(chamberMatch->id.rawId());
78 
80  iEvent.getByLabel("g4SimHits", "MuonCSCHits", simHits);
81  if (simHits.isValid()) {
82  for (edm::PSimHitContainer::const_iterator hit = simHits->begin(); hit != simHits->end(); hit++)
83  if (hit->trackId() == bestMatch)
84  checkSimHitForBestMatch(bestSegmentMatch, distance, *hit, detId, geometry);
85  } else
86  LogTrace("MuonIdentification") << "No CSC simulated hits are found";
87  }
88  if (distance < 9999) {
89  chamberMatch->truthMatches.push_back(bestSegmentMatch);
90  numberOfTruthMatchedChambers++;
91  LogTrace("MuonIdentification") << "Best truth matched hit:"
92  << "\tDetId: " << chamberMatch->id.rawId() << "\n"
93  << "\tprojection: ( " << bestSegmentMatch.x << ", " << bestSegmentMatch.y
94  << " )\n";
95  }
96  }
97  LogTrace("MuonIdentification") << "Truth matching summary:\n\tnumber of chambers: " << matches.size()
98  << "\n\tnumber of truth matched chambers: " << numberOfTruthMatchedChambers << "\n";
99 }
100 
102  double& distance,
103  const PSimHit& hit,
104  const DetId& chamberId,
106  printf("DONT FORGET TO CALL REGISTERCONSUMES()\n");
107 
108  // find the hit position projection at the reference surface of the chamber:
109  // first get entry and exit point of the hit in the global coordinates, then
110  // get local coordinates of these points wrt the chamber and then find the
111  // projected X-Y coordinates
112 
113  const GeomDet* chamberGeometry = geometry->idToDet(chamberId);
114  const GeomDet* simUnitGeometry = geometry->idToDet(DetId(hit.detUnitId()));
115 
116  if (chamberGeometry && simUnitGeometry) {
117  LocalPoint entryPoint = chamberGeometry->toLocal(simUnitGeometry->toGlobal(hit.entryPoint()));
118  LocalPoint exitPoint = chamberGeometry->toLocal(simUnitGeometry->toGlobal(hit.exitPoint()));
119  LocalVector direction = exitPoint - entryPoint;
120  if (fabs(direction.z()) > 0.001) {
121  LocalPoint projection = entryPoint - direction * (entryPoint.z() / direction.z());
122  if (fabs(projection.z()) > 0.001)
123  edm::LogWarning("MuonIdentification") << "z coordinate of the hit projection must be zero and it's not!\n";
124 
125  double new_distance = 99999;
126  if (entryPoint.z() * exitPoint.z() < -1) // changed sign, so the reference point is inside
127  new_distance = 0;
128  else {
129  if (fabs(entryPoint.z()) < fabs(exitPoint.z()))
130  new_distance = fabs(entryPoint.z());
131  else
132  new_distance = fabs(exitPoint.z());
133  }
134 
135  if (new_distance < distance) {
136  // find a SimHit closer to the reference surface, update segmentMatch
137  segmentMatch.x = projection.x();
138  segmentMatch.y = projection.y();
139  segmentMatch.xErr = 0;
140  segmentMatch.yErr = 0;
141  segmentMatch.dXdZ = direction.x() / direction.z();
142  segmentMatch.dYdZ = direction.y() / direction.z();
143  segmentMatch.dXdZErr = 0;
144  segmentMatch.dYdZErr = 0;
145  distance = new_distance;
146  LogTrace("MuonIdentificationVerbose")
147  << "Better truth matched segment found:\n"
148  << "\tDetId: " << chamberId.rawId() << "\n"
149  << "\tentry point: ( " << entryPoint.x() << ", " << entryPoint.y() << ", " << entryPoint.z() << " )\n"
150  << "\texit point: ( " << exitPoint.x() << ", " << exitPoint.y() << ", " << exitPoint.z() << " )\n"
151  << "\tprojection: ( " << projection.x() << ", " << projection.y() << ", " << projection.z() << " )\n";
152  }
153  }
154  } else {
155  if (!chamberGeometry)
156  edm::LogWarning("MuonIdentification") << "Cannot get chamber geomtry for DetId: " << chamberId.rawId();
157  if (!simUnitGeometry)
158  edm::LogWarning("MuonIdentification") << "Cannot get detector unit geomtry for DetId: " << hit.detUnitId();
159  }
160 }
161 
162 double MuonIdTruthInfo::matchChi2(const reco::Track& recoTrk, const SimTrack& simTrk) {
163  double deltaPhi = fabs(recoTrk.phi() - simTrk.momentum().phi());
164  if (deltaPhi > 1.8 * 3.1416)
165  deltaPhi = 2 * 3.1416 - deltaPhi; // take care of phi discontinuity
166  return pow((recoTrk.p() - simTrk.momentum().rho()) / simTrk.momentum().rho(), 2) + pow(deltaPhi / (2 * 3.1416), 2) +
167  pow((recoTrk.theta() - simTrk.momentum().theta()) / 3.1416, 2);
168 }
double p() const
momentum vector magnitude
Definition: TrackBase.h:599
static void truthMatchMuon(const edm::Event &iEvent, const edm::EventSetup &iSetup, reco::Muon &aMuon)
TrackRef track() const override
reference to a Track
Definition: Muon.h:46
EDGetTokenT< ProductType > mayConsume(edm::InputTag const &tag)
double theta() const
polar angle
Definition: TrackBase.h:581
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
T y() const
Definition: PV3DBase.h:60
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:58
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:614
static void checkSimHitForBestMatch(reco::MuonSegmentMatch &segmentMatch, double &distance, const PSimHit &hit, const DetId &chamberId, const edm::ESHandle< GlobalTrackingGeometry > &geometry)
int iEvent
Definition: GenABIO.cc:224
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:46
T z() const
Definition: PV3DBase.h:61
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:232
bool isValid() const
Definition: HandleBase.h:70
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:488
#define LogTrace(id)
Definition: DetId.h:17
def bestMatch(object, matchCollection)
Definition: deltar.py:138
void registerConsumes(edm::ConsumesCollector &iC)
std::vector< MuonChamberMatch > & matches()
get muon matching information
Definition: Muon.h:145
ESHandle< TrackerGeometry > geometry
const math::XYZTLorentzVectorD & momentum() const
Definition: CoreSimTrack.h:19
T get() const
Definition: EventSetup.h:73
const GeomDet * idToDet(DetId) const override
std::vector< PSimHit > PSimHitContainer
static constexpr int DT
Definition: MuonSubdetId.h:11
bool isValid() const
Definition: ESHandle.h:44
static constexpr int CSC
Definition: MuonSubdetId.h:12
T x() const
Definition: PV3DBase.h:59
std::vector< SimTrack > SimTrackContainer
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:43
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:30
unsigned int detUnitId() const
Definition: PSimHit.h:97
static double matchChi2(const reco::Track &recoTrk, const SimTrack &simTrk)