Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:18:45

0001 #ifndef AnalyzeTiming_h_
0002 #define AnalyzeTiming_h_
0003 
0004 #include <TH1F.h>
0005 #include <TFile.h>
0006 
0007 #include "DataFormats/HLTReco/interface/ModuleTiming.h"
0008 
0009 #include <map>
0010 #include <string>
0011 #include <iostream>
0012 
0013 // class for analyzing the timing of the trigger modules
0014 class AnalyzeTiming {
0015 private:
0016   // times are given in secs; convert to msecs here
0017   static const double secs2msecs;
0018 
0019 public:
0020   // create object with # of bins and range for histogram (in msecs)
0021   AnalyzeTiming(unsigned NBINS, double Xlow, double Xmax) : tot_time(0) {
0022     Nbins = NBINS;
0023     xlow = Xlow;
0024     xmax = Xmax;
0025     moduleTimes.clear();
0026     tot_time = createHistogram("full_event");
0027   }
0028   //
0029   ~AnalyzeTiming() {
0030     if (tot_time)
0031       delete tot_time;
0032 
0033     for (cIt it = moduleTimes.begin(); it != moduleTimes.end(); ++it)
0034       delete it->second;
0035     moduleTimes.clear();
0036   }
0037 
0038   // call this method for every event
0039   void analyze(const edm::EventTime &evtTime) {
0040     unsigned size = evtTime.size();
0041     for (unsigned int i = 0; i != size; ++i) {  // loop over all modules of event
0042 
0043       std::string module_name = evtTime.name(i);
0044       TH1F *tmp = 0;
0045       cIt it = moduleTimes.find(module_name);
0046       if (it != moduleTimes.end())
0047         tmp = it->second;
0048       else {
0049         // first time module is encountered; create new histogram
0050         tmp = createHistogram(module_name);
0051 
0052         // add entry to map
0053 
0054         /* this for some reason creates duplicate entries in moduleTimes... */
0055         //      moduleTimes[module_name] = tmp;
0056 
0057         moduleTimes[module_name.c_str()] = tmp;
0058       }
0059 
0060       tmp->Fill(evtTime.time(i) * secs2msecs);
0061 
0062     }  // loop over all modules of event
0063 
0064     tot_time->Fill(evtTime.tot_time() * secs2msecs);
0065   }
0066 
0067   // call this after all events have been processed
0068   void getResults() {
0069     TFile *out_file = new TFile("HLT_timing.root", "RECREATE");
0070     for (cIt it = moduleTimes.begin(); it != moduleTimes.end(); ++it)
0071       (it->second)->Write();
0072 
0073     tot_time->Write();
0074     //    out_file->Write();
0075     out_file->Close();
0076     delete out_file;
0077   }
0078 
0079 private:
0080   TH1F *createHistogram(std::string name) {
0081     std::cout << " Creating histogram for newly added module \"" << name << "\"" << std::endl;
0082     char title[1024];
0083     snprintf(title, 1024, "Processing times for %s", name.c_str());
0084     TH1F *newHistogram = new TH1F(name.c_str(), title, Nbins, xlow, xmax);
0085     newHistogram->GetXaxis()->SetTitle("t (ms)");
0086     // will do our own memory management
0087     newHistogram->SetDirectory(0);
0088     return newHistogram;
0089   }
0090 
0091   // time distributions
0092   // (a) per module
0093   std::map<std::string, TH1F *> moduleTimes;
0094   typedef std::map<std::string, TH1F *>::const_iterator cIt;
0095   // (b) per event
0096   TH1F *tot_time;
0097 
0098   // histogram parameters (common for all histograms)
0099   unsigned int Nbins;  // # of bins
0100   double xlow, xmax;   // histogram range
0101 };
0102 
0103 const double AnalyzeTiming::secs2msecs = 1000;
0104 
0105 #endif  // #define AnalyzeTiming_h_