00001
00009 #include "Validation/RecoMuon/src/MuonTrackAnalyzer.h"
00010
00011
00012 #include "FWCore/Framework/interface/Frameworkfwd.h"
00013 #include "FWCore/Framework/interface/Event.h"
00014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00015 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00016
00017 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
00018
00019 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
00020 #include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
00021 #include "DataFormats/TrackReco/interface/Track.h"
00022 #include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
00023 #include "TrackingTools/DetLayers/interface/DetLayer.h"
00024 #include "DataFormats/Math/interface/deltaR.h"
00025
00026 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
00027
00028 #include "CommonTools/Statistics/interface/ChiSquaredProbability.h"
00029
00030 #include "RecoMuon/TrackingTools/interface/MuonPatternRecoDumper.h"
00031 #include "RecoMuon/TrackingTools/interface/MuonServiceProxy.h"
00032 #include "RecoMuon/TrackingTools/interface/MuonUpdatorAtVertex.h"
00033
00034 #include "Validation/RecoMuon/src/Histograms.h"
00035 #include "Validation/RecoMuon/src/HTrack.h"
00036
00037 #include "SimDataFormats/TrackingHit/interface/PSimHit.h"
00038 #include "SimDataFormats/TrackingHit/interface/PSimHitContainer.h"
00039 #include "SimDataFormats/Track/interface/SimTrackContainer.h"
00040
00041 #include "TFile.h"
00042 #include "TH1F.h"
00043 #include "TH2F.h"
00044
00045 using namespace std;
00046 using namespace edm;
00047
00049 MuonTrackAnalyzer::MuonTrackAnalyzer(const ParameterSet& pset){
00050
00051
00052 ParameterSet serviceParameters = pset.getParameter<ParameterSet>("ServiceParameters");
00053
00054 theService = new MuonServiceProxy(serviceParameters);
00055
00056 theTracksLabel = pset.getParameter<InputTag>("Tracks");
00057 doTracksAnalysis = pset.getUntrackedParameter<bool>("DoTracksAnalysis",true);
00058
00059 doSeedsAnalysis = pset.getUntrackedParameter<bool>("DoSeedsAnalysis",false);
00060 if(doSeedsAnalysis){
00061 theSeedsLabel = pset.getParameter<InputTag>("MuonSeed");
00062 ParameterSet updatorPar = pset.getParameter<ParameterSet>("MuonUpdatorAtVertexParameters");
00063 theSeedPropagatorName = updatorPar.getParameter<string>("Propagator");
00064
00065 theUpdator = new MuonUpdatorAtVertex(updatorPar,theService);
00066 }
00067
00068 theCSCSimHitLabel = pset.getParameter<InputTag>("CSCSimHit");
00069 theDTSimHitLabel = pset.getParameter<InputTag>("DTSimHit");
00070 theRPCSimHitLabel = pset.getParameter<InputTag>("RPCSimHit");
00071
00072 theEtaRange = (EtaRange) pset.getParameter<int>("EtaRange");
00073
00074
00075 numberOfSimTracks=0;
00076
00077 numberOfRecTracks=0;
00078
00079 dbe_ = edm::Service<DQMStore>().operator->();
00080 out = pset.getUntrackedParameter<string>("rootFileName");
00081 dirName_ = pset.getUntrackedParameter<std::string>("dirName");
00082
00083 }
00084
00086 MuonTrackAnalyzer::~MuonTrackAnalyzer(){
00087 if (theService) delete theService;
00088 }
00089
00090 void MuonTrackAnalyzer::beginJob(){
00091 dbe_->showDirStructure();
00092
00093 dbe_->cd();
00094 InputTag algo = theTracksLabel;
00095 string dirName=dirName_;
00096 if (algo.process()!="")
00097 dirName+=algo.process()+"_";
00098 if(algo.label()!="")
00099 dirName+=algo.label()+"_";
00100 if(algo.instance()!="")
00101 dirName+=algo.instance()+"";
00102 if (dirName.find("Tracks")<dirName.length()){
00103 dirName.replace(dirName.find("Tracks"),6,"");
00104 }
00105 std::replace(dirName.begin(), dirName.end(), ':', '_');
00106 dbe_->setCurrentFolder(dirName.c_str());
00107
00108
00109 std::string simName = dirName;
00110 simName+="/SimTracks";
00111 hSimTracks = new HTrackVariables(simName.c_str(),"SimTracks");
00112
00113 dbe_->cd();
00114 dbe_->setCurrentFolder(dirName.c_str());
00115
00116
00117
00118
00119 if(doSeedsAnalysis){
00120 dbe_->cd();
00121 dbe_->setCurrentFolder(dirName.c_str());
00122 hRecoSeedInner = new HTrack(dirName.c_str(),"RecoSeed","Inner");
00123 hRecoSeedPCA = new HTrack(dirName.c_str(),"RecoSeed","PCA");
00124 }
00125
00126 if(doTracksAnalysis){
00127 dbe_->cd();
00128 dbe_->setCurrentFolder(dirName.c_str());
00129 hRecoTracksPCA = new HTrack(dirName.c_str(),"RecoTracks","PCA");
00130 hRecoTracksInner = new HTrack(dirName.c_str(),"RecoTracks","Inner");
00131 hRecoTracksOuter = new HTrack(dirName.c_str(),"RecoTracks","Outer");
00132
00133 dbe_->cd();
00134 dbe_->setCurrentFolder(dirName.c_str());
00135
00136
00137
00138
00139 hChi2 = dbe_->book1D("chi2","#chi^2",200,0,200);
00140 hChi2VsEta = dbe_->book2D("chi2VsEta","#chi^2 VS #eta",120,-3.,3.,200,0,200);
00141
00142 hChi2Norm = dbe_->book1D("chi2Norm","Normalized #chi^2",400,0,100);
00143 hChi2NormVsEta = dbe_->book2D("chi2NormVsEta","Normalized #chi^2 VS #eta",120,-3.,3.,400,0,100);
00144
00145 hHitsPerTrack = dbe_->book1D("HitsPerTrack","Number of hits per track",55,0,55);
00146 hHitsPerTrackVsEta = dbe_->book2D("HitsPerTrackVsEta","Number of hits per track VS #eta",
00147 120,-3.,3.,55,0,55);
00148
00149 hDof = dbe_->book1D("dof","Number of Degree of Freedom",55,0,55);
00150 hDofVsEta = dbe_->book2D("dofVsEta","Number of Degree of Freedom VS #eta",120,-3.,3.,55,0,55);
00151
00152 hChi2Prob = dbe_->book1D("chi2Prob","#chi^2 probability",200,0,1);
00153 hChi2ProbVsEta = dbe_->book2D("chi2ProbVsEta","#chi^2 probability VS #eta",120,-3.,3.,200,0,1);
00154
00155 hNumberOfTracks = dbe_->book1D("NumberOfTracks","Number of reconstructed tracks per event",200,0,200);
00156 hNumberOfTracksVsEta = dbe_->book2D("NumberOfTracksVsEta",
00157 "Number of reconstructed tracks per event VS #eta",
00158 120,-3.,3.,10,0,10);
00159
00160 hChargeVsEta = dbe_->book2D("ChargeVsEta","Charge vs #eta gen",120,-3.,3.,4,-2.,2.);
00161 hChargeVsPt = dbe_->book2D("ChargeVsPt","Charge vs P_{T} gen",250,0,200,4,-2.,2.);
00162 hPtRecVsPtGen = dbe_->book2D("PtRecVsPtGen","P_{T} rec vs P_{T} gen",250,0,200,250,0,200);
00163
00164 hDeltaPtVsEta = dbe_->book2D("DeltaPtVsEta","#Delta P_{t} vs #eta gen",120,-3.,3.,500,-250.,250.);
00165 hDeltaPt_In_Out_VsEta = dbe_->book2D("DeltaPt_In_Out_VsEta_","P^{in}_{t} - P^{out}_{t} vs #eta gen",120,-3.,3.,500,-250.,250.);
00166 }
00167
00168
00169 }
00170
00171 void MuonTrackAnalyzer::endJob(){
00172 LogInfo("MuonTrackAnalyzer")<< "Number of Sim tracks: " << numberOfSimTracks;
00173
00174 LogInfo("MuonTrackAnalyzer") << "Number of Reco tracks: " << numberOfRecTracks;
00175
00176
00177 if(doTracksAnalysis){
00178 double eff = hRecoTracksPCA->computeEfficiency(hSimTracks);
00179 LogInfo("MuonTrackAnalyzer") <<" *Track Efficiency* = "<< eff <<"%";
00180 }
00181
00182 if(doSeedsAnalysis){
00183 double eff = hRecoSeedInner->computeEfficiency(hSimTracks);
00184 LogInfo("MuonTrackAnalyzer")<<" *Seed Efficiency* = "<< eff <<"%";
00185 }
00186
00187 if ( out.size() != 0 && dbe_ ) dbe_->save(out);
00188 }
00189
00190 void MuonTrackAnalyzer::analyze(const Event & event, const EventSetup& eventSetup){
00191
00192 LogDebug("MuonTrackAnalyzer") << "Run: " << event.id().run() << " Event: " << event.id().event();
00193
00194
00195 theService->update(eventSetup);
00196
00197 Handle<SimTrackContainer> simTracks;
00198 event.getByLabel("g4SimHits",simTracks);
00199 fillPlots(event,simTracks);
00200
00201
00202 if(doTracksAnalysis)
00203 tracksAnalysis(event,eventSetup,simTracks);
00204
00205 if(doSeedsAnalysis)
00206 seedsAnalysis(event,eventSetup,simTracks);
00207
00208
00209 }
00210
00211 void MuonTrackAnalyzer::seedsAnalysis(const Event & event, const EventSetup& eventSetup,
00212 Handle<SimTrackContainer> simTracks){
00213
00214 MuonPatternRecoDumper debug;
00215
00216
00217 Handle<TrajectorySeedCollection> seeds;
00218 event.getByLabel(theSeedsLabel, seeds);
00219
00220 LogTrace("MuonTrackAnalyzer")<<"Number of reconstructed seeds: " << seeds->size()<<endl;
00221
00222 for(TrajectorySeedCollection::const_iterator seed = seeds->begin();
00223 seed != seeds->end(); ++seed){
00224 TrajectoryStateOnSurface seedTSOS = getSeedTSOS(*seed);
00225 pair<SimTrack,double> sim = getSimTrack(seedTSOS,simTracks);
00226 fillPlots(seedTSOS, sim.first,
00227 hRecoSeedInner, debug);
00228
00229 std::pair<bool,FreeTrajectoryState> propSeed =
00230 theUpdator->propagateToNominalLine(seedTSOS);
00231 if(propSeed.first)
00232 fillPlots(propSeed.second, sim.first,
00233 hRecoSeedPCA, debug);
00234 else
00235 LogTrace("MuonTrackAnalyzer")<<"Error in seed propagation"<<endl;
00236
00237 }
00238 }
00239
00240
00241 void MuonTrackAnalyzer::tracksAnalysis(const Event & event, const EventSetup& eventSetup,
00242 Handle<SimTrackContainer> simTracks){
00243 MuonPatternRecoDumper debug;
00244
00245
00246
00247 Handle<reco::TrackCollection> tracks;
00248 event.getByLabel(theTracksLabel, tracks);
00249
00250 LogTrace("MuonTrackAnalyzer")<<"Reconstructed tracks: " << tracks->size() << endl;
00251 hNumberOfTracks->Fill(tracks->size());
00252
00253 if(tracks->size()) numberOfRecTracks++;
00254
00255
00256 for(reco::TrackCollection::const_iterator t = tracks->begin(); t != tracks->end(); ++t) {
00257
00258 reco::TransientTrack track(*t,&*theService->magneticField(),theService->trackingGeometry());
00259
00260 TrajectoryStateOnSurface outerTSOS = track.outermostMeasurementState();
00261 TrajectoryStateOnSurface innerTSOS = track.innermostMeasurementState();
00262 TrajectoryStateOnSurface pcaTSOS = track.impactPointState();
00263
00264 pair<SimTrack,double> sim = getSimTrack(pcaTSOS,simTracks);
00265 SimTrack simTrack = sim.first;
00266 hNumberOfTracksVsEta->Fill(simTrack.momentum().eta(), tracks->size());
00267 fillPlots(track,simTrack);
00268
00269 LogTrace("MuonTrackAnalyzer") << "State at the outer surface: " << endl;
00270 fillPlots(outerTSOS,simTrack,hRecoTracksOuter,debug);
00271
00272 LogTrace("MuonTrackAnalyzer") << "State at the inner surface: " << endl;
00273 fillPlots(innerTSOS,simTrack,hRecoTracksInner,debug);
00274
00275 LogTrace("MuonTrackAnalyzer") << "State at PCA: " << endl;
00276 fillPlots(pcaTSOS,simTrack,hRecoTracksPCA,debug);
00277
00278 double deltaPt_in_out = innerTSOS.globalMomentum().perp()-outerTSOS.globalMomentum().perp();
00279 hDeltaPt_In_Out_VsEta->Fill(simTrack.momentum().eta(),deltaPt_in_out);
00280
00281 double deltaPt_pca_sim = pcaTSOS.globalMomentum().perp()-sqrt(simTrack.momentum().Perp2());
00282 hDeltaPtVsEta->Fill(simTrack.momentum().eta(),deltaPt_pca_sim);
00283
00284 hChargeVsEta->Fill(simTrack.momentum().eta(),pcaTSOS.charge());
00285
00286 hChargeVsPt->Fill(sqrt(simTrack.momentum().perp2()),pcaTSOS.charge());
00287
00288 hPtRecVsPtGen->Fill(sqrt(simTrack.momentum().perp2()),pcaTSOS.globalMomentum().perp());
00289 }
00290 LogTrace("MuonTrackAnalyzer")<<"--------------------------------------------"<<endl;
00291 }
00292
00293
00294
00295
00296 void MuonTrackAnalyzer::fillPlots(const Event &event, edm::Handle<edm::SimTrackContainer> &simTracks){
00297
00298 if(!checkMuonSimHitPresence(event,simTracks)) return;
00299
00300
00301 SimTrackContainer::const_iterator simTrack;
00302 LogTrace("MuonTrackAnalyzer")<<"Simulated tracks: "<<simTracks->size()<<endl;
00303
00304 for (simTrack = simTracks->begin(); simTrack != simTracks->end(); ++simTrack)
00305 if (abs((*simTrack).type()) == 13) {
00306
00307 if( !isInTheAcceptance( (*simTrack).momentum().eta()) ) continue;
00308
00309 numberOfSimTracks++;
00310
00311 LogTrace("MuonTrackAnalyzer")<<"Simualted muon:"<<endl;
00312 LogTrace("MuonTrackAnalyzer")<<"Sim pT: "<<sqrt((*simTrack).momentum().perp2())<<endl;
00313 LogTrace("MuonTrackAnalyzer")<<"Sim Eta: "<<(*simTrack).momentum().eta()<<endl;
00314
00315 hSimTracks->Fill((*simTrack).momentum().mag(),
00316 sqrt((*simTrack).momentum().perp2()),
00317 (*simTrack).momentum().eta(),
00318 (*simTrack).momentum().phi(),
00319 -(*simTrack).type()/ abs((*simTrack).type()) );
00320 LogTrace("MuonTrackAnalyzer") << "hSimTracks filled" << endl;
00321 }
00322
00323 LogTrace("MuonTrackAnalyzer") << endl;
00324 }
00325
00326
00327 void MuonTrackAnalyzer::fillPlots(reco::TransientTrack &track, SimTrack &simTrack){
00328
00329 LogTrace("MuonTrackAnalyzer")<<"Analizer: New track, chi2: "<<track.chi2()<<" dof: "<<track.ndof()<<endl;
00330 hChi2->Fill(track.chi2());
00331 hDof->Fill(track.ndof());
00332 hChi2Norm->Fill(track.normalizedChi2());
00333 hHitsPerTrack->Fill(track.recHitsSize());
00334
00335 hChi2Prob->Fill( ChiSquaredProbability(track.chi2(),track.ndof()) );
00336
00337 hChi2VsEta->Fill(simTrack.momentum().eta(),track.chi2());
00338 hChi2NormVsEta->Fill(simTrack.momentum().eta(),track.normalizedChi2());
00339 hChi2ProbVsEta->Fill(simTrack.momentum().eta(),ChiSquaredProbability(track.chi2(),track.ndof()));
00340 hHitsPerTrackVsEta->Fill(simTrack.momentum().eta(),track.recHitsSize());
00341 hDofVsEta->Fill(simTrack.momentum().eta(),track.ndof());
00342 }
00343
00344
00345 void MuonTrackAnalyzer::fillPlots(TrajectoryStateOnSurface &recoTSOS,SimTrack &simTrack,
00346 HTrack *histo, MuonPatternRecoDumper &debug){
00347
00348 LogTrace("MuonTrackAnalyzer") << debug.dumpTSOS(recoTSOS)<<endl;
00349 histo->Fill(recoTSOS);
00350
00351 GlobalVector tsosVect = recoTSOS.globalMomentum();
00352 math::XYZVectorD reco(tsosVect.x(), tsosVect.y(), tsosVect.z());
00353 double deltaRVal = deltaR<double>(reco.eta(),reco.phi(),
00354 simTrack.momentum().eta(),simTrack.momentum().phi());
00355 histo->FillDeltaR(deltaRVal);
00356
00357 histo->computeResolutionAndPull(recoTSOS,simTrack);
00358 }
00359
00360
00361 void MuonTrackAnalyzer::fillPlots(FreeTrajectoryState &recoFTS,SimTrack &simTrack,
00362 HTrack *histo, MuonPatternRecoDumper &debug){
00363
00364 LogTrace("MuonTrackAnalyzer") << debug.dumpFTS(recoFTS)<<endl;
00365 histo->Fill(recoFTS);
00366
00367 GlobalVector ftsVect = recoFTS.momentum();
00368 math::XYZVectorD reco(ftsVect.x(), ftsVect.y(), ftsVect.z());
00369 double deltaRVal = deltaR<double>(reco.eta(),reco.phi(),
00370 simTrack.momentum().eta(),simTrack.momentum().phi());
00371 histo->FillDeltaR(deltaRVal);
00372
00373 histo->computeResolutionAndPull(recoFTS,simTrack);
00374 }
00375
00376 pair<SimTrack,double> MuonTrackAnalyzer::getSimTrack(TrajectoryStateOnSurface &tsos,
00377 Handle<SimTrackContainer> simTracks){
00378
00379
00380
00381
00382
00383
00384
00385
00386
00387
00388
00389
00390
00391
00392
00393
00394
00395 SimTrackContainer::const_iterator simTrack;
00396
00397 SimTrack result;
00398
00399 double bestDeltaR = 10e5;
00400 for (simTrack = simTracks->begin(); simTrack != simTracks->end(); ++simTrack){
00401 if (abs((*simTrack).type()) != 13) continue;
00402
00403
00404 GlobalVector tsosVect = tsos.globalMomentum();
00405 math::XYZVectorD vect(tsosVect.x(), tsosVect.y(), tsosVect.z());
00406 double newDeltaR = deltaR<double>(vect.eta(),vect.phi(),
00407 simTrack->momentum().eta(),simTrack->momentum().phi());
00408
00409 if ( newDeltaR < bestDeltaR ) {
00410 LogTrace("MuonTrackAnalyzer") << "Matching Track with DeltaR = " << newDeltaR<<endl;
00411 bestDeltaR = newDeltaR;
00412 result = *simTrack;
00413 }
00414 }
00415 return pair<SimTrack,double>(result,bestDeltaR);
00416 }
00417
00418
00419 bool MuonTrackAnalyzer::isInTheAcceptance(double eta){
00420 switch(theEtaRange){
00421 case all:
00422 return ( abs(eta) <= 2.4 ) ? true : false;
00423 case barrel:
00424 return ( abs(eta) < 1.1 ) ? true : false;
00425 case endcap:
00426 return ( abs(eta) >= 1.1 && abs(eta) <= 2.4 ) ? true : false;
00427 default:
00428 {LogTrace("MuonTrackAnalyzer")<<"No correct Eta range selected!! "<<endl; return false;}
00429 }
00430 }
00431
00432 bool MuonTrackAnalyzer::checkMuonSimHitPresence(const Event & event,
00433 edm::Handle<edm::SimTrackContainer> simTracks){
00434
00435
00436 Handle<PSimHitContainer> dtSimHits;
00437 event.getByLabel(theDTSimHitLabel.instance(),theDTSimHitLabel.label(), dtSimHits);
00438
00439 Handle<PSimHitContainer> cscSimHits;
00440 event.getByLabel(theCSCSimHitLabel.instance(),theCSCSimHitLabel.label(), cscSimHits);
00441
00442 Handle<PSimHitContainer> rpcSimHits;
00443 event.getByLabel(theRPCSimHitLabel.instance(),theRPCSimHitLabel.label(), rpcSimHits);
00444
00445 map<unsigned int, vector<const PSimHit*> > mapOfMuonSimHits;
00446
00447 for(PSimHitContainer::const_iterator simhit = dtSimHits->begin();
00448 simhit != dtSimHits->end(); ++simhit) {
00449 if (abs(simhit->particleType()) != 13) continue;
00450 mapOfMuonSimHits[simhit->trackId()].push_back(&*simhit);
00451 }
00452
00453 for(PSimHitContainer::const_iterator simhit = cscSimHits->begin();
00454 simhit != cscSimHits->end(); ++simhit) {
00455 if (abs(simhit->particleType()) != 13) continue;
00456 mapOfMuonSimHits[simhit->trackId()].push_back(&*simhit);
00457 }
00458
00459 for(PSimHitContainer::const_iterator simhit = rpcSimHits->begin();
00460 simhit != rpcSimHits->end(); ++simhit) {
00461 if (abs(simhit->particleType()) != 13) continue;
00462 mapOfMuonSimHits[simhit->trackId()].push_back(&*simhit);
00463 }
00464
00465 bool presence = false;
00466
00467 for (SimTrackContainer::const_iterator simTrack = simTracks->begin();
00468 simTrack != simTracks->end(); ++simTrack){
00469
00470 if (abs(simTrack->type()) != 13) continue;
00471
00472 map<unsigned int, vector<const PSimHit*> >::const_iterator mapIterator =
00473 mapOfMuonSimHits.find(simTrack->trackId());
00474
00475 if (mapIterator != mapOfMuonSimHits.end())
00476 presence = true;
00477 }
00478
00479 return presence;
00480 }
00481
00482 TrajectoryStateOnSurface MuonTrackAnalyzer::getSeedTSOS(const TrajectorySeed& seed){
00483
00484
00485 PTrajectoryStateOnDet pTSOD = seed.startingState();
00486
00487
00488
00489
00490 DetId seedDetId(pTSOD.detId());
00491
00492 const GeomDet* gdet = theService->trackingGeometry()->idToDet( seedDetId );
00493
00494 TrajectoryStateOnSurface initialState = trajectoryStateTransform::transientState(pTSOD, &(gdet->surface()), &*theService->magneticField());
00495
00496
00497 const DetLayer *initialLayer = theService->detLayerGeometry()->idToLayer( seedDetId );
00498
00499 PropagationDirection detLayerOrder = oppositeToMomentum;
00500
00501
00502 vector<const DetLayer*> detLayers;
00503 detLayers = initialLayer->compatibleLayers( *initialState.freeState(),detLayerOrder);
00504
00505 TrajectoryStateOnSurface result = initialState;
00506 if(detLayers.size()){
00507 const DetLayer* finalLayer = detLayers.back();
00508 const TrajectoryStateOnSurface propagatedState = theService->propagator(theSeedPropagatorName)->propagate(initialState, finalLayer->surface());
00509 if(propagatedState.isValid())
00510 result = propagatedState;
00511 }
00512
00513 return result;
00514 }