CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/DQMOffline/Trigger/src/GeneralHLTOffline.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    GeneralHLTOffline
00004 // Class:      GeneralHLTOffline
00005 //
00012 //
00013 // Original Author:  Jason Michael Slaunwhite,512 1-008,`+41227670494,
00014 //         Created:  Fri Aug  5 10:34:47 CEST 2011
00015 // $Id: GeneralHLTOffline.cc,v 1.11 2013/02/02 16:24:44 rovere Exp $
00016 //
00017 //
00018 
00019 // system include files
00020 #include <memory>
00021 
00022 // user include files
00023 #include "DQMServices/Core/interface/DQMStore.h"
00024 #include "DQMServices/Core/interface/MonitorElement.h"
00025 
00026 #include "DataFormats/Common/interface/TriggerResults.h"
00027 #include "DataFormats/HLTReco/interface/TriggerObject.h"
00028 #include "DataFormats/HLTReco/interface/TriggerEvent.h"
00029 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
00030 
00031 #include "FWCore/Framework/interface/Frameworkfwd.h"
00032 #include "FWCore/Framework/interface/EDAnalyzer.h"
00033 #include "FWCore/Framework/interface/Event.h"
00034 #include "FWCore/Framework/interface/MakerMacros.h"
00035 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00036 #include "FWCore/ServiceRegistry/interface/Service.h"
00037 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00038 
00039 #include "HLTrigger/HLTcore/interface/HLTConfigProvider.h"
00040 
00041 #include "TMath.h"
00042 #include "TStyle.h"
00043 
00044 //
00045 // class declaration
00046 //
00047 
00048 class GeneralHLTOffline : public edm::EDAnalyzer {
00049  public:
00050   explicit GeneralHLTOffline(const edm::ParameterSet&);
00051   ~GeneralHLTOffline();
00052 
00053  private:
00054   virtual void beginJob();
00055   virtual void analyze(const edm::Event&, const edm::EventSetup&);
00056   virtual void endJob();
00057   virtual void beginRun(edm::Run const&, edm::EventSetup const&);
00058   virtual void endRun(edm::Run const&, edm::EventSetup const&);
00059   virtual void beginLuminosityBlock(edm::LuminosityBlock const&,
00060                                     edm::EventSetup const&);
00061   virtual void endLuminosityBlock(edm::LuminosityBlock const&,
00062                                   edm::EventSetup const&);
00063   virtual void setupHltMatrix(const std::string &, int);
00064   virtual void fillHltMatrix(const std::string &,
00065                              const std::string &,
00066                              double, double, bool);
00067 
00068   // ----------member data ---------------------------
00069 
00070 
00071   bool debugPrint;
00072   bool outputPrint;
00073   HLTConfigProvider hlt_config_;
00074 
00075   std::string plotDirectoryName;
00076   std::string hltTag;
00077   std::string hlt_menu_;
00078   std::vector< std::vector<std::string> > PDsVectorPathsVector;
00079   std::vector<std::string> AddedDatasets;
00080 
00081   DQMStore * dbe_;
00082   MonitorElement * cppath_;
00083 };
00084 
00085 //
00086 // constructors and destructor
00087 //
00088 GeneralHLTOffline::GeneralHLTOffline(const edm::ParameterSet& ps):hlt_menu_(""),
00089                                                                   dbe_(0),
00090                                                                   cppath_(0) {
00091   debugPrint  = false;
00092   outputPrint = false;
00093 
00094   plotDirectoryName = ps.getUntrackedParameter<std::string>("dirname",
00095                                                             "HLT/General");
00096 
00097   hltTag = ps.getParameter<std::string> ("HltProcessName");
00098 
00099   if (debugPrint) {
00100     std::cout << "Inside Constructor" << std::endl;
00101     std::cout << "Got plot dirname = " << plotDirectoryName << std::endl;
00102   }
00103 }
00104 
00105 
00106 GeneralHLTOffline::~GeneralHLTOffline() {
00107 }
00108 
00109 // ------------ method called for each event  ------------
00110 void
00111 GeneralHLTOffline::analyze(const edm::Event& iEvent,
00112                            const edm::EventSetup& iSetup) {
00113   if (debugPrint)
00114     std::cout << "Inside analyze - run, block, event "
00115               << iEvent.id().run() << " , " << iEvent.id().luminosityBlock()
00116               << " , " << iEvent.id() << " , " << std::endl;
00117 
00118 
00119   // Access Trigger Results
00120   edm::Handle<edm::TriggerResults> triggerResults;
00121   iEvent.getByLabel(edm::InputTag("TriggerResults", "", hltTag), triggerResults);
00122 
00123   if (!triggerResults.isValid()) {
00124     if (debugPrint)
00125       std::cout << "Trigger results not valid" << std::endl;
00126     return;
00127   }
00128 
00129   if (debugPrint)
00130     std::cout << "Found triggerResults" << std::endl;
00131 
00132   edm::Handle<trigger::TriggerEvent> aodTriggerEvent;
00133   iEvent.getByLabel(edm::InputTag("hltTriggerSummaryAOD", "", hltTag),
00134                     aodTriggerEvent);
00135 
00136   if (!aodTriggerEvent.isValid()) {
00137     if (debugPrint)
00138       std::cout << "No AOD trigger summary found! Returning...";
00139     return;
00140   }
00141 
00142   const std::vector<std::string> &nameStreams = hlt_config_.streamNames();
00143 
00144   const trigger::TriggerObjectCollection objects = aodTriggerEvent->getObjects();
00145 
00146   bool streamAfound = false;
00147   std::vector<std::string>::const_iterator si = nameStreams.begin();
00148   std::vector<std::string>::const_iterator se = nameStreams.end();
00149   for ( ; si != se; ++si) {
00150     if ((*si) == "A") {
00151       streamAfound = true;
00152       break;
00153     }
00154   }
00155 
00156   if (debugPrint) {
00157     std::cout << "Stream A "
00158               << (streamAfound ? "found" : "not found")
00159               << std::endl;
00160   }
00161 
00162   if (streamAfound) {
00163     const std::vector<std::string> &datasetNames =  hlt_config_.streamContent("A");
00164     // Loop over PDs
00165     for (unsigned int iPD = 0; iPD < datasetNames.size(); iPD++) {
00166       // Loop over Paths in each PD
00167       bool first_count = true;
00168       for (unsigned int iPath = 0;
00169            iPath < PDsVectorPathsVector[iPD].size(); iPath++) {
00170         std::string &pathName = PDsVectorPathsVector[iPD][iPath];
00171         unsigned int index = hlt_config_.triggerIndex(pathName);
00172         if (debugPrint) {
00173           std::cout << "Looking at path " << pathName << std::endl;
00174           std::cout << "Index = " << index
00175                     << " triggerResults->size() = " << triggerResults->size()
00176                     << std::endl;
00177         }
00178 
00179         // fill the histos with empty weights......
00180         const std::string &label = datasetNames[iPD];
00181         std::string fullPathToCPP = "HLT/GeneralHLTOffline/"
00182             + label + "/cppath_" + label + hlt_menu_;
00183         MonitorElement * ME_mini_cppath = dbe_->get(fullPathToCPP);
00184         TH1F * hist_mini_cppath = NULL;
00185         if (ME_mini_cppath)
00186           hist_mini_cppath = ME_mini_cppath->getTH1F();
00187 
00188         if (hist_mini_cppath) {
00189           TAxis * axis = hist_mini_cppath->GetXaxis();
00190           if (axis) {
00191             int bin_num = axis->FindBin(pathName.c_str());
00192             int bn = bin_num - 1;
00193             hist_mini_cppath->Fill(bn, 0);
00194             hist_mini_cppath->SetEntries(hist_mini_cppath->Integral());
00195           }
00196         }
00197 
00198         if (index < triggerResults->size()) {
00199           if (triggerResults->accept(index)) {
00200             cppath_->Fill(index, 1);
00201             if (debugPrint)
00202               std::cout << "Check Event " <<  iEvent.id()
00203                         << " Run " << iEvent.id().run()
00204                         << " fired path " << pathName << std::endl;
00205 
00206             // look up module labels for this path
00207             const std::vector<std::string> &modulesThisPath =
00208                 hlt_config_.moduleLabels(pathName);
00209 
00210             if (debugPrint)
00211               std::cout << "Looping over module labels " << std::endl;
00212 
00213             // Loop backward through module names
00214             for (int iModule = (modulesThisPath.size() - 1);
00215                  iModule >= 0; iModule--) {
00216               if (debugPrint)
00217                 std::cout << "Module name is "
00218                           << modulesThisPath[iModule] << std::endl;
00219               // check to see if you have savetags information
00220               if (hlt_config_.saveTags(modulesThisPath[iModule])) {
00221                 if (debugPrint)
00222                   std::cout << "For path " << pathName
00223                             << " this module " << modulesThisPath[iModule]
00224                             <<" is a saveTags module of type "
00225                             << hlt_config_.moduleType(modulesThisPath[iModule])
00226                             << std::endl;
00227                 if (hlt_config_.moduleType(modulesThisPath[iModule])
00228                     == "HLTLevel1GTSeed")
00229                   break;
00230                 edm::InputTag moduleWhoseResultsWeWant(modulesThisPath[iModule],
00231                                                        "",
00232                                                        hltTag);
00233                 unsigned int idx_module_aod_trg =
00234                     aodTriggerEvent->filterIndex(moduleWhoseResultsWeWant);
00235                 if (idx_module_aod_trg < aodTriggerEvent->sizeFilters()) {
00236                   const trigger::Keys &keys =
00237                       aodTriggerEvent->filterKeys(idx_module_aod_trg);
00238                   if (debugPrint)
00239                     std::cout << "Got Keys for index "
00240                               << idx_module_aod_trg
00241                               <<", size of keys is " << keys.size()
00242                               << std::endl;
00243                   if (keys.size() >= 1000)
00244                     edm::LogWarning("GeneralHLTOffline")
00245                         << "WARNING!! size of keys is " << keys.size()
00246                         << " for path " << pathName << " and module "
00247                         << modulesThisPath[iModule]<< std::endl;
00248 
00249                   // There can be > 100 keys (3-vectors) for some
00250                   // modules with no ID filled the first one has the
00251                   // highest value for single-object triggers for
00252                   // multi-object triggers, seems reasonable to use
00253                   // the first one as well So loop here has been
00254                   // commented out for ( size_t iKey = 0; iKey <
00255                   // keys.size(); iKey++ ) {
00256 
00257                   if (keys.size() > 0) {
00258                     trigger::TriggerObject foundObject = objects[keys[0]];
00259                     if (debugPrint || outputPrint)
00260                       std::cout << "This object has id (pt, eta, phi) = "
00261                                 << " " << foundObject.id() << " "
00262                                 << std::setw(10) << foundObject.pt()
00263                                 << ", " << std::setw(10) << foundObject.eta()
00264                                 << ", " << std::setw(10) << foundObject.phi()
00265                                 << "   for path = " << std::setw(20) << pathName
00266                                 << " module " << std::setw(40)
00267                                 << modulesThisPath[iModule] << std::endl;
00268                     if (debugPrint)
00269                       std::cout << "CHECK RUN " << iEvent.id().run() << " "
00270                                 << iEvent.id() << " " << pathName << " "
00271                                 << modulesThisPath[iModule] << " "
00272                                 << datasetNames[iPD] << " "
00273                                 << hlt_config_.moduleType(modulesThisPath[iModule])
00274                                 << " " << keys.size() << " "
00275                                 << std::setprecision(4) << foundObject.pt() << " "
00276                                 << foundObject.eta() << " "
00277                                 << foundObject.phi() << std::endl;
00278 
00279                     // first_count is to make sure that the top-level
00280                     // histograms of each dataset don't get filled
00281                     // more than once
00282                     fillHltMatrix(datasetNames[iPD], pathName,
00283                                   foundObject.eta(), foundObject.phi(),
00284                                   first_count);
00285                     first_count = false;
00286                   }  // at least one key
00287                 }  // end if filter in aodTriggerEvent
00288                 // OK, we found the last module. No need to look at
00289                 // the others.  get out of the loop
00290                 break;
00291               }  // end if saveTags
00292             }  // end Loop backward through module names
00293           }  // end if(triggerResults->accept(index))
00294         }  // end if (index < triggerResults->size())
00295       }  // end Loop over Paths in each PD
00296     }  // end Loop over PDs
00297   }
00298 }
00299 
00300 
00301 // ------------ method called once each job just before starting event loop  ------------
00302 void
00303 GeneralHLTOffline::beginJob() {
00304   if (debugPrint)
00305     std::cout << "Inside begin job" << std::endl;
00306 
00307   dbe_ = edm::Service<DQMStore>().operator->();
00308   if (dbe_)
00309     dbe_->setCurrentFolder(plotDirectoryName);
00310 }
00311 
00312 // ------------ method called once each job just after ending the event loop  ------------
00313 void
00314 GeneralHLTOffline::endJob() {
00315 }
00316 
00317 // ------------ method called when starting to processes a run  ------------
00318 void
00319 GeneralHLTOffline::beginRun(edm::Run const& iRun,
00320                             edm::EventSetup const& iSetup) {
00321   if (debugPrint)
00322     std::cout << "Inside beginRun" << std::endl;
00323 
00324   bool changed = true;
00325   hlt_menu_ = hlt_config_.tableName();
00326   for (unsigned int n = 0, e = hlt_menu_.length(); n != e; ++n)
00327     if (hlt_menu_[n] == '/' || hlt_menu_[n] == '.')
00328       hlt_menu_[n] = '_';
00329 
00330   if (!hlt_config_.init(iRun, iSetup, hltTag, changed)) {
00331     if (debugPrint)
00332       std::cout << "Warning, didn't find process HLT" << std::endl;
00333   } else {
00334     if (debugPrint)
00335       std::cout << " HLTConfig processName " << hlt_config_.processName()
00336                 << " tableName " << hlt_config_.tableName()
00337                 << " size " << hlt_config_.size() << std::endl;
00338   }
00339 
00341 
00342   dbe_->setCurrentFolder("HLT/GeneralHLTOffline/");
00343   cppath_ = dbe_->book1D("cppath" + hlt_menu_,
00344                          "Counts/Path",
00345                          hlt_config_.size(), 0, hlt_config_.size());
00346 
00347   const std::vector<std::string> &nameStreams = hlt_config_.streamNames();
00348   bool streamAfound = false;
00349   std::vector<std::string>::const_iterator si = nameStreams.begin();
00350   std::vector<std::string>::const_iterator se = nameStreams.end();
00351   for ( ; si != se; ++si) {
00352     if ((*si) == "A") {
00353       streamAfound = true;
00354       break;
00355     }
00356   }
00357 
00358   if (streamAfound) {
00359     const std::vector<std::string> &datasetNames =  hlt_config_.streamContent("A");
00360     if (debugPrint)
00361       std::cout << "Number of Stream A datasets "
00362                 << datasetNames.size() << std::endl;
00363 
00364     for (unsigned int i = 0; i < datasetNames.size(); i++) {
00365       const std::vector<std::string> &datasetPaths = hlt_config_.datasetContent(datasetNames[i]);
00366       if (debugPrint) {
00367         std::cout << "This is dataset " << datasetNames[i]
00368                   << "datasetPaths.size() = " << datasetPaths.size() << std::endl;
00369         for (unsigned int iPath = 0;
00370              iPath < datasetPaths.size(); iPath++) {
00371           std::cout << "Before setupHltMatrix -  MET dataset "
00372                     << datasetPaths[iPath] << std::endl;
00373         }
00374       }
00375       // Check if dataset has been added - if not add it
00376       // need to loop through AddedDatasets and compare
00377       bool foundDataset = false;
00378       int datasetNum = -1;
00379       for (unsigned int d = 0; d < AddedDatasets.size(); d++) {
00380         if (AddedDatasets[d].compare(datasetNames[i]) == 0) {
00381           foundDataset = true;
00382           datasetNum = d;
00383           if (debugPrint)
00384             std::cout << "Dataset " << datasetNames[i]
00385                       << " found in AddedDatasets at position " << d << std::endl;
00386           break;
00387         }
00388       }
00389 
00390       if (!foundDataset) {
00391         if (debugPrint)
00392           std::cout << " Fill trigger paths for dataset "
00393                     << datasetNames[i] << std::endl;
00394         PDsVectorPathsVector.push_back(datasetPaths);
00395         // store dataset pathname
00396         AddedDatasets.push_back(datasetNames[i]);
00397       } else {
00398         // This trigger path has already been added - this implies that
00399         // this is a new run What we want to do is check if there is a
00400         // new trigger that was not in the original dataset For a given
00401         // dataset, loop over the stored list of triggers, and compare
00402         // to the current list of triggers If any of the triggers are
00403         // missing, add them to the end of the appropriate dataset
00404         if (debugPrint)
00405           std::cout << " Additional runs : Check for additional"
00406                     << "trigger paths per dataset " << std::endl;
00407         // Loop over correct path of PDsVectorPathsVector
00408         bool found = false;
00409 
00410         // Loop over triggers in the path
00411         for (unsigned int iTrig = 0; iTrig < datasetPaths.size(); iTrig++) {
00412           if (debugPrint)
00413             std::cout << "Looping over trigger list in dataset "
00414                       <<  iTrig <<  "  "
00415                       << datasetPaths[iTrig] << std::endl;
00416           found = false;
00417           // Loop over triggers already on the list
00418           for (unsigned int od = 0; od < PDsVectorPathsVector[datasetNum].size(); od++) {
00419             if (debugPrint)
00420               std::cout << "Looping over existing trigger list " << od
00421                         <<  "  " << PDsVectorPathsVector[datasetNum][od] << std::endl;
00422             // Compare, see if match is found
00423             if (hlt_config_.removeVersion(datasetPaths[iTrig]).compare(
00424                     hlt_config_.removeVersion(PDsVectorPathsVector[datasetNum][od])) == 0) {
00425               found = true;
00426               if (debugPrint)
00427                 std::cout << " FOUND " << datasetPaths[iTrig] << std::endl;
00428               break;
00429             }
00430           }
00431           // If match is not found, add trigger to correct path of PDsVectorPathsVector
00432           if (!found)
00433             PDsVectorPathsVector[datasetNum].push_back(datasetPaths[iTrig]);
00434           if (debugPrint)
00435             std::cout << datasetPaths[iTrig]
00436                       << "  NOT FOUND - so we added it to the correct dataset "
00437                       << datasetNames[i] << std::endl;
00438         }
00439       }
00440       // Let's check this whole big structure
00441       if (debugPrint) {
00442         for (unsigned int is = 0; is < PDsVectorPathsVector.size(); is++) {
00443           std::cout << "   PDsVectorPathsVector[" << is << "] is "
00444                     << PDsVectorPathsVector[is].size() << std::endl;
00445           for (unsigned int ip = 0; ip < PDsVectorPathsVector[is].size(); ip++) {
00446             std::cout << "    trigger " << ip << " path "
00447                       << PDsVectorPathsVector[is][ip] << std::endl;
00448           }
00449         }
00450       }
00451 
00452       if (debugPrint)
00453         std::cout <<"Found PD: " << datasetNames[i] << std::endl;
00454 
00455       setupHltMatrix(datasetNames[i], i);
00456     }  // end of loop over dataset names
00457   }  // if stream A found
00458 }  // end of beginRun
00459 
00460 // ------------ method called when ending the processing of a run  ------------
00461 void GeneralHLTOffline::endRun(edm::Run const&, edm::EventSetup const&) {
00462   if (debugPrint)
00463     std::cout << " endRun called " << std::endl;
00464 }
00465 
00466 
00467 void GeneralHLTOffline::setupHltMatrix(const std::string & label, int iPD) {
00468   std::string h_name;
00469   std::string h_title;
00470   std::string h_name_1dEta;
00471   std::string h_name_1dPhi;
00472   std::string h_title_1dEta;
00473   std::string h_title_1dPhi;
00474   std::string h_name_1dEtaPath;
00475   std::string h_name_1dPhiPath;
00476   std::string h_title_1dEtaPath;
00477   std::string h_title_1dPhiPath;
00478   std::string pathName;
00479   std::string PD_Folder;
00480   std::string Path_Folder;
00481 
00482   PD_Folder = TString("HLT/GeneralHLTOffline");
00483   if (label != "SingleMu" && label != "SingleElectron" && label != "Jet")
00484     PD_Folder = TString("HLT/GeneralHLTOffline/"+label);
00485 
00486   dbe_->setCurrentFolder(PD_Folder.c_str());
00487   dbe_->bookString("hltMenuName", hlt_menu_.c_str());
00488 
00489   h_name = "HLT_" +label + "_EtaVsPhi";
00490   h_title = "HLT_" + label + "_EtaVsPhi";
00491   h_name_1dEta = "HLT_" + label + "_1dEta";
00492   h_name_1dPhi = "HLT_" + label + "_1dPhi";
00493   h_title_1dEta = label + " Occupancy Vs Eta";
00494   h_title_1dPhi = label + " Occupancy Vs Phi";
00495 
00496   Int_t numBinsEta = 30;
00497   Int_t numBinsPhi = 34;
00498   Int_t numBinsEtaFine = 60;
00499   Int_t numBinsPhiFine = 66;
00500   Double_t EtaMax = 2.610;
00501   Double_t PhiMax = 17.0*TMath::Pi()/16.0;
00502   Double_t PhiMaxFine = 33.0*TMath::Pi()/32.0;
00503   MonitorElement * service_me = NULL;
00504 
00505   service_me = dbe_->book2D(h_name.c_str(),
00506                             h_title.c_str(),
00507                             numBinsEta, -EtaMax, EtaMax,
00508                             numBinsPhi, -PhiMax, PhiMax);
00509   if (TH1 * service_histo = service_me->getTH2F())
00510     service_histo->SetMinimum(0);
00511 
00512   if (label != "MET" && label != "HT") {
00513     service_me = dbe_->book1D(h_name_1dEta.c_str(),
00514                               h_title_1dEta.c_str(),
00515                               numBinsEtaFine, -EtaMax, EtaMax);
00516     if (TH1 * service_histo = service_me->getTH1F())
00517       service_histo->SetMinimum(0);
00518   }
00519   if (label != "HT") {
00520     service_me = dbe_->book1D(h_name_1dPhi.c_str(),
00521                               h_title_1dPhi.c_str(),
00522                               numBinsPhiFine, -PhiMaxFine, PhiMaxFine);
00523     if (TH1 * service_histo = service_me->getTH1F())
00524       service_histo->SetMinimum(0);
00525   }
00526 
00527   // make it the top level directory, that is on the same dir level as
00528   // paths
00529   std::string folderz;
00530   folderz = TString("HLT/GeneralHLTOffline/"+label);
00531   dbe_->setCurrentFolder(folderz.c_str());
00532 
00533   std::string dnamez = "cppath_" + label + "_" + hlt_menu_;
00534   int sizez = PDsVectorPathsVector[iPD].size();
00535   TH1F * hist_mini_cppath = NULL;
00536   MonitorElement * hist_mini_cppath_me = dbe_->book1D(dnamez.c_str(),
00537                                                       dnamez.c_str(),
00538                                                       sizez,
00539                                                       0,
00540                                                       sizez);
00541   if (hist_mini_cppath_me)
00542     hist_mini_cppath = hist_mini_cppath_me->getTH1F();
00543 
00544   unsigned int jPath;
00545   for (unsigned int iPath = 0; iPath < PDsVectorPathsVector[iPD].size(); iPath++) {
00546     pathName = hlt_config_.removeVersion(PDsVectorPathsVector[iPD][iPath]);
00547     h_name_1dEtaPath = "HLT_" + pathName + "_1dEta";
00548     h_name_1dPhiPath = "HLT_" + pathName + "_1dPhi";
00549     h_title_1dEtaPath = pathName + " Occupancy Vs Eta";
00550     h_title_1dPhiPath = pathName + "Occupancy Vs Phi";
00551     jPath = iPath + 1;
00552 
00553     if (hist_mini_cppath) {
00554       TAxis * axis = hist_mini_cppath->GetXaxis();
00555       if (axis)
00556         axis->SetBinLabel(jPath, pathName.c_str());
00557     }
00558 
00559     Path_Folder = TString("HLT/GeneralHLTOffline/" + label + "/Paths");
00560     dbe_->setCurrentFolder(Path_Folder.c_str());
00561 
00562     dbe_->book1D(h_name_1dEtaPath.c_str(),
00563                  h_title_1dEtaPath.c_str(),
00564                  numBinsEtaFine, -EtaMax, EtaMax);
00565     dbe_->book1D(h_name_1dPhiPath.c_str(),
00566                  h_title_1dPhiPath.c_str(),
00567                  numBinsPhiFine, -PhiMaxFine, PhiMaxFine);
00568 
00569     if (debugPrint)
00570       std::cout << "book1D for " << pathName << std::endl;
00571   }
00572 
00573   if (debugPrint)
00574     std::cout << "Success setupHltMatrix( " << label << " , "
00575               << iPD << " )" << std::cout;
00576 }  // End setupHltMatrix
00577 
00578 
00579 void GeneralHLTOffline::fillHltMatrix(const std::string & label,
00580                                       const std::string & path,
00581                                       double Eta,
00582                                       double Phi,
00583                                       bool first_count) {
00584   if (debugPrint)
00585     std::cout << "Inside fillHltMatrix( " << label << " , "
00586               << path << " ) " << std::endl;
00587 
00588   std::string fullPathToME;
00589   std::string fullPathToME1dEta;
00590   std::string fullPathToME1dPhi;
00591   std::string fullPathToME1dEtaPath;
00592   std::string fullPathToME1dPhiPath;
00593   std::string fullPathToCPP;
00594 
00595 
00596   fullPathToME = "HLT/GeneralHLTOffline/HLT_" + label + "_EtaVsPhi";
00597   fullPathToME1dEta = "HLT/GeneralHLTOffline/HLT_" + label + "_1dEta";
00598   fullPathToME1dPhi = "HLT/GeneralHLTOffline/HLT_" + label + "_1dPhi";
00599   fullPathToCPP = "HLT/GeneralHLTOffline/" + label
00600       + "/cppath_" + label + "_" + hlt_menu_;
00601 
00602   if (label != "SingleMu" && label != "SingleElectron" && label != "Jet") {
00603     fullPathToME = "HLT/GeneralHLTOffline/"
00604         + label + "/HLT_" + label + "_EtaVsPhi";
00605     fullPathToME1dEta = "HLT/GeneralHLTOffline/"
00606         + label + "/HLT_" + label + "_1dEta";
00607     fullPathToME1dPhi = "HLT/GeneralHLTOffline/"
00608         + label + "/HLT_" + label + "_1dPhi";
00609   }
00610 
00611   fullPathToME1dEtaPath = "HLT/GeneralHLTOffline/"
00612       + label + "/Paths/HLT_"
00613       + hlt_config_.removeVersion(path) + "_1dEta";
00614   fullPathToME1dPhiPath = "HLT/GeneralHLTOffline/"
00615       + label + "/Paths/HLT_"
00616       + hlt_config_.removeVersion(path) + "_1dPhi";
00617 
00618   TH1F * hist_mini_cppath = NULL;
00619   MonitorElement * ME_mini_cppath = dbe_->get(fullPathToCPP);
00620   if (ME_mini_cppath)
00621     hist_mini_cppath = ME_mini_cppath->getTH1F();
00622 
00623   // fill top-level histograms
00624   if (first_count) {
00625     if (debugPrint)
00626       std::cout << " label " << label << " fullPathToME1dPhi "
00627                 << fullPathToME1dPhi << " path "  << path
00628                 << " Phi " << Phi << " Eta " << Eta << std::endl;
00629 
00630     if (label != "MET" && label != "HT") {
00631       MonitorElement * ME_1dEta = dbe_->get(fullPathToME1dEta);
00632       if (ME_1dEta) {
00633         TH1F * hist_1dEta = ME_1dEta->getTH1F();
00634         if (hist_1dEta)
00635           hist_1dEta->Fill(Eta);
00636       }
00637     }
00638     if (label != "HT") {
00639       MonitorElement * ME_1dPhi = dbe_->get(fullPathToME1dPhi);
00640       if (ME_1dPhi) {
00641         TH1F * hist_1dPhi = ME_1dPhi->getTH1F();
00642         if (hist_1dPhi)
00643           hist_1dPhi->Fill(Phi);
00644         if (debugPrint)
00645           std::cout << "  **FILLED** label " << label << " fullPathToME1dPhi "
00646                     << fullPathToME1dPhi << " path "  << path
00647                     << " Phi " << Phi << " Eta " << Eta << std::endl;
00648       }
00649     }
00650     if (label != "MET" && label != "HT") {
00651       MonitorElement * ME_2d = dbe_->get(fullPathToME);
00652       if (ME_2d) {
00653         TH2F * hist_2d = ME_2d->getTH2F();
00654         if (hist_2d)
00655           hist_2d->Fill(Eta, Phi);
00656       }
00657     }
00658   }  // end fill top-level histograms
00659 
00660   if (label != "MET" && label != "HT") {
00661     MonitorElement * ME_1dEtaPath = dbe_->get(fullPathToME1dEtaPath);
00662     if (ME_1dEtaPath) {
00663       TH1F * hist_1dEtaPath = ME_1dEtaPath->getTH1F();
00664       if (hist_1dEtaPath)
00665         hist_1dEtaPath->Fill(Eta);
00666     }
00667   }
00668   if (label != "HT") {
00669     MonitorElement * ME_1dPhiPath = dbe_->get(fullPathToME1dPhiPath);
00670     if (ME_1dPhiPath) {
00671       TH1F * hist_1dPhiPath = ME_1dPhiPath->getTH1F();
00672       if (hist_1dPhiPath)
00673         hist_1dPhiPath->Fill(Phi);
00674     }
00675   }
00676 
00677   if (debugPrint)
00678     if (label == "MET")
00679       std::cout << " MET Eta is " << Eta << std::endl;
00680 
00681   if (hist_mini_cppath) {
00682     TAxis * axis = hist_mini_cppath->GetXaxis();
00683     int bin_num = axis->FindBin(path.c_str());
00684     int bn = bin_num - 1;
00685     hist_mini_cppath->Fill(bn, 1);
00686   }
00687 
00688   if (debugPrint)
00689     std::cout << "hist->Fill" << std::endl;
00690 }  // End fillHltMatrix
00691 
00692 void GeneralHLTOffline::beginLuminosityBlock(edm::LuminosityBlock const&,
00693                                              edm::EventSetup const&) {
00694 }
00695 
00696 void
00697 GeneralHLTOffline::endLuminosityBlock(edm::LuminosityBlock const&,
00698                                       edm::EventSetup const&) {
00699 }
00700 
00701 DEFINE_FWK_MODULE(GeneralHLTOffline);