CMS 3D CMS Logo

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