CMS 3D CMS Logo

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