CMS 3D CMS Logo

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