CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_9_patch3/src/RecoMuon/MuonIdentification/src/MuonIdTruthInfo.cc

Go to the documentation of this file.
00001 #include "RecoMuon/MuonIdentification/interface/MuonIdTruthInfo.h"
00002 #include "SimDataFormats/Track/interface/SimTrackContainer.h"
00003 #include "SimDataFormats/TrackingHit/interface/PSimHitContainer.h"
00004 #include "DataFormats/MuonDetId/interface/CSCDetId.h"
00005 #include "DataFormats/MuonDetId/interface/DTChamberId.h"
00006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00007 #include "Geometry/Records/interface/GlobalTrackingGeometryRecord.h"
00008 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
00009 
00010 void MuonIdTruthInfo::truthMatchMuon(const edm::Event& iEvent, 
00011                                      const edm::EventSetup& iSetup,
00012                                      reco::Muon& aMuon)
00013 {
00014    // get a list of simulated track and find a track with the best match to
00015    // the muon.track(). Use its id and chamber id to localize hits
00016    // If a hit has non-zero local z coordinate, it's position wrt
00017    // to the center of a chamber is extrapolated by a straight line
00018    
00019    edm::Handle<edm::SimTrackContainer> simTracks;
00020    iEvent.getByType<edm::SimTrackContainer>(simTracks);
00021    if (! simTracks.isValid() ) {
00022       LogTrace("MuonIdentification") <<"No tracks found";
00023       return;
00024    }
00025    
00026    // get the tracking Geometry
00027    edm::ESHandle<GlobalTrackingGeometry> geometry;
00028    iSetup.get<GlobalTrackingGeometryRecord>().get(geometry);
00029    if ( ! geometry.isValid() )
00030      throw cms::Exception("FatalError") << "Unable to find GlobalTrackingGeometryRecord in event!\n";
00031    
00032    float bestMatchChi2 = 9999; //minimization creteria
00033    unsigned int bestMatch = 0;
00034    const unsigned int offset = 0; // kludge to fix a problem in trackId matching between tracks and hits.
00035    
00036    for( edm::SimTrackContainer::const_iterator simTrk = simTracks->begin();
00037         simTrk != simTracks->end(); simTrk++ )
00038      {
00039         float chi2 = matchChi2(*aMuon.track().get(),*simTrk);
00040         if (chi2>bestMatchChi2) continue;
00041         bestMatchChi2 = chi2;
00042         bestMatch = simTrk->trackId();
00043      }
00044    
00045    bestMatch -= offset;
00046    
00047    std::vector<reco::MuonChamberMatch>& matches = aMuon.matches();
00048    int numberOfTruthMatchedChambers = 0;
00049 
00050    // loop over chambers
00051    for(std::vector<reco::MuonChamberMatch>::iterator chamberMatch = matches.begin();
00052        chamberMatch != matches.end(); chamberMatch ++)
00053      {
00054         if (chamberMatch->id.det() != DetId::Muon ) {
00055            edm::LogWarning("MuonIdentification") << "Detector id of a muon chamber corresponds to not a muon detector";
00056            continue;
00057         }
00058         reco::MuonSegmentMatch bestSegmentMatch;
00059         double distance = 99999;
00060         
00061         if ( chamberMatch->id.subdetId() == MuonSubdetId::DT) {
00062            DTChamberId detId(chamberMatch->id.rawId());
00063 
00064            edm::Handle<edm::PSimHitContainer> simHits;
00065            iEvent.getByLabel("g4SimHits", "MuonDTHits", simHits);
00066            if ( simHits.isValid() ) {
00067               for( edm::PSimHitContainer::const_iterator hit = simHits->begin(); hit != simHits->end(); hit++)
00068                 if ( hit->trackId() == bestMatch ) checkSimHitForBestMatch(bestSegmentMatch, distance, *hit, detId, geometry );
00069            }else LogTrace("MuonIdentification") <<"No DT simulated hits are found";
00070         }
00071 
00072         if ( chamberMatch->id.subdetId() == MuonSubdetId::CSC) {
00073            CSCDetId detId(chamberMatch->id.rawId());
00074 
00075            edm::Handle<edm::PSimHitContainer> simHits;
00076            iEvent.getByLabel("g4SimHits", "MuonCSCHits", simHits);
00077            if ( simHits.isValid() ) {
00078               for( edm::PSimHitContainer::const_iterator hit = simHits->begin(); hit != simHits->end(); hit++)
00079                 if ( hit->trackId() == bestMatch ) checkSimHitForBestMatch(bestSegmentMatch, distance, *hit, detId, geometry );
00080            }else LogTrace("MuonIdentification") <<"No CSC simulated hits are found";
00081         }
00082         if (distance < 9999) {
00083            chamberMatch->truthMatches.push_back( bestSegmentMatch );
00084            numberOfTruthMatchedChambers++;
00085            LogTrace("MuonIdentification") << "Best truth matched hit:" <<
00086               "\tDetId: " << chamberMatch->id.rawId() << "\n" <<
00087               "\tprojection: ( " << bestSegmentMatch.x << ", " << bestSegmentMatch.y << " )\n";
00088         }
00089      }
00090    LogTrace("MuonIdentification") << "Truth matching summary:\n\tnumber of chambers: " << matches.size() <<
00091      "\n\tnumber of truth matched chambers: " << numberOfTruthMatchedChambers << "\n";
00092 }
00093 
00094 void MuonIdTruthInfo::checkSimHitForBestMatch(reco::MuonSegmentMatch& segmentMatch,
00095                                               double& distance,
00096                                               const PSimHit& hit, 
00097                                               const DetId& chamberId,
00098                                               const edm::ESHandle<GlobalTrackingGeometry>& geometry)
00099 {
00100    // find the hit position projection at the reference surface of the chamber:
00101    // first get entry and exit point of the hit in the global coordinates, then
00102    // get local coordinates of these points wrt the chamber and then find the
00103    // projected X-Y coordinates
00104    
00105    const GeomDet* chamberGeometry = geometry->idToDet( chamberId );
00106    const GeomDet* simUnitGeometry = geometry->idToDet( DetId(hit.detUnitId()) );
00107    
00108    if (chamberGeometry && simUnitGeometry ) {
00109       LocalPoint entryPoint = chamberGeometry->toLocal( simUnitGeometry->toGlobal( hit.entryPoint() ) );
00110       LocalPoint exitPoint =  chamberGeometry->toLocal( simUnitGeometry->toGlobal( hit.exitPoint()  ) );
00111       LocalVector direction = exitPoint - entryPoint;
00112       if ( fabs(direction.z()) > 0.001) {
00113          LocalPoint projection = entryPoint - direction*(entryPoint.z()/direction.z());
00114          if ( fabs(projection.z()) > 0.001 ) 
00115            edm::LogWarning("MuonIdentification") << "z coordinate of the hit projection must be zero and it's not!\n";
00116          
00117          double new_distance = 99999;
00118          if( entryPoint.z()*exitPoint.z() < -1 ) // changed sign, so the reference point is inside
00119            new_distance = 0;
00120          else {
00121             if ( fabs(entryPoint.z()) < fabs(exitPoint.z()) )
00122               new_distance = fabs(entryPoint.z());
00123             else
00124               new_distance = fabs(exitPoint.z());
00125          }
00126          
00127          if (new_distance < distance) { 
00128             // find a SimHit closer to the reference surface, update segmentMatch
00129             segmentMatch.x = projection.x();
00130             segmentMatch.y = projection.y();
00131             segmentMatch.xErr = 0;
00132             segmentMatch.yErr = 0;
00133             segmentMatch.dXdZ = direction.x()/direction.z();
00134             segmentMatch.dYdZ = direction.y()/direction.z();
00135             segmentMatch.dXdZErr = 0;
00136             segmentMatch.dYdZErr = 0;
00137             distance = new_distance;
00138             LogTrace("MuonIdentificationVerbose") << "Better truth matched segment found:\n" <<
00139               "\tDetId: " << chamberId.rawId() << "\n" <<
00140               "\tentry point: ( " << entryPoint.x() << ", " << entryPoint.y() << ", " << entryPoint.z() << " )\n" <<
00141               "\texit point: ( " << exitPoint.x() << ", " << exitPoint.y() << ", " << exitPoint.z() << " )\n" <<
00142               "\tprojection: ( " << projection.x() << ", " << projection.y() << ", " << projection.z() << " )\n";
00143          }
00144       }
00145    }else{
00146       if ( ! chamberGeometry ) edm::LogWarning("MuonIdentification") << "Cannot get chamber geomtry for DetId: " << chamberId.rawId();
00147       if ( ! simUnitGeometry ) edm::LogWarning("MuonIdentification") << "Cannot get detector unit geomtry for DetId: " << hit.detUnitId();
00148    }
00149 }
00150 
00151 double MuonIdTruthInfo::matchChi2( const reco::Track& recoTrk, const SimTrack& simTrk)
00152 {
00153    double deltaPhi = fabs(recoTrk.phi()-simTrk.momentum().phi());
00154    if (deltaPhi>1.8*3.1416) deltaPhi = 2*3.1416-deltaPhi; // take care of phi discontinuity
00155    return pow((recoTrk.p()-simTrk.momentum().rho())/simTrk.momentum().rho(),2)+
00156      pow(deltaPhi/(2*3.1416),2)+pow((recoTrk.theta()-simTrk.momentum().theta())/3.1416,2);
00157 }