CMS 3D CMS Logo

TrackingMonitor.cc

Go to the documentation of this file.
00001 /*
00002  *  See header file for a description of this class.
00003  *
00004  *  $Date: 2008/12/11 11:57:11 $
00005  *  $Revision: 1.1 $
00006  *  \author Suchandra Dutta , Giorgia Mila
00007  */
00008 
00009 #include "DataFormats/TrackReco/interface/Track.h"
00010 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00011 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
00012 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
00013 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
00014 #include "MagneticField/Engine/interface/MagneticField.h"
00015 
00016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00017 #include "FWCore/ParameterSet/interface/InputTag.h"
00018 #include "DQMServices/Core/interface/DQMStore.h"
00019 #include "DQM/TrackingMonitor/interface/TrackAnalyzer.h"
00020 #include "DQM/TrackingMonitor/plugins/TrackingMonitor.h"
00021 
00022 #include <string>
00023 
00024 TrackingMonitor::TrackingMonitor(const edm::ParameterSet& iConfig) {
00025   dqmStore_ = edm::Service<DQMStore>().operator->();
00026   conf_ = iConfig;
00027 
00028   // the track analyzer
00029   theTrackAnalyzer = new TrackAnalyzer(conf_);
00030 }
00031 
00032 TrackingMonitor::~TrackingMonitor() { 
00033   delete theTrackAnalyzer;
00034 }
00035 
00036 void TrackingMonitor::beginJob(edm::EventSetup const& iSetup) {
00037 
00038   using namespace edm;
00039 
00040   std::string AlgoName     = conf_.getParameter<std::string>("AlgoName");
00041   std::string MEFolderName = conf_.getParameter<std::string>("FolderName"); 
00042 
00043   dqmStore_->setCurrentFolder(MEFolderName);
00044 
00045   int    TKNoBin = conf_.getParameter<int>("TkSizeBin");
00046   double TKNoMin = conf_.getParameter<double>("TkSizeMin");
00047   double TKNoMax = conf_.getParameter<double>("TkSizeMax");
00048 
00049   int    TKHitBin = conf_.getParameter<int>("RecHitBin");
00050   double TKHitMin = conf_.getParameter<double>("RecHitMin");
00051   double TKHitMax = conf_.getParameter<double>("RecHitMax");
00052 
00053   int    TKLayBin = conf_.getParameter<int>("RecLayBin");
00054   double TKLayMin = conf_.getParameter<double>("RecLayMin");
00055   double TKLayMax = conf_.getParameter<double>("RecLayMax");
00056 
00057   histname = "NumberOfTracks_";
00058   NumberOfTracks = dqmStore_->book1D(histname+AlgoName, histname+AlgoName, TKNoBin, TKNoMin, TKNoMax);
00059 
00060   histname = "NumberOfMeanRecHitsPerTrack_";
00061   NumberOfMeanRecHitsPerTrack = dqmStore_->book1D(histname+AlgoName, histname+AlgoName, TKHitBin, TKHitMin, TKHitMax);
00062   NumberOfMeanRecHitsPerTrack->setAxisTitle("Mean number of RecHits per track");
00063 
00064   histname = "NumberOfMeanLayersPerTrack_";
00065   NumberOfMeanLayersPerTrack = dqmStore_->book1D(histname+AlgoName, histname+AlgoName, TKLayBin, TKLayMin, TKLayMax);
00066   NumberOfMeanLayersPerTrack->setAxisTitle("Mean number of Layers per track");
00067 
00068   theTrackAnalyzer->beginJob(iSetup, dqmStore_);
00069  
00070 }
00071 
00072 //
00073 // -- Analyse
00074 //
00075 void TrackingMonitor::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
00076 
00077   using namespace edm;
00078   
00079   InputTag trackProducer = conf_.getParameter<edm::InputTag>("TrackProducer");
00080   
00081   Handle<reco::TrackCollection> trackCollection;
00082   iEvent.getByLabel(trackProducer, trackCollection);
00083   if (!trackCollection.isValid()) return;
00084 
00085   NumberOfTracks->Fill(trackCollection->size());
00086   
00087   int totalRecHits = 0, totalLayers = 0;
00088   for (reco::TrackCollection::const_iterator track = trackCollection->begin(); track!=trackCollection->end(); ++track) {
00089   
00090     totalRecHits += track->found();
00091     totalLayers += track->hitPattern().trackerLayersWithMeasurement();
00092 
00093     theTrackAnalyzer->analyze(iEvent, iSetup, *track);
00094   }
00095 
00096   double meanrechits = 0, meanlayers = 0;
00097   // check that track size to avoid division by zero.
00098   if (trackCollection->size()) {
00099     meanrechits = static_cast<double>(totalRecHits)/static_cast<double>(trackCollection->size());
00100     meanlayers = static_cast<double>(totalLayers)/static_cast<double>(trackCollection->size());
00101   }
00102   NumberOfMeanRecHitsPerTrack->Fill(meanrechits);
00103   NumberOfMeanLayersPerTrack->Fill(meanlayers);
00104 }
00105 
00106 
00107 void TrackingMonitor::endJob(void) {
00108   bool outputMEsInRootFile = conf_.getParameter<bool>("OutputMEsInRootFile");
00109   std::string outputFileName = conf_.getParameter<std::string>("OutputFileName");
00110   if(outputMEsInRootFile){
00111     dqmStore_->showDirStructure();
00112     dqmStore_->save(outputFileName);
00113   }
00114 
00115   
00116 }
00117 
00118 DEFINE_FWK_MODULE(TrackingMonitor);

Generated on Tue Jun 9 17:33:43 2009 for CMSSW by  doxygen 1.5.4