Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:20:13

0001 // -*- C++ -*-
0002 //
0003 // Package:     GctErrorAnalyzer
0004 // Class:   GctErrorAnalyzer
0005 //
0006 /**\class GctErrorAnalyzer GctErrorAnalyzer.cc L1Trigger/L1GctAnalyzer/src/GctErrorAnalyzer.cc
0007 
0008 Description: Tool to debug the GCT with useful output
0009 
0010 Implementation:
0011 <Notes on implementation>
0012 */
0013 //
0014 // Original Author:  Jad Marrouche
0015 //         Created:  Wed May 20 14:19:23 CEST 2009
0016 //
0017 //
0018 
0019 // system include files
0020 #include <memory>
0021 
0022 // user include files
0023 #include "FWCore/Framework/interface/Frameworkfwd.h"
0024 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0025 #include "FWCore/Framework/interface/Event.h"
0026 #include "FWCore/Framework/interface/MakerMacros.h"
0027 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0028 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0029 #include "FWCore/Utilities/interface/InputTag.h"
0030 #include "FWCore/ServiceRegistry/interface/Service.h"
0031 //TFile maker include
0032 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0033 //ROOT includes
0034 #include "TH1.h"
0035 #include "TH2.h"
0036 #include "TAxis.h"
0037 //RCT and GCT DataFormat Collections
0038 #include "DataFormats/L1GlobalCaloTrigger/interface/L1GctCollections.h"
0039 #include "DataFormats/L1CaloTrigger/interface/L1CaloCollections.h"
0040 //GctErrorAnalyzer includes
0041 #include "L1Trigger/L1GctAnalyzer/interface/GctErrorAnalyzerDefinitions.h"
0042 #include "L1Trigger/L1GctAnalyzer/interface/compareCands.h"
0043 #include "L1Trigger/L1GctAnalyzer/interface/compareRingSums.h"
0044 #include "L1Trigger/L1GctAnalyzer/interface/compareBitCounts.h"
0045 #include "L1Trigger/L1GctAnalyzer/interface/compareTotalEnergySums.h"
0046 #include "L1Trigger/L1GctAnalyzer/interface/compareMissingEnergySums.h"
0047 //STL includes
0048 #include <string>
0049 #include <vector>
0050 #include <sstream>
0051 #include <algorithm>
0052 
0053 //
0054 // class declaration
0055 //
0056 
0057 class GctErrorAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources> {
0058 public:
0059   GctErrorAnalyzer() = delete;
0060   GctErrorAnalyzer(const GctErrorAnalyzer &) = delete;
0061   GctErrorAnalyzer operator=(const GctErrorAnalyzer &) = delete;
0062 
0063 private:
0064   void plotRCTRegions(const edm::Handle<L1CaloRegionCollection> &caloRegions);
0065   void plotIsoEm(const edm::Handle<L1GctEmCandCollection> &data, const edm::Handle<L1GctEmCandCollection> &emu);
0066   void plotNonIsoEm(const edm::Handle<L1GctEmCandCollection> &data, const edm::Handle<L1GctEmCandCollection> &emu);
0067   void plotEGErrors(const edm::Handle<L1GctEmCandCollection> &dataiso,
0068                     const edm::Handle<L1GctEmCandCollection> &emuiso,
0069                     const edm::Handle<L1GctEmCandCollection> &datanoniso,
0070                     const edm::Handle<L1GctEmCandCollection> &emunoniso,
0071                     const edm::Handle<L1CaloEmCollection> &regions);
0072   void plotCenJets(const edm::Handle<L1GctJetCandCollection> &data, const edm::Handle<L1GctJetCandCollection> &emu);
0073   void plotTauJets(const edm::Handle<L1GctJetCandCollection> &data, const edm::Handle<L1GctJetCandCollection> &emu);
0074   void plotForJets(const edm::Handle<L1GctJetCandCollection> &data, const edm::Handle<L1GctJetCandCollection> &emu);
0075   void plotIntJets(const edm::Handle<L1GctInternJetDataCollection> &emu);
0076   static bool sortJets(const jetData &jet1,
0077                        const jetData &jet2);  //define this as static as it doesn't need a GctErrorAnalyzer object
0078   void plotJetErrors(const edm::Handle<L1GctJetCandCollection> &cendata,
0079                      const edm::Handle<L1GctJetCandCollection> &cenemu,
0080                      const edm::Handle<L1GctJetCandCollection> &taudata,
0081                      const edm::Handle<L1GctJetCandCollection> &tauemu,
0082                      const edm::Handle<L1GctJetCandCollection> &fordata,
0083                      const edm::Handle<L1GctJetCandCollection> &foremu,
0084                      const edm::Handle<L1CaloRegionCollection> &regions);
0085   void plotHFRingSums(const edm::Handle<L1GctHFRingEtSumsCollection> &data,
0086                       const edm::Handle<L1GctHFRingEtSumsCollection> &emu);
0087   void plotHFBitCounts(const edm::Handle<L1GctHFBitCountsCollection> &hfBitCountsD,
0088                        const edm::Handle<L1GctHFBitCountsCollection> &hfBitCountsE);
0089   void plotHFErrors(const edm::Handle<L1GctHFRingEtSumsCollection> &hfRingSumsD,
0090                     const edm::Handle<L1GctHFRingEtSumsCollection> &hfRingSumsE,
0091                     const edm::Handle<L1GctHFBitCountsCollection> &hfBitCountsD,
0092                     const edm::Handle<L1GctHFBitCountsCollection> &hfBitCountsE,
0093                     const edm::Handle<L1CaloRegionCollection> &caloRegions);
0094   void plotTotalE(const edm::Handle<L1GctEtTotalCollection> &totalEtD,
0095                   const edm::Handle<L1GctEtTotalCollection> &totalEtE);
0096   void plotTotalH(const edm::Handle<L1GctEtHadCollection> &totalHtD, const edm::Handle<L1GctEtHadCollection> &totalHtE);
0097   void plotTotalEErrors(const edm::Handle<L1GctEtTotalCollection> &totalEtD,
0098                         const edm::Handle<L1GctEtTotalCollection> &totalEtE,
0099                         const edm::Handle<L1GctEtHadCollection> &totalHtD,
0100                         const edm::Handle<L1GctEtHadCollection> &totalHtE,
0101                         const edm::Handle<L1CaloRegionCollection> &caloRegions);
0102   void plotMissingEt(const edm::Handle<L1GctEtMissCollection> &missingEtD,
0103                      const edm::Handle<L1GctEtMissCollection> &missingEtE);
0104   void plotMissingHt(const edm::Handle<L1GctHtMissCollection> &missingHtD,
0105                      const edm::Handle<L1GctHtMissCollection> &missingHtE);
0106   void plotMissingEErrors(const edm::Handle<L1GctEtMissCollection> &missingEtD,
0107                           const edm::Handle<L1GctEtMissCollection> &missingEtE,
0108                           const edm::Handle<L1GctHtMissCollection> &missingHtD,
0109                           const edm::Handle<L1GctHtMissCollection> &missingHtE,
0110                           edm::Handle<L1CaloRegionCollection> &caloRegions,
0111                           const edm::Handle<L1GctInternJetDataCollection> &intjetsemu,
0112                           const edm::Handle<L1GctInternHtMissCollection> intMissingHtD);
0113   template <class T>
0114   bool checkCollections(const T &collection, const unsigned int &constraint, const std::string &label);
0115 
0116 public:
0117   explicit GctErrorAnalyzer(const edm::ParameterSet &);
0118   ~GctErrorAnalyzer() override;
0119 
0120 private:
0121   void beginJob() override;
0122   void analyze(const edm::Event &, const edm::EventSetup &) override;
0123   void endJob() override;
0124 
0125   // ----------member data ---------------------------
0126   //the following flags select what we'd like to plot and whether or not we want error information
0127   bool doRCT_;
0128   bool doEg_;
0129   bool doIsoDebug_;
0130   bool doNonIsoDebug_;
0131   bool doJets_;
0132   bool doCenJetsDebug_;
0133   bool doTauJetsDebug_;
0134   bool doForJetsDebug_;
0135   bool doHF_;
0136   bool doRingSumDebug_;
0137   bool doBitCountDebug_;
0138   bool doTotalEnergySums_;
0139   bool doTotalEtDebug_;
0140   bool doTotalHtDebug_;
0141   bool doMissingEnergySums_;
0142   bool doMissingETDebug_;
0143   bool doMissingHTDebug_;
0144   bool doExtraMissingHTDebug_;
0145   //the following flags configure whether or not we want multiple BX behaviour for
0146   //1. RCT regions
0147   //2. Emulator output
0148   //3. Hardware output
0149   bool doRCTMBx_;
0150   bool doEmuMBx_;
0151   bool doGCTMBx_;
0152   //the following values select the definition of the "triggered" Bx i.e. where to define Bx=0
0153   int RCTTrigBx_;
0154   int EmuTrigBx_;
0155   int GCTTrigBx_;
0156   //the following flags contain the location of the hardware and emulator digis
0157   edm::InputTag dataTag_;
0158   edm::InputTag emuTag_;
0159   //the following is a string which dictates whether or not we want to use the lab or full system parameters
0160   std::string useSys_;
0161 
0162   //the following declares a struct to hold the MBX Info to make it easy to pass the information around
0163   GctErrorAnalyzerMBxInfo MBxInfo;
0164 
0165   // histograms
0166   //RCT Regions
0167   TH2I *RCT_EtEtaPhi_, *RCT_TvEtaPhi_, *RCT_FgEtaPhi_, *RCT_OfEtaPhi_;
0168   //isoEg
0169   TH1I *isoEgD_Rank_, *isoEgE_Rank_;
0170   TH2I *isoEgD_EtEtaPhi_, *isoEgE_EtEtaPhi_;
0171   TH2I *isoEgD_OccEtaPhi_, *isoEgE_OccEtaPhi_;
0172   TH1I *isoEg_errorFlag_;
0173   //Global Error Histograms
0174   TH1I *isoEgD_GlobalError_Rank_;
0175   TH1I *isoEgE_GlobalError_Rank_;
0176   TH2I *isoEgD_GlobalError_EtEtaPhi_;
0177   TH2I *isoEgE_GlobalError_EtEtaPhi_;
0178   //nonIsoEg
0179   TH1I *nonIsoEgD_Rank_, *nonIsoEgE_Rank_;
0180   TH2I *nonIsoEgD_EtEtaPhi_, *nonIsoEgE_EtEtaPhi_;
0181   TH2I *nonIsoEgD_OccEtaPhi_, *nonIsoEgE_OccEtaPhi_;
0182   TH1I *nonIsoEg_errorFlag_;
0183   //Global Error Histograms
0184   TH1I *nonIsoEgD_GlobalError_Rank_;
0185   TH1I *nonIsoEgE_GlobalError_Rank_;
0186   TH2I *nonIsoEgD_GlobalError_EtEtaPhi_;
0187   TH2I *nonIsoEgE_GlobalError_EtEtaPhi_;
0188   //cenJet
0189   TH1I *cenJetD_Rank_, *cenJetE_Rank_;
0190   TH2I *cenJetD_EtEtaPhi_, *cenJetE_EtEtaPhi_;
0191   TH2I *cenJetD_OccEtaPhi_, *cenJetE_OccEtaPhi_;
0192   TH1I *cenJet_errorFlag_;
0193   //Global Error Histograms
0194   TH1I *cenJetD_GlobalError_Rank_;
0195   TH1I *cenJetE_GlobalError_Rank_;
0196   TH2I *cenJetD_GlobalError_EtEtaPhi_;
0197   TH2I *cenJetE_GlobalError_EtEtaPhi_;
0198   //tauJet
0199   TH1I *tauJetD_Rank_, *tauJetE_Rank_;
0200   TH2I *tauJetD_EtEtaPhi_, *tauJetE_EtEtaPhi_;
0201   TH2I *tauJetD_OccEtaPhi_, *tauJetE_OccEtaPhi_;
0202   TH1I *tauJet_errorFlag_;
0203   //Global Error Histograms
0204   TH1I *tauJetD_GlobalError_Rank_;
0205   TH1I *tauJetE_GlobalError_Rank_;
0206   TH2I *tauJetD_GlobalError_EtEtaPhi_;
0207   TH2I *tauJetE_GlobalError_EtEtaPhi_;
0208   //forJet
0209   TH1I *forJetD_Rank_, *forJetE_Rank_;
0210   TH2I *forJetD_EtEtaPhi_, *forJetE_EtEtaPhi_;
0211   TH2I *forJetD_OccEtaPhi_, *forJetE_OccEtaPhi_;
0212   TH1I *forJet_errorFlag_;
0213   //Global Error Histograms
0214   TH1I *forJetD_GlobalError_Rank_;
0215   TH1I *forJetE_GlobalError_Rank_;
0216   TH2I *forJetD_GlobalError_EtEtaPhi_;
0217   TH2I *forJetE_GlobalError_EtEtaPhi_;
0218   //intJet
0219   TH2I *intJetEtEtaPhiE_;
0220   TH1I *intJetE_Et_;
0221   TH1I *intJetE_Of_;
0222   TH1I *intJetE_Jet1Et_;
0223   TH1I *intJetE_Jet2Et_;
0224   TH1I *intJetE_Jet3Et_;
0225   TH1I *intJetE_Jet4Et_;
0226   //ringSums
0227   TH1I *hfRingSumD_1pos_, *hfRingSumD_1neg_, *hfRingSumD_2pos_, *hfRingSumD_2neg_;
0228   TH1I *hfRingSumE_1pos_, *hfRingSumE_1neg_, *hfRingSumE_2pos_, *hfRingSumE_2neg_;
0229   TH1I *hfRingSum_errorFlag_;
0230   //bitcounts
0231   TH1I *hfBitCountD_1pos_, *hfBitCountD_1neg_, *hfBitCountD_2pos_, *hfBitCountD_2neg_;
0232   TH1I *hfBitCountE_1pos_, *hfBitCountE_1neg_, *hfBitCountE_2pos_, *hfBitCountE_2neg_;
0233   TH1I *hfBitCount_errorFlag_;
0234   //totalEt
0235   TH1I *totalEtD_, *totalEtE_;
0236   TH1I *totalEtD_Of_, *totalEtE_Of_;
0237   TH1I *totalEt_errorFlag_;
0238   //ET GlobalError Histograms
0239   //TH1I *totalEtD_GlobalError_, *totalEtE_GlobalError_;
0240   //TH1I *totalEtD_GlobalError_Of_, *totalEtE_GlobalError_Of_;
0241   //totalHt
0242   TH1I *totalHtD_, *totalHtE_;
0243   TH1I *totalHtD_Of_, *totalHtE_Of_;
0244   TH1I *totalHt_errorFlag_;
0245   //HT GlobalError Histograms
0246   //TH1I *totalHtD_GlobalError_, *totalHtE_GlobalError_;
0247   //TH1I *totalHtD_GlobalError_Of_, *totalHtE_GlobalError_Of_;
0248   //missingET
0249   TH1I *missingEtD_, *missingEtE_;
0250   TH1I *missingEtD_Of_, *missingEtE_Of_;
0251   TH1I *missingEtD_Phi_, *missingEtE_Phi_;
0252   TH1I *missingEt_errorFlag_;
0253   //missingHT
0254   TH1I *missingHtD_, *missingHtE_;
0255   TH1I *missingHtD_Of_, *missingHtE_Of_;
0256   TH1I *missingHtD_Phi_, *missingHtE_Phi_;
0257   TH1I *missingHt_errorFlag_;
0258   TH1I *missingHtD_HtXPosLeaf1, *missingHtD_HtXPosLeaf2, *missingHtD_HtXPosLeaf3, *missingHtD_HtXNegLeaf1,
0259       *missingHtD_HtXNegLeaf2, *missingHtD_HtXNegLeaf3;
0260   TH1I *missingHtD_HtYPosLeaf1, *missingHtD_HtYPosLeaf2, *missingHtD_HtYPosLeaf3, *missingHtD_HtYNegLeaf1,
0261       *missingHtD_HtYNegLeaf2, *missingHtD_HtYNegLeaf3;
0262 
0263   //error flags to decide whether or not to print debug info
0264   bool isIsoError;
0265   bool isNonIsoError;
0266   bool isCenJetError;
0267   bool isTauJetError;
0268   bool isForJetError;
0269   bool isRingSumError;
0270   bool isBitCountError;
0271   bool isTotalEError;
0272   bool isTotalHError;
0273   bool isMissingEError;
0274   bool isMissingHError;
0275 
0276   //Directories - put this here because we want to
0277   //add directories dynamically to this folder
0278   //depending on the errors we find
0279   std::vector<TFileDirectory> errorHistCat;
0280 
0281   //the event number
0282   unsigned int eventNumber;
0283 
0284   const unsigned int *RCT_REGION_QUANTA;
0285 };
0286 
0287 //
0288 // constants, enums and typedefs
0289 // use in conjunction with the templated bits
0290 typedef compareCands<edm::Handle<L1GctEmCandCollection> > compareEG;
0291 typedef compareCands<edm::Handle<L1GctJetCandCollection> > compareJets;
0292 typedef compareTotalEnergySums<edm::Handle<L1GctEtTotalCollection> > compareTotalE;
0293 typedef compareTotalEnergySums<edm::Handle<L1GctEtHadCollection> > compareTotalH;
0294 typedef compareMissingEnergySums<edm::Handle<L1GctEtMissCollection> > compareMissingE;
0295 typedef compareMissingEnergySums<edm::Handle<L1GctHtMissCollection> > compareMissingH;
0296 //
0297 // static data member definitions
0298 //
0299 
0300 //
0301 // constructors and destructor
0302 //
0303 GctErrorAnalyzer::GctErrorAnalyzer(const edm::ParameterSet &iConfig)
0304     : doRCT_(iConfig.getUntrackedParameter<bool>("doRCT", true)),
0305       doEg_(iConfig.getUntrackedParameter<bool>("doEg", true)),
0306       doIsoDebug_(iConfig.getUntrackedParameter<bool>("doIsoDebug", true)),
0307       doNonIsoDebug_(iConfig.getUntrackedParameter<bool>("doNonIsoDebug", true)),
0308       doJets_(iConfig.getUntrackedParameter<bool>("doJets", true)),
0309       doCenJetsDebug_(iConfig.getUntrackedParameter<bool>("doCenJetsDebug", true)),
0310       doTauJetsDebug_(iConfig.getUntrackedParameter<bool>("doTauJetsDebug", true)),
0311       doForJetsDebug_(iConfig.getUntrackedParameter<bool>("doForJetsDebug", true)),
0312       doHF_(iConfig.getUntrackedParameter<bool>("doHF", true)),
0313       doRingSumDebug_(iConfig.getUntrackedParameter<bool>("doRingSumDebug", true)),
0314       doBitCountDebug_(iConfig.getUntrackedParameter<bool>("doBitCountDebug", true)),
0315       doTotalEnergySums_(iConfig.getUntrackedParameter<bool>("doTotalEnergySums", true)),
0316       doTotalEtDebug_(iConfig.getUntrackedParameter<bool>("doTotalEtDebug", true)),
0317       doTotalHtDebug_(iConfig.getUntrackedParameter<bool>("doTotalHtDebug", true)),
0318       doMissingEnergySums_(iConfig.getUntrackedParameter<bool>("doMissingEnergySums", true)),
0319       doMissingETDebug_(iConfig.getUntrackedParameter<bool>("doMissingETDebug", true)),
0320       doMissingHTDebug_(iConfig.getUntrackedParameter<bool>("doMissingHTDebug", true)),
0321       doExtraMissingHTDebug_(iConfig.getUntrackedParameter<bool>("doExtraMissingHTDebug", false)),
0322       doRCTMBx_(iConfig.getUntrackedParameter<bool>("doRCTMBx", false)),
0323       doEmuMBx_(iConfig.getUntrackedParameter<bool>("doEmuMBx", false)),
0324       doGCTMBx_(iConfig.getUntrackedParameter<bool>("doGCTMBx", false)),
0325       RCTTrigBx_(iConfig.getUntrackedParameter<int>("RCTTrigBx", 0)),
0326       EmuTrigBx_(iConfig.getUntrackedParameter<int>("EmuTrigBx", 0)),
0327       GCTTrigBx_(iConfig.getUntrackedParameter<int>("GCTTrigBx", 0)),
0328       dataTag_(iConfig.getUntrackedParameter<edm::InputTag>("dataTag", edm::InputTag("gctDigis"))),
0329       emuTag_(iConfig.getUntrackedParameter<edm::InputTag>("emuTag", edm::InputTag("gctEmuDigis"))),
0330       useSys_(iConfig.getUntrackedParameter<std::string>("useSys", "P5")) {
0331   //now do what ever initialization is needed
0332   //make the root file
0333   usesResource(TFileService::kSharedResource);
0334   edm::Service<TFileService> fs;
0335 
0336   //to try to make this look more elegant
0337   //make a string for each folder we'd like for the Data and Emulator Histograms
0338   std::vector<std::string> quantities;
0339   quantities.push_back("IsoEm");
0340   quantities.push_back("NonIsoEM");
0341   quantities.push_back("CenJets");
0342   quantities.push_back("TauJets");
0343   quantities.push_back("ForJets");
0344   quantities.push_back("HFRingSums");
0345   quantities.push_back("HFBitCounts");
0346   quantities.push_back("TotalESums");
0347   quantities.push_back("MissingESums");
0348 
0349   //make the Emulator Histogram directory
0350   TFileDirectory emuHist = fs->mkdir("EmulatorHistograms");
0351   std::vector<TFileDirectory> emuHistCat;
0352 
0353   //make the Data Histogram directory
0354   TFileDirectory dataHist = fs->mkdir("DataHistograms");
0355   std::vector<TFileDirectory> dataHistCat;
0356 
0357   //make the ErrorFlags directory
0358   TFileDirectory errorHistFlags = fs->mkdir("ErrorHistograms_Flags");
0359 
0360   //make the ErrorDebug directory
0361   TFileDirectory errorHistDetails = fs->mkdir("ErrorHistograms_Details");
0362 
0363   for (unsigned int i = 0; i < quantities.size(); i++) {
0364     //fill the data and emulator folders with the directories
0365     emuHistCat.push_back(emuHist.mkdir(quantities.at(i)));
0366     dataHistCat.push_back(dataHist.mkdir(quantities.at(i)));
0367   }
0368 
0369   //add a folder for RCT Regions - which only exist in data
0370   dataHistCat.push_back(dataHist.mkdir("RCTRegions"));
0371   //and add a folder for the Intermediate Jets - which only exist in emulator
0372   emuHistCat.push_back(emuHist.mkdir("IntJets"));
0373 
0374   //Fill the ErrorDebug folder with the directories
0375   errorHistCat.push_back(errorHistDetails.mkdir("EM"));
0376   errorHistCat.push_back(errorHistDetails.mkdir("Jets"));
0377   errorHistCat.push_back(errorHistDetails.mkdir("HF"));
0378   errorHistCat.push_back(errorHistDetails.mkdir("TotalE"));
0379   errorHistCat.push_back(errorHistDetails.mkdir("MissingE"));
0380 
0381   //BOOK HISTOGRAMS
0382   RCT_EtEtaPhi_ = dataHistCat.at(9).make<TH2I>(
0383       "RCT_EtEtaPhi", "RCT_EtEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
0384   RCT_TvEtaPhi_ = dataHistCat.at(9).make<TH2I>(
0385       "RCT_TvEtaPhi", "RCT_TvEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
0386   RCT_FgEtaPhi_ = dataHistCat.at(9).make<TH2I>(
0387       "RCT_FgEtaPhi", "RCT_FgEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
0388   RCT_OfEtaPhi_ = dataHistCat.at(9).make<TH2I>(
0389       "RCT_OfEtEtaPhi", "RCT_OfEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
0390   //isoEg
0391   isoEgD_Rank_ = dataHistCat.at(0).make<TH1I>("isoEgD_Rank", "isoEgD_Rank;Rank;Number of Events", 64, -0.5, 63.5);
0392   isoEgE_Rank_ = emuHistCat.at(0).make<TH1I>("isoEgE_Rank", "isoEgE_Rank;Rank;Number of Events", 64, -0.5, 63.5);
0393   isoEgD_EtEtaPhi_ = dataHistCat.at(0).make<TH2I>(
0394       "isoEgD_EtEtaPhi", "isoEgD_EtEtaPhi;#eta (GCT Units);#phi(GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
0395   isoEgE_EtEtaPhi_ = emuHistCat.at(0).make<TH2I>(
0396       "isoEgE_EtEtaPhi", "isoEgE_EtEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
0397   isoEgD_OccEtaPhi_ = dataHistCat.at(0).make<TH2I>(
0398       "isoEgD_OccEtaPhi", "isoEgD_OccEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
0399   isoEgE_OccEtaPhi_ = emuHistCat.at(0).make<TH2I>(
0400       "isoEgE_OccEtaPhi", "isoEgE_OccEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
0401   isoEg_errorFlag_ =
0402       errorHistFlags.make<TH1I>("isoEg_errorFlag", "isoEg_errorFlag;Status;Number of Candidates", 3, -0.5, 2.5);
0403   //Global isoEg Error
0404   isoEgD_GlobalError_Rank_ = errorHistCat.at(0).make<TH1I>(
0405       "isoEgD_GlobalError_Rank", "isoEgD_GlobalError_Rank;Rank;Number of Events", 64, -0.5, 63.5);
0406   isoEgE_GlobalError_Rank_ = errorHistCat.at(0).make<TH1I>(
0407       "isoEgE_GlobalError_Rank", "isoEgE_GlobalError_Rank;Rank;Number of Events", 64, -0.5, 63.5);
0408   isoEgD_GlobalError_EtEtaPhi_ = errorHistCat.at(0).make<TH2I>(
0409       "isoEgD_GlobalError_EtEtaPhi", "isoEgD_GlobalError_EtEtaPhi", 22, -0.5, 21.5, 18, -0.5, 17.5);
0410   isoEgE_GlobalError_EtEtaPhi_ = errorHistCat.at(0).make<TH2I>(
0411       "isoEgE_GlobalError_EtEtaPhi", "isoEgE_GlobalError_EtEtaPhi", 22, -0.5, 21.5, 18, -0.5, 17.5);
0412   //nonIsoEg
0413   nonIsoEgD_Rank_ =
0414       dataHistCat.at(1).make<TH1I>("nonIsoEgD_Rank", "nonIsoEgD_Rank;Rank;Number of Events", 64, -0.5, 63.5);
0415   nonIsoEgE_Rank_ =
0416       emuHistCat.at(1).make<TH1I>("nonIsoEgE_Rank", "nonIsoEgE_Rank;Rank;Number of Events", 64, -0.5, 63.5);
0417   nonIsoEgD_EtEtaPhi_ = dataHistCat.at(1).make<TH2I>(
0418       "nonIsoEgD_EtEtaPhi", "nonIsoEgD_EtEtaPhi;#eta (GCT Units);#phi(GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
0419   nonIsoEgE_EtEtaPhi_ = emuHistCat.at(1).make<TH2I>(
0420       "nonIsoEgE_EtEtaPhi", "nonIsoEgE_EtEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
0421   nonIsoEgD_OccEtaPhi_ = dataHistCat.at(1).make<TH2I>(
0422       "nonIsoEgD_OccEtaPhi", "nonIsoEgD_OccEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
0423   nonIsoEgE_OccEtaPhi_ = emuHistCat.at(1).make<TH2I>(
0424       "nonIsoEgE_OccEtaPhi", "nonIsoEgE_OccEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
0425   nonIsoEg_errorFlag_ =
0426       errorHistFlags.make<TH1I>("nonIsoEg_errorFlag", "nonIsoEg_errorFlag;Status;Number of Candidates", 3, -0.5, 2.5);
0427   //Global nonIsoEg Error
0428   nonIsoEgD_GlobalError_Rank_ = errorHistCat.at(0).make<TH1I>(
0429       "nonIsoEgD_GlobalError_Rank", "nonIsoEgD_GlobalError_Rank;Rank;Number of Events", 64, -0.5, 63.5);
0430   nonIsoEgE_GlobalError_Rank_ = errorHistCat.at(0).make<TH1I>(
0431       "nonIsoEgE_GlobalError_Rank", "nonIsoEgE_GlobalError_Rank;Rank;Number of Events", 64, -0.5, 63.5);
0432   nonIsoEgD_GlobalError_EtEtaPhi_ = errorHistCat.at(0).make<TH2I>(
0433       "nonIsoEgD_GlobalError_EtEtaPhi", "nonIsoEgD_GlobalError_EtEtaPhi", 22, -0.5, 21.5, 18, -0.5, 17.5);
0434   nonIsoEgE_GlobalError_EtEtaPhi_ = errorHistCat.at(0).make<TH2I>(
0435       "nonIsoEgE_GlobalError_EtEtaPhi", "nonIsoEgE_GlobalError_EtEtaPhi", 22, -0.5, 21.5, 18, -0.5, 17.5);
0436   //CenJets
0437   cenJetD_Rank_ = dataHistCat.at(2).make<TH1I>("cenJetD_Rank", "cenJetD_Rank;Rank;Number of Events", 64, -0.5, 63.5);
0438   cenJetE_Rank_ = emuHistCat.at(2).make<TH1I>("cenJetE_Rank", "cenJetE_Rank;Rank;Number of Events", 64, -0.5, 63.5);
0439   cenJetD_EtEtaPhi_ = dataHistCat.at(2).make<TH2I>(
0440       "cenJetD_EtEtaPhi", "cenJetD_EtEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
0441   cenJetE_EtEtaPhi_ = emuHistCat.at(2).make<TH2I>(
0442       "cenJetE_EtEtaPhi", "cenJetE_EtEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
0443   cenJetD_OccEtaPhi_ = dataHistCat.at(2).make<TH2I>(
0444       "cenJetD_OccEtaPhi", "cenJetD_OccEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
0445   cenJetE_OccEtaPhi_ = emuHistCat.at(2).make<TH2I>(
0446       "cenJetE_OccEtaPhi", "cenJetE_OccEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
0447   cenJet_errorFlag_ =
0448       errorHistFlags.make<TH1I>("cenJet_errorFlag", "cenJet_errorFlag;Status;Number of Candidates", 3, -0.5, 2.5);
0449   //Global CenJet Error
0450   cenJetD_GlobalError_Rank_ =
0451       errorHistCat.at(1).make<TH1I>("cenJetD_GlobalError_Rank", "cenJetD_GlobalError_Rank", 64, -0.5, 63.5);
0452   cenJetE_GlobalError_Rank_ =
0453       errorHistCat.at(1).make<TH1I>("cenJetE_GlobalError_Rank", "cenJetE_GlobalError_Rank", 64, -0.5, 63.5);
0454   cenJetD_GlobalError_EtEtaPhi_ = errorHistCat.at(1).make<TH2I>(
0455       "cenJetD_GlobalError_EtEtaPhi", "cenJetD_GlobalError_EtEtaPhi", 22, -0.5, 21.5, 18, -0.5, 17.5);
0456   cenJetE_GlobalError_EtEtaPhi_ = errorHistCat.at(1).make<TH2I>(
0457       "cenJetE_GlobalError_EtEtaPhi", "cenJetE_GlobalError_EtEtaPhi", 22, -0.5, 21.5, 18, -0.5, 17.5);
0458   //TauJets
0459   tauJetD_Rank_ = dataHistCat.at(3).make<TH1I>("tauJetD_Rank", "tauJetD_Rank;Rank;Number of Events", 64, -0.5, 63.5);
0460   tauJetE_Rank_ = emuHistCat.at(3).make<TH1I>("tauJetE_Rank", "tauJetE_Rank;Rank;Number of Events", 64, -0.5, 63.5);
0461   tauJetD_EtEtaPhi_ = dataHistCat.at(3).make<TH2I>(
0462       "tauJetD_EtEtaPhi", "tauJetD_EtEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
0463   tauJetE_EtEtaPhi_ = emuHistCat.at(3).make<TH2I>(
0464       "tauJetE_EtEtaPhi", "tauJetE_EtEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
0465   tauJetD_OccEtaPhi_ = dataHistCat.at(3).make<TH2I>(
0466       "tauJetD_OccEtaPhi", "tauJetD_OccEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
0467   tauJetE_OccEtaPhi_ = emuHistCat.at(3).make<TH2I>(
0468       "tauJetE_OccEtaPhi", "tauJetE_OccEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
0469   tauJet_errorFlag_ =
0470       errorHistFlags.make<TH1I>("tauJet_errorFlag", "tauJet_errorFlag;Status;Number of Candidates", 3, -0.5, 2.5);
0471   //Global TauJet Error
0472   tauJetD_GlobalError_Rank_ =
0473       errorHistCat.at(1).make<TH1I>("tauJetD_GlobalError_Rank", "tauJetD_GlobalError_Rank", 64, -0.5, 63.5);
0474   tauJetE_GlobalError_Rank_ =
0475       errorHistCat.at(1).make<TH1I>("tauJetE_GlobalError_Rank", "tauJetE_GlobalError_Rank", 64, -0.5, 63.5);
0476   tauJetD_GlobalError_EtEtaPhi_ = errorHistCat.at(1).make<TH2I>(
0477       "tauJetD_GlobalError_EtEtaPhi", "tauJetD_GlobalError_EtEtaPhi", 22, -0.5, 21.5, 18, -0.5, 17.5);
0478   tauJetE_GlobalError_EtEtaPhi_ = errorHistCat.at(1).make<TH2I>(
0479       "tauJetE_GlobalError_EtEtaPhi", "tauJetE_GlobalError_EtEtaPhi", 22, -0.5, 21.5, 18, -0.5, 17.5);
0480   //ForJets
0481   forJetD_Rank_ = dataHistCat.at(4).make<TH1I>("forJetD_Rank", "forJetD_Rank;Rank;Number of Events", 64, -0.5, 63.5);
0482   forJetE_Rank_ = emuHistCat.at(4).make<TH1I>("forJetE_Rank", "forJetE_Rank;Rank;Number of Events", 64, -0.5, 63.5);
0483   forJetD_EtEtaPhi_ = dataHistCat.at(4).make<TH2I>(
0484       "forJetD_EtEtaPhi", "forJetD_EtEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
0485   forJetE_EtEtaPhi_ = emuHistCat.at(4).make<TH2I>(
0486       "forJetE_EtEtaPhi", "forJetE_EtEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
0487   forJetD_OccEtaPhi_ = dataHistCat.at(4).make<TH2I>(
0488       "forJetD_OccEtaPhi", "forJetD_OccEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
0489   forJetE_OccEtaPhi_ = emuHistCat.at(4).make<TH2I>(
0490       "forJetE_OccEtaPhi", "forJetE_OccEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
0491   forJet_errorFlag_ =
0492       errorHistFlags.make<TH1I>("forJet_errorFlag", "forJet_errorFlag;Status;Number of Candidates", 3, -0.5, 2.5);
0493   //Global ForJet Error
0494   forJetD_GlobalError_Rank_ =
0495       errorHistCat.at(1).make<TH1I>("forJetD_GlobalError_Rank", "forJetD_GlobalError_Rank", 64, -0.5, 63.5);
0496   forJetE_GlobalError_Rank_ =
0497       errorHistCat.at(1).make<TH1I>("forJetE_GlobalError_Rank", "forJetE_GlobalError_Rank", 64, -0.5, 63.5);
0498   forJetD_GlobalError_EtEtaPhi_ = errorHistCat.at(1).make<TH2I>(
0499       "forJetD_GlobalError_EtEtaPhi", "forJetD_GlobalError_EtEtaPhi", 22, -0.5, 21.5, 18, -0.5, 17.5);
0500   forJetE_GlobalError_EtEtaPhi_ = errorHistCat.at(1).make<TH2I>(
0501       "forJetE_GlobalError_EtEtaPhi", "forJetE_GlobalError_EtEtaPhi", 22, -0.5, 21.5, 18, -0.5, 17.5);
0502   //IntJets
0503   intJetEtEtaPhiE_ = emuHistCat.at(9).make<TH2I>(
0504       "intJetEtEtaPhiE_", "intJetEtEtaPhiE_;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
0505   intJetE_Et_ = emuHistCat.at(9).make<TH1I>("intJetE_Et", "intJetE_Et;E_{T};Number of Events", 1024, -0.5, 1023.5);
0506   intJetE_Of_ =
0507       emuHistCat.at(9).make<TH1I>("intJetE_Of", "intJetE_Of;Overflow Bit Status;Number of Events", 2, -0.5, 1.5);
0508   intJetE_Jet1Et_ =
0509       emuHistCat.at(9).make<TH1I>("intJetE_Jet1Et", "intJetE_Jet1Et;E_{T};Number of Events", 1024, -0.5, 1023.5);
0510   intJetE_Jet2Et_ =
0511       emuHistCat.at(9).make<TH1I>("intJetE_Jet2Et", "intJetE_Jet2Et;E_{T};Number of Events", 1024, -0.5, 1023.5);
0512   intJetE_Jet3Et_ =
0513       emuHistCat.at(9).make<TH1I>("intJetE_Jet3Et", "intJetE_Jet3Et;E_{T};Number of Events", 1024, -0.5, 1023.5);
0514   intJetE_Jet4Et_ =
0515       emuHistCat.at(9).make<TH1I>("intJetE_Jet4Et", "intJetE_Jet4Et;E_{T};Number of Events", 1024, -0.5, 1023.5);
0516   //HFRing Sums
0517   hfRingSumD_1pos_ = dataHistCat.at(5).make<TH1I>("hfRingSumD_1+", "hfRingSumD_1+;Rank;Number of Events", 8, -0.5, 7.5);
0518   hfRingSumD_1neg_ = dataHistCat.at(5).make<TH1I>("hfRingSumD_1-", "hfRingSumD_1-;Rank;Number of Events", 8, -0.5, 7.5);
0519   hfRingSumD_2pos_ = dataHistCat.at(5).make<TH1I>("hfRingSumD_2+", "hfRingSumD_2+;Rank;Number of Events", 8, -0.5, 7.5);
0520   hfRingSumD_2neg_ = dataHistCat.at(5).make<TH1I>("hfRingSumD_2-", "hfRingSumD_2-;Rank;Number of Events", 8, -0.5, 7.5);
0521   hfRingSumE_1pos_ = emuHistCat.at(5).make<TH1I>("hfRingSumE_1+", "hfRingSumE_1+;Rank;Number of Events", 8, -0.5, 7.5);
0522   hfRingSumE_1neg_ = emuHistCat.at(5).make<TH1I>("hfRingSumE_1-", "hfRingSumE_1-;Rank;Number of Events", 8, -0.5, 7.5);
0523   hfRingSumE_2pos_ = emuHistCat.at(5).make<TH1I>("hfRingSumE_2+", "hfRingSumE_2+;Rank;Number of Events", 8, -0.5, 7.5);
0524   hfRingSumE_2neg_ = emuHistCat.at(5).make<TH1I>("hfRingSumE_2-", "hfRingSumE_2-;Rank;Number of Events", 8, -0.5, 7.5);
0525   hfRingSum_errorFlag_ =
0526       errorHistFlags.make<TH1I>("hfRingSum_errorFlag", "hfRingSum_errorFlag;Status;Number of Candidates", 2, -0.5, 1.5);
0527   //HFRing BitCounts
0528   hfBitCountD_1pos_ =
0529       dataHistCat.at(6).make<TH1I>("hfBitCountD_1+", "hfBitCountD_1+;Rank;Number of Events", 8, -0.5, 7.5);
0530   hfBitCountD_1neg_ =
0531       dataHistCat.at(6).make<TH1I>("hfBitCountD_1-", "hfBitCountD_1-;Rank;Number of Events", 8, -0.5, 7.5);
0532   hfBitCountD_2pos_ =
0533       dataHistCat.at(6).make<TH1I>("hfBitCountD_2+", "hfBitCountD_2+;Rank;Number of Events", 8, -0.5, 7.5);
0534   hfBitCountD_2neg_ =
0535       dataHistCat.at(6).make<TH1I>("hfBitCountD_2-", "hfBitCountD_2-;Rank;Number of Events", 8, -0.5, 7.5);
0536   hfBitCountE_1pos_ =
0537       emuHistCat.at(6).make<TH1I>("hfBitCountE_1+", "hfBitCountE_1+;Rank;Number of Events", 8, -0.5, 7.5);
0538   hfBitCountE_1neg_ =
0539       emuHistCat.at(6).make<TH1I>("hfBitCountE_1-", "hfBitCountE_1-;Rank;Number of Events", 8, -0.5, 7.5);
0540   hfBitCountE_2pos_ =
0541       emuHistCat.at(6).make<TH1I>("hfBitCountE_2+", "hfBitCountE_2+;Rank;Number of Events", 8, -0.5, 7.5);
0542   hfBitCountE_2neg_ =
0543       emuHistCat.at(6).make<TH1I>("hfBitCountE_2-", "hfBitCountE_2-;Rank;Number of Events", 8, -0.5, 7.5);
0544   hfBitCount_errorFlag_ = errorHistFlags.make<TH1I>(
0545       "hfBitCount_errorFlag", "hfBitCount_errorFlag;Status;Number of Candidates", 2, -0.5, 1.5);
0546   //Total ET
0547   totalEtD_ = dataHistCat.at(7).make<TH1I>("totalEtD", "totalEtD;E_{T};Number of Events", 2048, -0.5, 2047.5);
0548   totalEtD_Of_ =
0549       dataHistCat.at(7).make<TH1I>("totalEtD_Of", "totalEtD_Of;Overflow Bit Status;Number of Events", 2, -0.5, 1.5);
0550   totalEtE_ = emuHistCat.at(7).make<TH1I>("totalEtE", "totalEtE;E_{T};Number of Events", 2048, -0.5, 2047.5);
0551   totalEtE_Of_ =
0552       emuHistCat.at(7).make<TH1I>("totalEtE_Of", "totalEtE_Of;Overflow Bit Status;Number of Events", 2, -0.5, 1.5);
0553   totalEt_errorFlag_ =
0554       errorHistFlags.make<TH1I>("totalEt_errorFlag", "totalEt_errorFlag;Status;Number of Candidates", 2, -0.5, 1.5);
0555   //Book the Global ET Error histograms in the errorHistCat
0556   //totalEtD_GlobalError_ = errorHistCat.at(3).make<TH1I>("totalEtD_GlobalError", "totalEtD_GlobalError;E_{T};Number of Events", 1024, -0.5, 1023.5);
0557   //totalEtE_GlobalError_ = errorHistCat.at(3).make<TH1I>("totalEtE_GlobalError", "totalEtE_GlobalError;E_{T};Number of Events", 1024, -0.5, 1023.5);
0558   //totalEtD_GlobalError_Of_ = errorHistCat.at(3).make<TH1I>("totalEtD_GlobalError_Of", "totalEtD_GlobalError_Of;Overflow Bit Status;Number of Events", 2, -0.5, 1.5);
0559   //totalEtE_GlobalError_Of_ = errorHistCat.at(3).make<TH1I>("totalEtE_GlobalError_Of", "totalEtE_GlobalError_Of;Overflow Bit Status;Number of Events", 2, -0.5, 1.5);
0560   //Total HT
0561   totalHtD_ = dataHistCat.at(7).make<TH1I>("totalHtD", "totalHtD;H_{T};Number of Events", 2048, -0.5, 2047.5);
0562   totalHtD_Of_ =
0563       dataHistCat.at(7).make<TH1I>("totalHtD_Of", "totalHtD_Of;Overflow Bit Status;Number of Events", 2, -0.5, 1.5);
0564   totalHtE_ = emuHistCat.at(7).make<TH1I>("totalHtE", "totalHtE;H_{T};Number of Events", 2048, -0.5, 2047.5);
0565   totalHtE_Of_ =
0566       emuHistCat.at(7).make<TH1I>("totalHtE_Of", "totalHtE_Of;Overflow Bit Status;Number of Events", 2, -0.5, 1.5);
0567   totalHt_errorFlag_ =
0568       errorHistFlags.make<TH1I>("totalHt_errorFlag", "totalHt_errorFlag;Status;Number of Candidates", 2, -0.5, 1.5);
0569   //Book the Global HT Error histograms in the errorHistCat
0570   //totalHtD_GlobalError_ = errorHistCat.at(3).make<TH1I>("totalHtD_GlobalError", "totalHtD_GlobalError;E_{T};Number of Events", 1024, -0.5, 1023.5);
0571   //totalHtE_GlobalError_ = errorHistCat.at(3).make<TH1I>("totalHtE_GlobalError", "totalHtE_GlobalError;E_{T};Number of Events", 1024, -0.5, 1023.5);
0572   //totalHtD_GlobalError_Of_ = errorHistCat.at(3).make<TH1I>("totalHtD_GlobalError_Of", "totalHtD_GlobalError_Of;Overflow Bit Status;Number of Events", 2, -0.5, 1.5);
0573   //totalHtE_GlobalError_Of_ = errorHistCat.at(3).make<TH1I>("totalHtE_GlobalError_Of", "totalHtE_GlobalError_Of;Overflow Bit Status;Number of Events", 2, -0.5, 1.5);
0574   //MissingEt
0575   missingEtD_ =
0576       dataHistCat.at(8).make<TH1I>("missingEtD", "missingEtD;Missing E_{T};Number of Events", 1024, -0.5, 1023.5);
0577   missingEtD_Of_ =
0578       dataHistCat.at(8).make<TH1I>("missingEtD_Of", "missingEtD_Of;Overflow Bit Status;Number of Events", 2, -0.5, 1.5);
0579   missingEtD_Phi_ = dataHistCat.at(8).make<TH1I>(
0580       "missingEtD_Phi", "missingEtD_Phi;Missing E_{T} #phi;Number of Events", 72, -0.5, 71.5);
0581   missingEtE_ =
0582       emuHistCat.at(8).make<TH1I>("missingEtE", "missingEtE;Missing E_{T};Number of Events", 1024, -0.5, 1023.5);
0583   missingEtE_Of_ =
0584       emuHistCat.at(8).make<TH1I>("missingEtE_Of", "missingEtE_Of;Overflow Bit Status;Number of Events", 2, -0.5, 1.5);
0585   missingEtE_Phi_ = emuHistCat.at(8).make<TH1I>(
0586       "missingEtE_Phi", "missingEtE_Phi;Missing E_{T} #phi;Number of Events", 72, -0.5, 71.5);
0587   missingEt_errorFlag_ =
0588       errorHistFlags.make<TH1I>("missingEt_errorFlag", "missingEt_errorFlag;Status;Number of Candidates", 4, -0.5, 3.5);
0589   //MissingHt
0590   missingHtD_ =
0591       dataHistCat.at(8).make<TH1I>("missingHtD", "missingHtD;Missing H_{T};Number of Events", 1024, -0.5, 1023.5);
0592   missingHtD_Of_ =
0593       dataHistCat.at(8).make<TH1I>("missingHtD_Of", "missingHtD_Of;Overflow Bit Status;Number of Events", 2, -0.5, 1.5);
0594   missingHtD_Phi_ = dataHistCat.at(8).make<TH1I>(
0595       "missingHtD_Phi", "missingHtD_Phi;Missing H_{T} #phi;Number of Events", 72, -0.5, 71.5);
0596   missingHtE_ =
0597       emuHistCat.at(8).make<TH1I>("missingHtE", "missingHtE;Missing H_{T};Number of Events", 1024, -0.5, 1023.5);
0598   missingHtE_Of_ =
0599       emuHistCat.at(8).make<TH1I>("missingHtE_Of", "missingHtE_Of;Overflow Bit Status;Number of Events", 2, -0.5, 1.5);
0600   missingHtE_Phi_ = emuHistCat.at(8).make<TH1I>(
0601       "missingHtE_Phi", "missingHtE_Phi;Missing H_{T} #phi;Number of Events", 72, -0.5, 71.5);
0602   missingHt_errorFlag_ =
0603       errorHistFlags.make<TH1I>("missingHt_errorFlag", "missingHt_errorFlag;Status;Number of Candidates", 4, -0.5, 3.5);
0604   //Additional MissingHt Debug histograms
0605   missingHtD_HtXPosLeaf1 = dataHistCat.at(8).make<TH1I>(
0606       "missingHtD_HtXPosLeaf1", "missingHtD;Missing H_{T} X PosLeaf1;Number of Events", 4096, -2048.5, 2047.5);
0607   missingHtD_HtXPosLeaf2 = dataHistCat.at(8).make<TH1I>(
0608       "missingHtD_HtXPosLeaf2", "missingHtD;Missing H_{T} X PosLeaf2;Number of Events", 4096, -2048.5, 2047.5);
0609   missingHtD_HtXPosLeaf3 = dataHistCat.at(8).make<TH1I>(
0610       "missingHtD_HtXPosLeaf3", "missingHtD;Missing H_{T} X PosLeaf3;Number of Events", 4096, -2048.5, 2047.5);
0611   missingHtD_HtXNegLeaf1 = dataHistCat.at(8).make<TH1I>(
0612       "missingHtD_HtXNegLeaf1", "missingHtD;Missing H_{T} X NegLeaf1;Number of Events", 4096, -2048.5, 2047.5);
0613   missingHtD_HtXNegLeaf2 = dataHistCat.at(8).make<TH1I>(
0614       "missingHtD_HtXNegLeaf2", "missingHtD;Missing H_{T} X NegLeaf2;Number of Events", 4096, -2048.5, 2047.5);
0615   missingHtD_HtXNegLeaf3 = dataHistCat.at(8).make<TH1I>(
0616       "missingHtD_HtXNegLeaf3", "missingHtD;Missing H_{T} X NegLeaf3;Number of Events", 4096, -2048.5, 2047.5);
0617 
0618   missingHtD_HtYPosLeaf1 = dataHistCat.at(8).make<TH1I>(
0619       "missingHtD_HtYPosLeaf1", "missingHtD;Missing H_{T} Y PosLeaf1;Number of Events", 4096, -2048.5, 2047.5);
0620   missingHtD_HtYPosLeaf2 = dataHistCat.at(8).make<TH1I>(
0621       "missingHtD_HtYPosLeaf2", "missingHtD;Missing H_{T} Y PosLeaf2;Number of Events", 4096, -2048.5, 2047.5);
0622   missingHtD_HtYPosLeaf3 = dataHistCat.at(8).make<TH1I>(
0623       "missingHtD_HtYPosLeaf3", "missingHtD;Missing H_{T} Y PosLeaf3;Number of Events", 4096, -2048.5, 2047.5);
0624   missingHtD_HtYNegLeaf1 = dataHistCat.at(8).make<TH1I>(
0625       "missingHtD_HtYNegLeaf1", "missingHtD;Missing H_{T} Y NegLeaf1;Number of Events", 4096, -2048.5, 2047.5);
0626   missingHtD_HtYNegLeaf2 = dataHistCat.at(8).make<TH1I>(
0627       "missingHtD_HtYNegLeaf2", "missingHtD;Missing H_{T} Y NegLeaf2;Number of Events", 4096, -2048.5, 2047.5);
0628   missingHtD_HtYNegLeaf3 = dataHistCat.at(8).make<TH1I>(
0629       "missingHtD_HtYNegLeaf3", "missingHtD;Missing H_{T} Y NegLeaf3;Number of Events", 4096, -2048.5, 2047.5);
0630 
0631   //Annotate the labels of the error flags
0632   //For the electrons and jets
0633   std::vector<std::string> errorFlagLabels;
0634   errorFlagLabels.push_back("Matched");
0635   errorFlagLabels.push_back("Unmatched Data Cand");
0636   errorFlagLabels.push_back("Unmatched Emul Cand");
0637 
0638   for (unsigned int i = 0; i < errorFlagLabels.size(); i++) {
0639     isoEg_errorFlag_->GetXaxis()->SetBinLabel(i + 1, errorFlagLabels.at(i).c_str());
0640     nonIsoEg_errorFlag_->GetXaxis()->SetBinLabel(i + 1, errorFlagLabels.at(i).c_str());
0641     cenJet_errorFlag_->GetXaxis()->SetBinLabel(i + 1, errorFlagLabels.at(i).c_str());
0642     tauJet_errorFlag_->GetXaxis()->SetBinLabel(i + 1, errorFlagLabels.at(i).c_str());
0643     forJet_errorFlag_->GetXaxis()->SetBinLabel(i + 1, errorFlagLabels.at(i).c_str());
0644   }
0645   errorFlagLabels.clear();
0646 
0647   //For the Total Energy Sums and HF
0648   errorFlagLabels.push_back("Matched");
0649   errorFlagLabels.push_back("Unmatched");
0650 
0651   for (unsigned int i = 0; i < errorFlagLabels.size(); i++) {
0652     hfRingSum_errorFlag_->GetXaxis()->SetBinLabel(i + 1, errorFlagLabels.at(i).c_str());
0653     hfBitCount_errorFlag_->GetXaxis()->SetBinLabel(i + 1, errorFlagLabels.at(i).c_str());
0654     totalEt_errorFlag_->GetXaxis()->SetBinLabel(i + 1, errorFlagLabels.at(i).c_str());
0655     totalHt_errorFlag_->GetXaxis()->SetBinLabel(i + 1, errorFlagLabels.at(i).c_str());
0656   }
0657   errorFlagLabels.clear();
0658 
0659   //For the Missing Energy Sums
0660   errorFlagLabels.push_back("Matched");
0661   errorFlagLabels.push_back("Matched Mag");
0662   errorFlagLabels.push_back("Matched Phi");
0663   errorFlagLabels.push_back("Unmatched");
0664 
0665   for (unsigned int i = 0; i < errorFlagLabels.size(); i++) {
0666     missingEt_errorFlag_->GetXaxis()->SetBinLabel(i + 1, errorFlagLabels.at(i).c_str());
0667     missingHt_errorFlag_->GetXaxis()->SetBinLabel(i + 1, errorFlagLabels.at(i).c_str());
0668   }
0669 
0670   //initialise - set all flags to false as they will be set on an event-by-event basis
0671   isIsoError = false;
0672   isNonIsoError = false;
0673   isCenJetError = false;
0674   isTauJetError = false;
0675   isForJetError = false;
0676   isRingSumError = false;
0677   isBitCountError = false;
0678   isTotalEError = false;
0679   isTotalHError = false;
0680   isMissingEError = false;
0681   isMissingHError = false;
0682 
0683   //fill the struct of MBXinformation. It is easier to pass this information to the respective functions as used below this way
0684   MBxInfo.RCTTrigBx = RCTTrigBx_;
0685   MBxInfo.EmuTrigBx = EmuTrigBx_;
0686   MBxInfo.GCTTrigBx = GCTTrigBx_;
0687 
0688   //set the parameters according to the system chosen
0689   if (useSys_ == "P5") {
0690     RCT_REGION_QUANTA = &RCT_REGION_QUANTA_P5;
0691   } else if (useSys_ == "Lab") {
0692     RCT_REGION_QUANTA = &RCT_REGION_QUANTA_LAB;
0693   } else {
0694     edm::LogWarning("ChosenSystem") << " "
0695                                     << "The system you chose to use (" << useSys_
0696                                     << ") was not recognised. Defaulting to the full system geometry";
0697     RCT_REGION_QUANTA = &RCT_REGION_QUANTA_P5;
0698   }
0699 }
0700 
0701 GctErrorAnalyzer::~GctErrorAnalyzer() {
0702   // do anything here that needs to be done at desctruction time
0703   // (e.g. close files, deallocate resources etc.)
0704 }
0705 
0706 //
0707 // member functions
0708 //
0709 
0710 // ------------ method called to for each event  ------------
0711 void GctErrorAnalyzer::analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) {
0712   using namespace edm;
0713   using namespace std;
0714 
0715   Handle<L1CaloRegionCollection> caloRegions;
0716   Handle<L1CaloEmCollection> emRegions;
0717 
0718   Handle<L1GctEmCandCollection> nonIsoEgD;
0719   Handle<L1GctEmCandCollection> nonIsoEgE;
0720   Handle<L1GctEmCandCollection> isoEgD;
0721   Handle<L1GctEmCandCollection> isoEgE;
0722 
0723   Handle<L1GctJetCandCollection> cenJetsD;
0724   Handle<L1GctJetCandCollection> cenJetsE;
0725   Handle<L1GctJetCandCollection> forJetsD;
0726   Handle<L1GctJetCandCollection> forJetsE;
0727   Handle<L1GctJetCandCollection> tauJetsD;
0728   Handle<L1GctJetCandCollection> tauJetsE;
0729 
0730   Handle<L1GctInternJetDataCollection> intJetsE;
0731 
0732   Handle<L1GctHFRingEtSumsCollection> hfRingSumsD;
0733   Handle<L1GctHFRingEtSumsCollection> hfRingSumsE;
0734 
0735   Handle<L1GctHFBitCountsCollection> hfBitCountsD;
0736   Handle<L1GctHFBitCountsCollection> hfBitCountsE;
0737 
0738   Handle<L1GctEtTotalCollection> totalEtD;
0739   Handle<L1GctEtTotalCollection> totalEtE;
0740 
0741   Handle<L1GctEtHadCollection> totalHtD;
0742   Handle<L1GctEtHadCollection> totalHtE;
0743 
0744   Handle<L1GctEtMissCollection> missingEtD;
0745   Handle<L1GctEtMissCollection> missingEtE;
0746 
0747   Handle<L1GctHtMissCollection> missingHtD;
0748   Handle<L1GctHtMissCollection> missingHtE;
0749 
0750   Handle<L1GctInternHtMissCollection> intHtMissD;
0751 
0752   //we need this for all user cases...
0753   iEvent.getByLabel(dataTag_.label(), caloRegions);
0754 
0755   //in order to allow the debug folders to have a unique name (so that when jobs are split in crab, we can merge)
0756   //use the eventnum in the folder name
0757   eventNumber = iEvent.id().event();
0758 
0759   if (doRCT_) {
0760     if (checkCollections(caloRegions, *RCT_REGION_QUANTA, "RCT CaloRegions"))
0761       plotRCTRegions(caloRegions);
0762   }
0763 
0764   if (doEg_) {
0765     iEvent.getByLabel(dataTag_.label(), "nonIsoEm", nonIsoEgD);
0766     iEvent.getByLabel(emuTag_.label(), "nonIsoEm", nonIsoEgE);
0767 
0768     iEvent.getByLabel(dataTag_.label(), "isoEm", isoEgD);
0769     iEvent.getByLabel(emuTag_.label(), "isoEm", isoEgE);
0770 
0771     isIsoError = false;
0772     isNonIsoError = false;
0773 
0774     if (checkCollections(isoEgD, GCT_OBJECT_QUANTA, "Iso e/g Data") &&
0775         checkCollections(isoEgE, GCT_OBJECT_QUANTA, "Iso e/g Emulator")) {
0776       plotIsoEm(isoEgD, isoEgE);
0777       compareEG isoCompare(isoEgD, isoEgE, MBxInfo);
0778       isIsoError = isoCompare.doCompare(isoEg_errorFlag_,
0779                                         isoEgD_GlobalError_Rank_,
0780                                         isoEgD_GlobalError_EtEtaPhi_,
0781                                         isoEgE_GlobalError_Rank_,
0782                                         isoEgE_GlobalError_EtEtaPhi_);
0783     }
0784 
0785     if (checkCollections(nonIsoEgD, GCT_OBJECT_QUANTA, "NonIso e/g Data") &&
0786         checkCollections(nonIsoEgE, GCT_OBJECT_QUANTA, "NonIso e/g Emulator")) {
0787       plotNonIsoEm(nonIsoEgD, nonIsoEgE);
0788       compareEG nonIsoCompare(nonIsoEgD, nonIsoEgE, MBxInfo);
0789       isNonIsoError = nonIsoCompare.doCompare(nonIsoEg_errorFlag_,
0790                                               nonIsoEgD_GlobalError_Rank_,
0791                                               nonIsoEgD_GlobalError_EtEtaPhi_,
0792                                               nonIsoEgE_GlobalError_Rank_,
0793                                               nonIsoEgE_GlobalError_EtEtaPhi_);
0794     }
0795 
0796     if ((isIsoError && doIsoDebug_) || (isNonIsoError && doNonIsoDebug_)) {
0797       iEvent.getByLabel(dataTag_.label(), emRegions);
0798       if (checkCollections(emRegions, RCT_EM_OBJECT_QUANTA, "RCT EMRegions"))
0799         plotEGErrors(isoEgD, isoEgE, nonIsoEgD, nonIsoEgE, emRegions);
0800     }
0801   }
0802 
0803   if (doJets_) {
0804     iEvent.getByLabel(emuTag_.label(), "cenJets", cenJetsE);
0805     iEvent.getByLabel(dataTag_.label(), "cenJets", cenJetsD);
0806 
0807     iEvent.getByLabel(emuTag_.label(), "forJets", forJetsE);
0808     iEvent.getByLabel(dataTag_.label(), "forJets", forJetsD);
0809 
0810     iEvent.getByLabel(emuTag_.label(), "tauJets", tauJetsE);
0811     iEvent.getByLabel(dataTag_.label(), "tauJets", tauJetsD);
0812 
0813     iEvent.getByLabel(emuTag_.label(), intJetsE);
0814 
0815     isCenJetError = false;
0816     isTauJetError = false;
0817     isForJetError = false;
0818 
0819     //Central Jets
0820     if (checkCollections(cenJetsD, GCT_OBJECT_QUANTA, "Central Jets Data") &&
0821         checkCollections(cenJetsE, GCT_OBJECT_QUANTA, "Central Jets Emulator")) {
0822       plotCenJets(cenJetsD, cenJetsE);
0823       compareJets cenJetsCompare(cenJetsD, cenJetsE, MBxInfo);
0824       isCenJetError = cenJetsCompare.doCompare(cenJet_errorFlag_,
0825                                                cenJetD_GlobalError_Rank_,
0826                                                cenJetD_GlobalError_EtEtaPhi_,
0827                                                cenJetE_GlobalError_Rank_,
0828                                                cenJetE_GlobalError_EtEtaPhi_);
0829     }
0830 
0831     //Tau Jets
0832     if (checkCollections(tauJetsD, GCT_OBJECT_QUANTA, "Tau Jets Data") &&
0833         checkCollections(tauJetsE, GCT_OBJECT_QUANTA, "Tau Jets Emulator")) {
0834       plotTauJets(tauJetsD, tauJetsE);
0835       compareJets tauJetsCompare(tauJetsD, tauJetsE, MBxInfo);
0836       isTauJetError = tauJetsCompare.doCompare(tauJet_errorFlag_,
0837                                                tauJetD_GlobalError_Rank_,
0838                                                tauJetD_GlobalError_EtEtaPhi_,
0839                                                tauJetE_GlobalError_Rank_,
0840                                                tauJetE_GlobalError_EtEtaPhi_);
0841     }
0842 
0843     //For Jets
0844     if (checkCollections(forJetsD, GCT_OBJECT_QUANTA, "Forward Jets Data") &&
0845         checkCollections(forJetsE, GCT_OBJECT_QUANTA, "Forward Jets Emulator")) {
0846       plotForJets(forJetsD, forJetsE);
0847       compareJets forJetsCompare(forJetsD, forJetsE, MBxInfo);
0848       isForJetError = forJetsCompare.doCompare(forJet_errorFlag_,
0849                                                forJetD_GlobalError_Rank_,
0850                                                forJetD_GlobalError_EtEtaPhi_,
0851                                                forJetE_GlobalError_Rank_,
0852                                                forJetE_GlobalError_EtEtaPhi_);
0853     }
0854 
0855     //Emulator Intermediate Jets
0856     if (checkCollections(intJetsE, NUM_INT_JETS, "Intermediate Jets Emulator"))
0857       plotIntJets(intJetsE);
0858 
0859     if ((isCenJetError && doCenJetsDebug_) || (isTauJetError && doTauJetsDebug_) ||
0860         (isForJetError && doForJetsDebug_)) {
0861       plotJetErrors(cenJetsD, cenJetsE, tauJetsD, tauJetsE, forJetsD, forJetsE, caloRegions);
0862     }
0863   }
0864 
0865   if (doHF_) {
0866     iEvent.getByLabel(dataTag_.label(), hfRingSumsD);
0867     iEvent.getByLabel(emuTag_.label(), hfRingSumsE);
0868 
0869     iEvent.getByLabel(dataTag_.label(), hfBitCountsD);
0870     iEvent.getByLabel(emuTag_.label(), hfBitCountsE);
0871 
0872     isRingSumError = false;
0873     isBitCountError = false;
0874 
0875     if (checkCollections(hfRingSumsD, GCT_SUMS_QUANTA, "HF Ring Sums Data") &&
0876         checkCollections(hfRingSumsE, GCT_SUMS_QUANTA, "HF Ring Sums Emulator")) {
0877       plotHFRingSums(hfRingSumsD, hfRingSumsE);
0878       compareRingSums HFRingSums(hfRingSumsD, hfRingSumsE, MBxInfo);
0879       isRingSumError = HFRingSums.doCompare(hfRingSum_errorFlag_);
0880     }
0881 
0882     if (checkCollections(hfBitCountsD, GCT_SUMS_QUANTA, "HF Bit Counts Data") &&
0883         checkCollections(hfBitCountsE, GCT_SUMS_QUANTA, "HF Bit Counts Emulator")) {
0884       plotHFBitCounts(hfBitCountsD, hfBitCountsE);
0885       compareBitCounts HFBitCounts(hfBitCountsD, hfBitCountsE, MBxInfo);
0886       isBitCountError = HFBitCounts.doCompare(hfBitCount_errorFlag_);
0887     }
0888 
0889     if ((isRingSumError && doRingSumDebug_) || (isBitCountError && doBitCountDebug_)) {
0890       plotHFErrors(hfRingSumsD, hfRingSumsE, hfBitCountsD, hfBitCountsE, caloRegions);
0891     }
0892   }
0893 
0894   if (doTotalEnergySums_) {
0895     iEvent.getByLabel(dataTag_.label(), totalEtD);
0896     iEvent.getByLabel(emuTag_.label(), totalEtE);
0897 
0898     iEvent.getByLabel(dataTag_.label(), totalHtD);
0899     iEvent.getByLabel(emuTag_.label(), totalHtE);
0900 
0901     isTotalEError = false;
0902     isTotalHError = false;
0903 
0904     if (checkCollections(totalEtD, GCT_SUMS_QUANTA, "Total Et Data") &&
0905         checkCollections(totalEtE, GCT_SUMS_QUANTA, "Total Et Emulator")) {
0906       plotTotalE(totalEtD, totalEtE);
0907       compareTotalE compareET(totalEtD, totalEtE, MBxInfo);
0908       isTotalEError = compareET.doCompare(totalEt_errorFlag_);
0909     }
0910 
0911     if (checkCollections(totalHtD, GCT_SUMS_QUANTA, "Total Ht Data") &&
0912         checkCollections(totalHtE, GCT_SUMS_QUANTA, "Total Ht Emulator")) {
0913       plotTotalH(totalHtD, totalHtE);
0914       compareTotalH compareHT(totalHtD, totalHtE, MBxInfo);
0915       isTotalHError = compareHT.doCompare(totalHt_errorFlag_);
0916     }
0917 
0918     if ((isTotalEError && doTotalEtDebug_) || (isTotalHError && doTotalHtDebug_)) {
0919       plotTotalEErrors(totalEtD, totalEtE, totalHtD, totalHtE, caloRegions);
0920     }
0921   }
0922 
0923   if (doMissingEnergySums_) {
0924     iEvent.getByLabel(dataTag_.label(), missingEtD);
0925     iEvent.getByLabel(emuTag_.label(), missingEtE);
0926 
0927     iEvent.getByLabel(dataTag_.label(), missingHtD);
0928     iEvent.getByLabel(emuTag_.label(), missingHtE);
0929 
0930     isMissingEError = false;
0931     isMissingHError = false;
0932 
0933     if (checkCollections(missingEtD, GCT_SUMS_QUANTA, "Missing Et Data") &&
0934         checkCollections(missingEtE, GCT_SUMS_QUANTA, "Missing Et Emulator")) {
0935       plotMissingEt(missingEtD, missingEtE);
0936       compareMissingE compareMET(missingEtD, missingEtE, MBxInfo);
0937       isMissingEError = compareMET.doCompare(missingEt_errorFlag_);
0938     }
0939 
0940     if (checkCollections(missingHtD, GCT_SUMS_QUANTA, "Missing Ht Data") &&
0941         checkCollections(missingHtE, GCT_SUMS_QUANTA, "Missing Ht Emulator")) {
0942       plotMissingHt(missingHtD, missingHtE);
0943       compareMissingH compareMHT(missingHtD, missingHtE, MBxInfo);
0944       isMissingHError = compareMHT.doCompare(missingHt_errorFlag_);
0945 
0946       //added 19/03/2010 for intermediate information on MissingHt quantities in the data
0947       if (doExtraMissingHTDebug_) {
0948         iEvent.getByLabel(dataTag_.label(), "", intHtMissD);
0949         if (checkCollections(intHtMissD, GCT_INT_HTMISS_QUANTA, "Internal Missing Ht Data")) {
0950           for (unsigned int i = 0; i < intHtMissD->size(); i++) {
0951             if (doGCTMBx_ || intHtMissD->at(i).bx() == GCTTrigBx_) {
0952               if (!intHtMissD->at(i).overflow()) {
0953                 //the capBlock 0x301 is the input pipeline at the wheel for positive eta, whereas 0x701 is for negative eta
0954                 if (intHtMissD->at(i).capBlock() == 0x301 && intHtMissD->at(i).capIndex() == 0 &&
0955                     intHtMissD->at(i).isThereHtx())
0956                   missingHtD_HtXPosLeaf1->Fill(intHtMissD->at(i).htx());
0957                 if (intHtMissD->at(i).capBlock() == 0x301 && intHtMissD->at(i).capIndex() == 1 &&
0958                     intHtMissD->at(i).isThereHtx())
0959                   missingHtD_HtXPosLeaf2->Fill(intHtMissD->at(i).htx());
0960                 if (intHtMissD->at(i).capBlock() == 0x301 && intHtMissD->at(i).capIndex() == 2 &&
0961                     intHtMissD->at(i).isThereHtx())
0962                   missingHtD_HtXPosLeaf3->Fill(intHtMissD->at(i).htx());
0963                 if (intHtMissD->at(i).capBlock() == 0x701 && intHtMissD->at(i).capIndex() == 0 &&
0964                     intHtMissD->at(i).isThereHtx())
0965                   missingHtD_HtXNegLeaf1->Fill(intHtMissD->at(i).htx());
0966                 if (intHtMissD->at(i).capBlock() == 0x701 && intHtMissD->at(i).capIndex() == 1 &&
0967                     intHtMissD->at(i).isThereHtx())
0968                   missingHtD_HtXNegLeaf2->Fill(intHtMissD->at(i).htx());
0969                 if (intHtMissD->at(i).capBlock() == 0x701 && intHtMissD->at(i).capIndex() == 2 &&
0970                     intHtMissD->at(i).isThereHtx())
0971                   missingHtD_HtXNegLeaf3->Fill(intHtMissD->at(i).htx());
0972 
0973                 if (intHtMissD->at(i).capBlock() == 0x301 && intHtMissD->at(i).capIndex() == 0 &&
0974                     intHtMissD->at(i).isThereHty())
0975                   missingHtD_HtYPosLeaf1->Fill(intHtMissD->at(i).hty());
0976                 if (intHtMissD->at(i).capBlock() == 0x301 && intHtMissD->at(i).capIndex() == 1 &&
0977                     intHtMissD->at(i).isThereHty())
0978                   missingHtD_HtYPosLeaf2->Fill(intHtMissD->at(i).hty());
0979                 if (intHtMissD->at(i).capBlock() == 0x301 && intHtMissD->at(i).capIndex() == 2 &&
0980                     intHtMissD->at(i).isThereHty())
0981                   missingHtD_HtYPosLeaf3->Fill(intHtMissD->at(i).hty());
0982                 if (intHtMissD->at(i).capBlock() == 0x701 && intHtMissD->at(i).capIndex() == 0 &&
0983                     intHtMissD->at(i).isThereHty())
0984                   missingHtD_HtYNegLeaf1->Fill(intHtMissD->at(i).hty());
0985                 if (intHtMissD->at(i).capBlock() == 0x701 && intHtMissD->at(i).capIndex() == 1 &&
0986                     intHtMissD->at(i).isThereHty())
0987                   missingHtD_HtYNegLeaf2->Fill(intHtMissD->at(i).hty());
0988                 if (intHtMissD->at(i).capBlock() == 0x701 && intHtMissD->at(i).capIndex() == 2 &&
0989                     intHtMissD->at(i).isThereHty())
0990                   missingHtD_HtYNegLeaf3->Fill(intHtMissD->at(i).hty());
0991               }
0992             }
0993           }
0994         }
0995       }
0996     }
0997 
0998     if ((isMissingEError && doMissingETDebug_) || (isMissingHError && doMissingHTDebug_)) {
0999       plotMissingEErrors(missingEtD, missingEtE, missingHtD, missingHtE, caloRegions, intJetsE, intHtMissD);
1000     }
1001   }
1002 }
1003 
1004 // ------------ method called once each job just before starting event loop  ------------
1005 void GctErrorAnalyzer::beginJob() {}
1006 
1007 // ------------ method called once each job just after ending the event loop  ------------
1008 void GctErrorAnalyzer::endJob() {}
1009 
1010 void GctErrorAnalyzer::plotRCTRegions(const edm::Handle<L1CaloRegionCollection> &caloRegions) {
1011   //if more than one Bx is readout per event, then caloRegions->size() will be some multiple of 396
1012   for (unsigned int i = 0; i < caloRegions->size(); i++) {
1013     //if the RCTMBx flag is set to true, write out all the info into the same histogram
1014     //otherwise only the RCTTrigBx will be written out - could skip (RCT_REGION_QUANTA-1) events here to speed things up...
1015     if (doRCTMBx_ || caloRegions->at(i).bx() == RCTTrigBx_) {
1016       if (caloRegions->at(i).et() > 0)
1017         RCT_EtEtaPhi_->Fill(caloRegions->at(i).gctEta(), caloRegions->at(i).gctPhi(), caloRegions->at(i).et());
1018       if (caloRegions->at(i).tauVeto())
1019         RCT_TvEtaPhi_->Fill(caloRegions->at(i).gctEta(), caloRegions->at(i).gctPhi());
1020       if (caloRegions->at(i).fineGrain())
1021         RCT_FgEtaPhi_->Fill(caloRegions->at(i).gctEta(), caloRegions->at(i).gctPhi());
1022       if (caloRegions->at(i).overFlow())
1023         RCT_OfEtaPhi_->Fill(caloRegions->at(i).gctEta(), caloRegions->at(i).gctPhi());
1024     }
1025   }
1026 }
1027 
1028 void GctErrorAnalyzer::plotIsoEm(const edm::Handle<L1GctEmCandCollection> &isoEgD,
1029                                  const edm::Handle<L1GctEmCandCollection> &isoEgE) {
1030   //loop over all the data candidates - if multiple bx, then this should be a multiple of GCT_OBJECT_QUANTA
1031   for (unsigned int i = 0; i < isoEgD->size(); i++) {
1032     //if the GCTMBx flag is set, plot all Bx for this quantity, otherwise only plot Bx = GCTTrigBx_
1033     if (doGCTMBx_ || isoEgD->at(i).bx() == GCTTrigBx_) {
1034       isoEgD_Rank_->Fill(isoEgD->at(i).rank());
1035       if (isoEgD->at(i).rank() > 0) {
1036         isoEgD_EtEtaPhi_->Fill(isoEgD->at(i).regionId().ieta(), isoEgD->at(i).regionId().iphi(), isoEgD->at(i).rank());
1037         isoEgD_OccEtaPhi_->Fill(isoEgD->at(i).regionId().ieta(), isoEgD->at(i).regionId().iphi());
1038       }
1039     }
1040   }
1041   //now repeat for the emulator candidates
1042   for (unsigned int i = 0; i < isoEgE->size(); i++) {
1043     if (doEmuMBx_ || isoEgE->at(i).bx() == EmuTrigBx_) {
1044       isoEgE_Rank_->Fill(isoEgE->at(i).rank());
1045       if (isoEgE->at(i).rank() > 0) {
1046         isoEgE_EtEtaPhi_->Fill(isoEgE->at(i).regionId().ieta(), isoEgE->at(i).regionId().iphi(), isoEgE->at(i).rank());
1047         isoEgE_OccEtaPhi_->Fill(isoEgE->at(i).regionId().ieta(), isoEgE->at(i).regionId().iphi());
1048       }
1049     }
1050   }
1051 }
1052 
1053 void GctErrorAnalyzer::plotNonIsoEm(const edm::Handle<L1GctEmCandCollection> &nonIsoEgD,
1054                                     const edm::Handle<L1GctEmCandCollection> &nonIsoEgE) {
1055   //loop over all the data candidates - if multiple bx, then this should be a multiple of GCT_OBJECT_QUANTA
1056   for (unsigned int i = 0; i < nonIsoEgD->size(); i++) {
1057     //if the GCTMBx flag is set, plot all Bx for this quantity, otherwise only plot Bx = GCTTrigBx_
1058     if (doGCTMBx_ || nonIsoEgD->at(i).bx() == GCTTrigBx_) {
1059       nonIsoEgD_Rank_->Fill(nonIsoEgD->at(i).rank());
1060       if (nonIsoEgD->at(i).rank() > 0) {
1061         nonIsoEgD_EtEtaPhi_->Fill(
1062             nonIsoEgD->at(i).regionId().ieta(), nonIsoEgD->at(i).regionId().iphi(), nonIsoEgD->at(i).rank());
1063         nonIsoEgD_OccEtaPhi_->Fill(nonIsoEgD->at(i).regionId().ieta(), nonIsoEgD->at(i).regionId().iphi());
1064       }
1065     }
1066   }
1067   //now repeat for the emulator candidates
1068   for (unsigned int i = 0; i < nonIsoEgE->size(); i++) {
1069     if (doEmuMBx_ || nonIsoEgE->at(i).bx() == EmuTrigBx_) {
1070       nonIsoEgE_Rank_->Fill(nonIsoEgE->at(i).rank());
1071       if (nonIsoEgE->at(i).rank() > 0) {
1072         nonIsoEgE_EtEtaPhi_->Fill(
1073             nonIsoEgE->at(i).regionId().ieta(), nonIsoEgE->at(i).regionId().iphi(), nonIsoEgE->at(i).rank());
1074         nonIsoEgE_OccEtaPhi_->Fill(nonIsoEgE->at(i).regionId().ieta(), nonIsoEgE->at(i).regionId().iphi());
1075       }
1076     }
1077   }
1078 }
1079 
1080 void GctErrorAnalyzer::plotEGErrors(const edm::Handle<L1GctEmCandCollection> &isoEgD,
1081                                     const edm::Handle<L1GctEmCandCollection> &isoEgE,
1082                                     const edm::Handle<L1GctEmCandCollection> &nonIsoEgD,
1083                                     const edm::Handle<L1GctEmCandCollection> &nonIsoEgE,
1084                                     const edm::Handle<L1CaloEmCollection> &emRegions) {
1085   std::string errorDirName = "err_";
1086   if (isIsoError)
1087     errorDirName.append("I");
1088   if (isNonIsoError)
1089     errorDirName.append("N");
1090   std::stringstream caseNumber;
1091   caseNumber << eventNumber;
1092   errorDirName.append(caseNumber.str());
1093   TFileDirectory errorDir = errorHistCat.at(0).mkdir(errorDirName);
1094 
1095   TH2I *errorEmRegionIsoEtEtaPhi_ = errorDir.make<TH2I>("errorEmRegionIsoEtEtaPhi",
1096                                                         "errorEmRegionIsoEtEtaPhi;#eta (GCT Units);#phi (GCT Units)",
1097                                                         22,
1098                                                         -0.5,
1099                                                         21.5,
1100                                                         18,
1101                                                         -0.5,
1102                                                         17.5);
1103   TH2I *errorEmRegionNonIsoEtEtaPhi_ =
1104       errorDir.make<TH2I>("errorEmRegionNonIsoEtEtaPhi",
1105                           "errorEmRegionNonIsoEtEtaPhi;#eta (GCT Units);#phi (GCT Units)",
1106                           22,
1107                           -0.5,
1108                           21.5,
1109                           18,
1110                           -0.5,
1111                           17.5);
1112   TH2I *errorIsoEtEtaPhiD_ = errorDir.make<TH2I>(
1113       "errorIsoEtEtaPhiD", "errorIsoEtEtaPhiD;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
1114   TH2I *errorNonIsoEtEtaPhiD_ = errorDir.make<TH2I>(
1115       "errorNonIsoEtEtaPhiD", "errorNonIsoEtEtaPhiD;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
1116   TH2I *errorIsoEtEtaPhiE_ = errorDir.make<TH2I>(
1117       "errorIsoEtEtaPhiE", "errorIsoEtEtaPhiE;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
1118   TH2I *errorNonIsoEtEtaPhiE_ = errorDir.make<TH2I>(
1119       "errorNonIsoEtEtaPhiE", "errorNonIsoEtEtaPhiE;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
1120 
1121   //fill the EM input collection
1122   //should only fill the correct bx for emRegions - and since this is showing an error in the comparison, we should plot the input to this comparison i.e. Bx=RCTTrigBx
1123   //this assumes that comparison is done on the central Bx i.e. RctBx=0 corresponds to GctBx=0, and EmuBx=0 takes RctBx=0
1124   for (unsigned int i = 0; i < emRegions->size(); i++) {
1125     if (emRegions->at(i).bx() == RCTTrigBx_) {
1126       if (emRegions->at(i).isolated()) {
1127         if (emRegions->at(i).rank() > 0)
1128           errorEmRegionIsoEtEtaPhi_->Fill(
1129               emRegions->at(i).regionId().ieta(), emRegions->at(i).regionId().iphi(), emRegions->at(i).rank());
1130       } else {
1131         if (emRegions->at(i).rank() > 0)
1132           errorEmRegionNonIsoEtEtaPhi_->Fill(
1133               emRegions->at(i).regionId().ieta(), emRegions->at(i).regionId().iphi(), emRegions->at(i).rank());
1134       }
1135     }
1136   }
1137 
1138   //no need to have the rank plot, because you can't have two electrons in the same place (eta,phi), in the same event...
1139   //in this case, since we're actually comparing the GCTTrigBx_ with the EmuTrigBx_, we plot these individually
1140   for (unsigned int i = 0; i < isoEgD->size(); i++) {
1141     if (isoEgD->at(i).bx() == GCTTrigBx_) {
1142       if (isoEgD->at(i).rank() > 0)
1143         errorIsoEtEtaPhiD_->Fill(
1144             isoEgD->at(i).regionId().ieta(), isoEgD->at(i).regionId().iphi(), isoEgD->at(i).rank());
1145     }
1146   }
1147   for (unsigned int i = 0; i < nonIsoEgD->size(); i++) {
1148     if (nonIsoEgD->at(i).bx() == GCTTrigBx_) {
1149       if (nonIsoEgD->at(i).rank() > 0)
1150         errorNonIsoEtEtaPhiD_->Fill(
1151             nonIsoEgD->at(i).regionId().ieta(), nonIsoEgD->at(i).regionId().iphi(), nonIsoEgD->at(i).rank());
1152     }
1153   }
1154 
1155   //now for the emulator candidates
1156   for (unsigned int i = 0; i < isoEgE->size(); i++) {
1157     if (isoEgE->at(i).bx() == EmuTrigBx_) {
1158       if (isoEgE->at(i).rank() > 0)
1159         errorIsoEtEtaPhiE_->Fill(
1160             isoEgE->at(i).regionId().ieta(), isoEgE->at(i).regionId().iphi(), isoEgE->at(i).rank());
1161     }
1162   }
1163   for (unsigned int i = 0; i < nonIsoEgE->size(); i++) {
1164     if (nonIsoEgE->at(i).bx() == EmuTrigBx_) {
1165       if (nonIsoEgE->at(i).rank() > 0)
1166         errorNonIsoEtEtaPhiE_->Fill(
1167             nonIsoEgE->at(i).regionId().ieta(), nonIsoEgE->at(i).regionId().iphi(), nonIsoEgE->at(i).rank());
1168     }
1169   }
1170 }
1171 
1172 void GctErrorAnalyzer::plotCenJets(const edm::Handle<L1GctJetCandCollection> &cenJetsD,
1173                                    const edm::Handle<L1GctJetCandCollection> &cenJetsE) {
1174   for (unsigned int i = 0; i < cenJetsD->size(); i++) {
1175     if (doGCTMBx_ || cenJetsD->at(i).bx() == GCTTrigBx_) {
1176       cenJetD_Rank_->Fill(cenJetsD->at(i).rank());
1177       if (cenJetsD->at(i).rank() > 0) {
1178         cenJetD_EtEtaPhi_->Fill(
1179             cenJetsD->at(i).regionId().ieta(), cenJetsD->at(i).regionId().iphi(), cenJetsD->at(i).rank());
1180         cenJetD_OccEtaPhi_->Fill(cenJetsD->at(i).regionId().ieta(), cenJetsD->at(i).regionId().iphi());
1181       }
1182     }
1183   }
1184 
1185   for (unsigned int i = 0; i < cenJetsE->size(); i++) {
1186     if (doEmuMBx_ || cenJetsE->at(i).bx() == EmuTrigBx_) {
1187       cenJetE_Rank_->Fill(cenJetsE->at(i).rank());
1188       if (cenJetsE->at(i).rank() > 0) {
1189         cenJetE_EtEtaPhi_->Fill(
1190             cenJetsE->at(i).regionId().ieta(), cenJetsE->at(i).regionId().iphi(), cenJetsE->at(i).rank());
1191         cenJetE_OccEtaPhi_->Fill(cenJetsE->at(i).regionId().ieta(), cenJetsE->at(i).regionId().iphi());
1192       }
1193     }
1194   }
1195 }
1196 
1197 void GctErrorAnalyzer::plotTauJets(const edm::Handle<L1GctJetCandCollection> &tauJetsD,
1198                                    const edm::Handle<L1GctJetCandCollection> &tauJetsE) {
1199   for (unsigned int i = 0; i < tauJetsD->size(); i++) {
1200     if (doGCTMBx_ || tauJetsD->at(i).bx() == GCTTrigBx_) {
1201       tauJetD_Rank_->Fill(tauJetsD->at(i).rank());
1202       if (tauJetsD->at(i).rank() > 0) {
1203         tauJetD_EtEtaPhi_->Fill(
1204             tauJetsD->at(i).regionId().ieta(), tauJetsD->at(i).regionId().iphi(), tauJetsD->at(i).rank());
1205         tauJetD_OccEtaPhi_->Fill(tauJetsD->at(i).regionId().ieta(), tauJetsD->at(i).regionId().iphi());
1206       }
1207     }
1208   }
1209 
1210   for (unsigned int i = 0; i < tauJetsE->size(); i++) {
1211     if (doEmuMBx_ || tauJetsE->at(i).bx() == EmuTrigBx_) {
1212       tauJetE_Rank_->Fill(tauJetsE->at(i).rank());
1213       if (tauJetsE->at(i).rank() > 0) {
1214         tauJetE_EtEtaPhi_->Fill(
1215             tauJetsE->at(i).regionId().ieta(), tauJetsE->at(i).regionId().iphi(), tauJetsE->at(i).rank());
1216         tauJetE_OccEtaPhi_->Fill(tauJetsE->at(i).regionId().ieta(), tauJetsE->at(i).regionId().iphi());
1217       }
1218     }
1219   }
1220 }
1221 
1222 void GctErrorAnalyzer::plotForJets(const edm::Handle<L1GctJetCandCollection> &forJetsD,
1223                                    const edm::Handle<L1GctJetCandCollection> &forJetsE) {
1224   for (unsigned int i = 0; i < forJetsD->size(); i++) {
1225     if (doGCTMBx_ || forJetsD->at(i).bx() == GCTTrigBx_) {
1226       forJetD_Rank_->Fill(forJetsD->at(i).rank());
1227       if (forJetsD->at(i).rank() > 0) {
1228         forJetD_EtEtaPhi_->Fill(
1229             forJetsD->at(i).regionId().ieta(), forJetsD->at(i).regionId().iphi(), forJetsD->at(i).rank());
1230         forJetD_OccEtaPhi_->Fill(forJetsD->at(i).regionId().ieta(), forJetsD->at(i).regionId().iphi());
1231       }
1232     }
1233   }
1234 
1235   for (unsigned int i = 0; i < forJetsE->size(); i++) {
1236     if (doEmuMBx_ || forJetsE->at(i).bx() == EmuTrigBx_) {
1237       forJetE_Rank_->Fill(forJetsE->at(i).rank());
1238       if (forJetsE->at(i).rank() > 0) {
1239         forJetE_EtEtaPhi_->Fill(
1240             forJetsE->at(i).regionId().ieta(), forJetsE->at(i).regionId().iphi(), forJetsE->at(i).rank());
1241         forJetE_OccEtaPhi_->Fill(forJetsE->at(i).regionId().ieta(), forJetsE->at(i).regionId().iphi());
1242       }
1243     }
1244   }
1245 }
1246 
1247 void GctErrorAnalyzer::plotIntJets(const edm::Handle<L1GctInternJetDataCollection> &intJetsE) {
1248   jetData intJet;
1249   std::vector<jetData> intJetCollection(
1250       NUM_INT_JETS);  //define fixed size for the vector to avoid reallocation (i.e. max size possible)
1251 
1252   //since we don't read out the intermediate (i.e. leaf card) jets, we can only plot the emulator distributions
1253   //the 1st-4th jet Et will prove useful in understanding and motivating cuts on individual jets in HT and MHT.
1254   for (unsigned int i = 0; i < intJetsE->size(); i++) {
1255     if (doEmuMBx_ || intJetsE->at(i).bx() == EmuTrigBx_) {
1256       //the intermediate jets are not sorted in terms of Et so
1257       //in order to do this independently of the data format,
1258       //copy to a user defined struct and sort that way
1259       intJet.et = intJetsE->at(i).et();
1260       intJet.phi = intJetsE->at(i).phi();
1261       intJet.eta = intJetsE->at(i).eta();
1262       intJetCollection.at(i % NUM_INT_JETS) = intJet;
1263 
1264       //remember, if the event has 1 overflowed jet, then we fill the internal jet dist overflow histogram
1265       //and skip the event - there is no point looking at the leading jet distributions etc for an event
1266       //with an overflowed jet - this will imply HT, ET, MET and MHT all overflow too.
1267       if (intJetsE->at(i).oflow()) {
1268         intJetE_Of_->Fill(intJetsE->at(i).oflow());
1269         return;
1270       }
1271 
1272       //plot the (et,eta,phi) distribution of the intermediate jets (for non-zero et)
1273       if (intJetsE->at(i).et())
1274         intJetEtEtaPhiE_->Fill(
1275             intJetsE->at(i).regionId().ieta(), intJetsE->at(i).regionId().iphi(), intJetsE->at(i).et());
1276     }
1277   }
1278 
1279   //if we get this far, there are no jets with an overflow bit so fill the overflow histogram and
1280   //sort the intJetCollection according to the rule defined in sortJets (i.e. largest et first)
1281   intJetE_Of_->Fill(0);
1282   std::sort(intJetCollection.begin(), intJetCollection.end(), sortJets);
1283 
1284   std::vector<TH1I *> leadingJetDist(4);
1285   leadingJetDist.at(0) = intJetE_Jet1Et_;
1286   leadingJetDist.at(1) = intJetE_Jet2Et_;
1287   leadingJetDist.at(2) = intJetE_Jet3Et_;
1288   leadingJetDist.at(3) = intJetE_Jet4Et_;
1289 
1290   unsigned int i = 0;
1291   unsigned int j = 0;
1292   unsigned int currentEt = 0;
1293   while (intJetCollection.at(i).et > 0) {
1294     if (j < leadingJetDist.size()) {
1295       if (i == 0) {
1296         leadingJetDist.at(j)->Fill(intJetCollection.at(i).et);
1297         currentEt = intJetCollection.at(i).et;
1298         j++;
1299       } else {
1300         if (intJetCollection.at(i).et < currentEt) {
1301           leadingJetDist.at(j)->Fill(intJetCollection.at(i).et);
1302           currentEt = intJetCollection.at(i).et;
1303           j++;
1304         }
1305       }
1306     }
1307 
1308     intJetE_Et_->Fill(intJetCollection.at(i).et);
1309     i++;
1310   }
1311   return;
1312 }
1313 
1314 bool GctErrorAnalyzer::sortJets(const jetData &jet1, const jetData &jet2) { return jet1.et > jet2.et; }
1315 
1316 template <class T>
1317 bool GctErrorAnalyzer::checkCollections(const T &collection, const unsigned int &constraint, const std::string &label) {
1318   //unfortunately, the dataformats are not consistent with the name() method (i.e. some have it, others don't)
1319   //and a typeof() function doesn't exist in ANSI C++, so to identify the templated type, we pass a std::string
1320 
1321   if (!collection.isValid()) {
1322     edm::LogWarning("DataNotFound") << " Could not find " << label << " label";
1323     return false;
1324   }
1325   if (collection->size() % constraint != 0 || collection->empty()) {
1326     edm::LogWarning("CollectionSizeError")
1327         << " " << label << " collection size is " << collection->size() << ", expected multiple of " << constraint;
1328     return false;
1329   }
1330 
1331   return true;
1332 }
1333 
1334 void GctErrorAnalyzer::plotJetErrors(const edm::Handle<L1GctJetCandCollection> &cenJetsD,
1335                                      const edm::Handle<L1GctJetCandCollection> &cenJetsE,
1336                                      const edm::Handle<L1GctJetCandCollection> &tauJetsD,
1337                                      const edm::Handle<L1GctJetCandCollection> &tauJetsE,
1338                                      const edm::Handle<L1GctJetCandCollection> &forJetsD,
1339                                      const edm::Handle<L1GctJetCandCollection> &forJetsE,
1340                                      const edm::Handle<L1CaloRegionCollection> &caloRegions) {
1341   std::string errorDirName = "err_";
1342   if (isCenJetError)
1343     errorDirName.append("C");
1344   if (isTauJetError)
1345     errorDirName.append("T");
1346   if (isForJetError)
1347     errorDirName.append("F");
1348   std::stringstream caseNumber;
1349   caseNumber << eventNumber;
1350   errorDirName.append(caseNumber.str());
1351   TFileDirectory errorDir = errorHistCat.at(1).mkdir(errorDirName);
1352 
1353   TH2I *errorRegionEtEtaPhi_ = errorDir.make<TH2I>(
1354       "errorRegionEtEtaPhi", "errorRegionEtEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
1355   TH2I *errorRegionTvEtaPhi_ = errorDir.make<TH2I>(
1356       "errorRegionTvEtaPhi", "errorRegionTvEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
1357   TH2I *errorRegionOfEtaPhi_ = errorDir.make<TH2I>(
1358       "errorRegionOfEtaPhi", "errorRegionOfEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
1359 
1360   //make sure to only plot the caloRegion bx which corresponds to the data vs emulator comparison
1361   for (unsigned int i = 0; i < caloRegions->size(); i++) {
1362     if (caloRegions->at(i).bx() == RCTTrigBx_) {
1363       if (caloRegions->at(i).et() > 0)
1364         errorRegionEtEtaPhi_->Fill(caloRegions->at(i).gctEta(), caloRegions->at(i).gctPhi(), caloRegions->at(i).et());
1365       if (caloRegions->at(i).tauVeto())
1366         errorRegionTvEtaPhi_->Fill(caloRegions->at(i).gctEta(), caloRegions->at(i).gctPhi());
1367       if (caloRegions->at(i).overFlow())
1368         errorRegionOfEtaPhi_->Fill(caloRegions->at(i).gctEta(), caloRegions->at(i).gctPhi());
1369     }
1370   }
1371 
1372   TH2I *cenJet_errorEtEtaPhiData_ = errorDir.make<TH2I>("cenJet_errorEtEtaPhiData",
1373                                                         "cenJet_errorEtEtaPhiData;#eta (GCT Units);#phi (GCT Units)",
1374                                                         22,
1375                                                         -0.5,
1376                                                         21.5,
1377                                                         18,
1378                                                         -0.5,
1379                                                         17.5);
1380   TH2I *cenJet_errorEtEtaPhiEmu_ = errorDir.make<TH2I>("cenJet_errorEtEtaPhiEmu",
1381                                                        "cenJet_errorEtEtaPhiEmu;#eta (GCT Units);#phi (GCT Units)",
1382                                                        22,
1383                                                        -0.5,
1384                                                        21.5,
1385                                                        18,
1386                                                        -0.5,
1387                                                        17.5);
1388   TH2I *tauJet_errorEtEtaPhiData_ = errorDir.make<TH2I>("tauJet_errorEtEtaPhiData",
1389                                                         "tauJet_errorEtEtaPhiData;#eta (GCT Units);#phi (GCT Units)",
1390                                                         22,
1391                                                         -0.5,
1392                                                         21.5,
1393                                                         18,
1394                                                         -0.5,
1395                                                         17.5);
1396   TH2I *tauJet_errorEtEtaPhiEmu_ = errorDir.make<TH2I>("tauJet_errorEtEtaPhiEmu",
1397                                                        "tauJet_errorEtEtaPhiEmu;#eta (GCT Units);#phi (GCT Units)",
1398                                                        22,
1399                                                        -0.5,
1400                                                        21.5,
1401                                                        18,
1402                                                        -0.5,
1403                                                        17.5);
1404   TH2I *forJet_errorEtEtaPhiData_ = errorDir.make<TH2I>("forJet_errorEtEtaPhiData",
1405                                                         "forJet_errorEtEtaPhiData;#eta (GCT Units);#phi (GCT Units)",
1406                                                         22,
1407                                                         -0.5,
1408                                                         21.5,
1409                                                         18,
1410                                                         -0.5,
1411                                                         17.5);
1412   TH2I *forJet_errorEtEtaPhiEmu_ = errorDir.make<TH2I>("forJet_errorEtEtaPhiEmu",
1413                                                        "forJet_errorEtEtaPhiEmu;#eta (GCT Units);#phi (GCT Units)",
1414                                                        22,
1415                                                        -0.5,
1416                                                        21.5,
1417                                                        18,
1418                                                        -0.5,
1419                                                        17.5);
1420 
1421   //first plot the data candiates for the Trigger Bx that this error corresponds to
1422   for (unsigned int i = 0; i < cenJetsD->size(); i++) {
1423     if (cenJetsD->at(i).bx() == GCTTrigBx_) {
1424       if (cenJetsD->at(i).rank() > 0)
1425         cenJet_errorEtEtaPhiData_->Fill(
1426             cenJetsD->at(i).regionId().ieta(), cenJetsD->at(i).regionId().iphi(), cenJetsD->at(i).rank());
1427     }
1428   }
1429   for (unsigned int i = 0; i < tauJetsD->size(); i++) {
1430     if (tauJetsD->at(i).bx() == GCTTrigBx_) {
1431       if (tauJetsD->at(i).rank() > 0)
1432         tauJet_errorEtEtaPhiData_->Fill(
1433             tauJetsD->at(i).regionId().ieta(), tauJetsD->at(i).regionId().iphi(), tauJetsD->at(i).rank());
1434     }
1435   }
1436   for (unsigned int i = 0; i < forJetsD->size(); i++) {
1437     if (forJetsD->at(i).bx() == GCTTrigBx_) {
1438       if (forJetsD->at(i).rank() > 0)
1439         forJet_errorEtEtaPhiData_->Fill(
1440             forJetsD->at(i).regionId().ieta(), forJetsD->at(i).regionId().iphi(), forJetsD->at(i).rank());
1441     }
1442   }
1443 
1444   //now the emulator candidates
1445   for (unsigned int i = 0; i < cenJetsE->size(); i++) {
1446     if (cenJetsE->at(i).bx() == EmuTrigBx_) {
1447       if (cenJetsE->at(i).rank() > 0)
1448         cenJet_errorEtEtaPhiEmu_->Fill(
1449             cenJetsE->at(i).regionId().ieta(), cenJetsE->at(i).regionId().iphi(), cenJetsE->at(i).rank());
1450     }
1451   }
1452   for (unsigned int i = 0; i < tauJetsE->size(); i++) {
1453     if (tauJetsE->at(i).bx() == EmuTrigBx_) {
1454       if (tauJetsE->at(i).rank() > 0)
1455         tauJet_errorEtEtaPhiEmu_->Fill(
1456             tauJetsE->at(i).regionId().ieta(), tauJetsE->at(i).regionId().iphi(), tauJetsE->at(i).rank());
1457     }
1458   }
1459   for (unsigned int i = 0; i < forJetsE->size(); i++) {
1460     if (forJetsE->at(i).bx() == EmuTrigBx_) {
1461       if (forJetsE->at(i).rank() > 0)
1462         forJet_errorEtEtaPhiEmu_->Fill(
1463             forJetsE->at(i).regionId().ieta(), forJetsE->at(i).regionId().iphi(), forJetsE->at(i).rank());
1464     }
1465   }
1466 }
1467 
1468 void GctErrorAnalyzer::plotHFRingSums(const edm::Handle<L1GctHFRingEtSumsCollection> &hfRingSumsD,
1469                                       const edm::Handle<L1GctHFRingEtSumsCollection> &hfRingSumsE) {
1470   for (unsigned int i = 0; i < hfRingSumsD->size(); i++) {
1471     if (doGCTMBx_ || hfRingSumsD->at(i).bx() == GCTTrigBx_) {
1472       //there are 4 rings - just fill the histograms
1473       hfRingSumD_1pos_->Fill(hfRingSumsD->at(i).etSum(0));
1474       hfRingSumD_1neg_->Fill(hfRingSumsD->at(i).etSum(1));
1475       hfRingSumD_2pos_->Fill(hfRingSumsD->at(i).etSum(2));
1476       hfRingSumD_2neg_->Fill(hfRingSumsD->at(i).etSum(3));
1477     }
1478   }
1479 
1480   for (unsigned int i = 0; i < hfRingSumsE->size(); i++) {
1481     if (doEmuMBx_ || hfRingSumsE->at(i).bx() == EmuTrigBx_) {
1482       hfRingSumE_1pos_->Fill(hfRingSumsE->at(i).etSum(0));
1483       hfRingSumE_1neg_->Fill(hfRingSumsE->at(i).etSum(1));
1484       hfRingSumE_2pos_->Fill(hfRingSumsE->at(i).etSum(2));
1485       hfRingSumE_2neg_->Fill(hfRingSumsE->at(i).etSum(3));
1486     }
1487   }
1488 }
1489 
1490 void GctErrorAnalyzer::plotHFBitCounts(const edm::Handle<L1GctHFBitCountsCollection> &hfBitCountsD,
1491                                        const edm::Handle<L1GctHFBitCountsCollection> &hfBitCountsE) {
1492   for (unsigned int i = 0; i < hfBitCountsD->size(); i++) {
1493     if (doGCTMBx_ || hfBitCountsD->at(i).bx() == GCTTrigBx_) {
1494       //there are 4 rings - just fill the histograms
1495       hfBitCountD_1pos_->Fill(hfBitCountsD->at(i).bitCount(0));
1496       hfBitCountD_1neg_->Fill(hfBitCountsD->at(i).bitCount(1));
1497       hfBitCountD_2pos_->Fill(hfBitCountsD->at(i).bitCount(2));
1498       hfBitCountD_2neg_->Fill(hfBitCountsD->at(i).bitCount(3));
1499     }
1500   }
1501   for (unsigned int i = 0; i < hfBitCountsE->size(); i++) {
1502     if (doEmuMBx_ || hfBitCountsE->at(i).bx() == EmuTrigBx_) {
1503       hfBitCountE_1pos_->Fill(hfBitCountsE->at(i).bitCount(0));
1504       hfBitCountE_1neg_->Fill(hfBitCountsE->at(i).bitCount(1));
1505       hfBitCountE_2pos_->Fill(hfBitCountsE->at(i).bitCount(2));
1506       hfBitCountE_2neg_->Fill(hfBitCountsE->at(i).bitCount(3));
1507     }
1508   }
1509 }
1510 
1511 void GctErrorAnalyzer::plotHFErrors(const edm::Handle<L1GctHFRingEtSumsCollection> &hfRingSumsD,
1512                                     const edm::Handle<L1GctHFRingEtSumsCollection> &hfRingSumsE,
1513                                     const edm::Handle<L1GctHFBitCountsCollection> &hfBitCountsD,
1514                                     const edm::Handle<L1GctHFBitCountsCollection> &hfBitCountsE,
1515                                     const edm::Handle<L1CaloRegionCollection> &caloRegions) {
1516   std::string errorDirName = "err_";
1517   if (isRingSumError)
1518     errorDirName.append("R");
1519   if (isBitCountError)
1520     errorDirName.append("B");
1521   std::stringstream caseNumber;
1522   caseNumber << eventNumber;
1523   errorDirName.append(caseNumber.str());
1524   TFileDirectory errorDir = errorHistCat.at(2).mkdir(errorDirName);
1525 
1526   TH2I *errorRegionEtEtaPhi_ = errorDir.make<TH2I>(
1527       "errorRegionEtEtaPhi", "errorRegionEtEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
1528   TH2I *errorRegionFgEtaPhi_ = errorDir.make<TH2I>(
1529       "errorRegionFgEtaPhi", "errorRegionFgEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
1530   TH2I *errorRegionOfEtaPhi_ = errorDir.make<TH2I>(
1531       "errorRegionOfEtaPhi", "errorRegionOfEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
1532 
1533   TH1I *errorHFRingSumD_1pos_ =
1534       errorDir.make<TH1I>("errorHFRingSumD_1+", "errorHFRingSumD_1+;Rank;Number of Events", 8, -0.5, 7.5);
1535   TH1I *errorHFRingSumD_2pos_ =
1536       errorDir.make<TH1I>("errorHFRingSumD_2+", "errorHFRingSumD_2+;Rank;Number of Events", 8, -0.5, 7.5);
1537   TH1I *errorHFRingSumD_1neg_ =
1538       errorDir.make<TH1I>("errorHFRingSumD_1-", "errorHFRingSumD_1-;Rank;Number of Events", 8, -0.5, 7.5);
1539   TH1I *errorHFRingSumD_2neg_ =
1540       errorDir.make<TH1I>("errorHFRingSumD_2-", "errorHFRingSumD_2-;Rank;Number of Events", 8, -0.5, 7.5);
1541   TH1I *errorHFRingSumE_1pos_ =
1542       errorDir.make<TH1I>("errorHFRingSumE_1+", "errorHFRingSumE_1+;Rank;Number of Events", 8, -0.5, 7.5);
1543   TH1I *errorHFRingSumE_2pos_ =
1544       errorDir.make<TH1I>("errorHFRingSumE_2+", "errorHFRingSumE_2+;Rank;Number of Events", 8, -0.5, 7.5);
1545   TH1I *errorHFRingSumE_1neg_ =
1546       errorDir.make<TH1I>("errorHFRingSumE_1-", "errorHFRingSumE_1-;Rank;Number of Events", 8, -0.5, 7.5);
1547   TH1I *errorHFRingSumE_2neg_ =
1548       errorDir.make<TH1I>("errorHFRingSumE_2-", "errorHFRingSumE_2-;Rank;Number of Events", 8, -0.5, 7.5);
1549 
1550   TH1I *errorHFBitCountD_1pos_ =
1551       errorDir.make<TH1I>("errorHFBitCountD_1+", "errorHFBitCountD_1+;Rank;Number of Events", 8, -0.5, 7.5);
1552   TH1I *errorHFBitCountD_2pos_ =
1553       errorDir.make<TH1I>("errorHFBitCountD_2+", "errorHFBitCountD_2+;Rank;Number of Events", 8, -0.5, 7.5);
1554   TH1I *errorHFBitCountD_1neg_ =
1555       errorDir.make<TH1I>("errorHFBitCountD_1-", "errorHFBitCountD_1-;Rank;Number of Events", 8, -0.5, 7.5);
1556   TH1I *errorHFBitCountD_2neg_ =
1557       errorDir.make<TH1I>("errorHFBitCountD_2-", "errorHFBitCountD_2-;Rank;Number of Events", 8, -0.5, 7.5);
1558   TH1I *errorHFBitCountE_1pos_ =
1559       errorDir.make<TH1I>("errorHFBitCountE_1+", "errorHFBitCountE_1+;Rank;Number of Events", 8, -0.5, 7.5);
1560   TH1I *errorHFBitCountE_2pos_ =
1561       errorDir.make<TH1I>("errorHFBitCountE_2+", "errorHFBitCountE_2+;Rank;Number of Events", 8, -0.5, 7.5);
1562   TH1I *errorHFBitCountE_1neg_ =
1563       errorDir.make<TH1I>("errorHFBitCountE_1-", "errorHFBitCountE_1-;Rank;Number of Events", 8, -0.5, 7.5);
1564   TH1I *errorHFBitCountE_2neg_ =
1565       errorDir.make<TH1I>("errorHFBitCountE_2-", "errorHFBitCountE_2-;Rank;Number of Events", 8, -0.5, 7.5);
1566 
1567   for (unsigned int i = 0; i < caloRegions->size(); i++) {
1568     if (caloRegions->at(i).bx() == RCTTrigBx_) {
1569       if (caloRegions->at(i).et() > 0)
1570         errorRegionEtEtaPhi_->Fill(caloRegions->at(i).gctEta(), caloRegions->at(i).gctPhi(), caloRegions->at(i).et());
1571       if (caloRegions->at(i).fineGrain())
1572         errorRegionFgEtaPhi_->Fill(caloRegions->at(i).gctEta(), caloRegions->at(i).gctPhi());
1573       if (caloRegions->at(i).overFlow())
1574         errorRegionOfEtaPhi_->Fill(caloRegions->at(i).gctEta(), caloRegions->at(i).gctPhi());
1575     }
1576   }
1577 
1578   for (unsigned int i = 0; i < hfRingSumsD->size(); i++) {
1579     if (hfRingSumsD->at(i).bx() == GCTTrigBx_) {
1580       errorHFRingSumD_1pos_->Fill(hfRingSumsD->at(i).etSum(0));
1581       errorHFRingSumD_1neg_->Fill(hfRingSumsD->at(i).etSum(1));
1582       errorHFRingSumD_2pos_->Fill(hfRingSumsD->at(i).etSum(2));
1583       errorHFRingSumD_2neg_->Fill(hfRingSumsD->at(i).etSum(3));
1584     }
1585   }
1586   for (unsigned int i = 0; i < hfRingSumsE->size(); i++) {
1587     if (hfRingSumsE->at(i).bx() == EmuTrigBx_) {
1588       errorHFRingSumE_1pos_->Fill(hfRingSumsE->at(i).etSum(0));
1589       errorHFRingSumE_1neg_->Fill(hfRingSumsE->at(i).etSum(1));
1590       errorHFRingSumE_2pos_->Fill(hfRingSumsE->at(i).etSum(2));
1591       errorHFRingSumE_2neg_->Fill(hfRingSumsE->at(i).etSum(3));
1592     }
1593   }
1594 
1595   for (unsigned int i = 0; i < hfBitCountsD->size(); i++) {
1596     if (hfBitCountsD->at(i).bx() == GCTTrigBx_) {
1597       errorHFBitCountD_1pos_->Fill(hfBitCountsD->at(i).bitCount(0));
1598       errorHFBitCountD_1neg_->Fill(hfBitCountsD->at(i).bitCount(1));
1599       errorHFBitCountD_2pos_->Fill(hfBitCountsD->at(i).bitCount(2));
1600       errorHFBitCountD_2neg_->Fill(hfBitCountsD->at(i).bitCount(3));
1601     }
1602   }
1603   for (unsigned int i = 0; i < hfBitCountsE->size(); i++) {
1604     if (hfBitCountsE->at(i).bx() == EmuTrigBx_) {
1605       errorHFBitCountE_1pos_->Fill(hfBitCountsE->at(i).bitCount(0));
1606       errorHFBitCountE_1neg_->Fill(hfBitCountsE->at(i).bitCount(1));
1607       errorHFBitCountE_2pos_->Fill(hfBitCountsE->at(i).bitCount(2));
1608       errorHFBitCountE_2neg_->Fill(hfBitCountsE->at(i).bitCount(3));
1609     }
1610   }
1611 }
1612 
1613 void GctErrorAnalyzer::plotTotalE(const edm::Handle<L1GctEtTotalCollection> &totalEtD,
1614                                   const edm::Handle<L1GctEtTotalCollection> &totalEtE) {
1615   for (unsigned int i = 0; i < totalEtD->size(); i++) {
1616     if (doGCTMBx_ || totalEtD->at(i).bx() == GCTTrigBx_) {
1617       totalEtD_Of_->Fill(totalEtD->at(i).overFlow());
1618       if (!totalEtD->at(i).overFlow())
1619         totalEtD_->Fill(totalEtD->at(i).et());
1620     }
1621   }
1622   for (unsigned int i = 0; i < totalEtE->size(); i++) {
1623     if (doEmuMBx_ || totalEtE->at(i).bx() == EmuTrigBx_) {
1624       totalEtE_Of_->Fill(totalEtE->at(i).overFlow());
1625       if (!totalEtE->at(i).overFlow())
1626         totalEtE_->Fill(totalEtE->at(i).et());
1627     }
1628   }
1629 }
1630 
1631 void GctErrorAnalyzer::plotTotalH(const edm::Handle<L1GctEtHadCollection> &totalHtD,
1632                                   const edm::Handle<L1GctEtHadCollection> &totalHtE) {
1633   for (unsigned int i = 0; i < totalHtD->size(); i++) {
1634     if (doGCTMBx_ || totalHtD->at(i).bx() == GCTTrigBx_) {
1635       totalHtD_Of_->Fill(totalHtD->at(i).overFlow());
1636       if (!totalHtD->at(i).overFlow())
1637         totalHtD_->Fill(totalHtD->at(i).et());
1638     }
1639   }
1640   for (unsigned int i = 0; i < totalHtE->size(); i++) {
1641     if (doEmuMBx_ || totalHtE->at(i).bx() == EmuTrigBx_) {
1642       totalHtE_Of_->Fill(totalHtE->at(i).overFlow());
1643       if (!totalHtE->at(i).overFlow())
1644         totalHtE_->Fill(totalHtE->at(i).et());
1645     }
1646   }
1647 }
1648 
1649 void GctErrorAnalyzer::plotTotalEErrors(const edm::Handle<L1GctEtTotalCollection> &totalEtD,
1650                                         const edm::Handle<L1GctEtTotalCollection> &totalEtE,
1651                                         const edm::Handle<L1GctEtHadCollection> &totalHtD,
1652                                         const edm::Handle<L1GctEtHadCollection> &totalHtE,
1653                                         const edm::Handle<L1CaloRegionCollection> &caloRegions) {
1654   std::string errorDirName = "err_";
1655   if (isTotalEError)
1656     errorDirName.append("E");
1657   if (isTotalHError)
1658     errorDirName.append("H");
1659   std::stringstream caseNumber;
1660   caseNumber << eventNumber;
1661   errorDirName.append(caseNumber.str());
1662   TFileDirectory errorDir = errorHistCat.at(3).mkdir(errorDirName);
1663 
1664   TH2I *errorRegionEtEtaPhi_ = errorDir.make<TH2I>(
1665       "errorRegionEtEtaPhi", "errorRegionEtEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
1666   TH2I *errorRegionOfEtaPhi_ = errorDir.make<TH2I>(
1667       "errorRegionOfEtaPhi", "errorRegionOfEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
1668   TH1I *errorTotalEtD_ =
1669       errorDir.make<TH1I>("errorTotalEtD", "errorTotalEtD;E_{T};Number of Events", 1024, -0.5, 1023.5);
1670   TH1I *errorTotalEtD_Of_ =
1671       errorDir.make<TH1I>("errorTotalEtD_Of", "errorTotalEtD_Of;Overflow Bit Status;Number of Events", 2, -0.5, 1.5);
1672   TH1I *errorTotalEtE_ =
1673       errorDir.make<TH1I>("errorTotalEtE", "errorTotalEtE;E_{T};Number of Events", 1024, -0.5, 1023.5);
1674   TH1I *errorTotalEtE_Of_ =
1675       errorDir.make<TH1I>("errorTotalEtE_Of", "errorTotalEtE_Of;Overflow Bit Status;Number of Events", 2, -0.5, 1.5);
1676   TH1I *errorTotalHtD_ =
1677       errorDir.make<TH1I>("errorTotalHtD", "errorTotalHtD;E_{T};Number of Events", 1024, -0.5, 1023.5);
1678   TH1I *errorTotalHtD_Of_ =
1679       errorDir.make<TH1I>("errorTotalHtD_Of", "errorTotalHtD_Of;Overflow Bit Status;Number of Events", 2, -0.5, 1.5);
1680   TH1I *errorTotalHtE_ =
1681       errorDir.make<TH1I>("errorTotalHtE", "errorTotalHtE;E_{T};Number of Events", 1024, -0.5, 1023.5);
1682   TH1I *errorTotalHtE_Of_ =
1683       errorDir.make<TH1I>("errorTotalHtE_Of", "errorTotalHtE_Of;Overflow Bit Status;Number of Events", 2, -0.5, 1.5);
1684 
1685   //plot the region ET and OF bits
1686   for (unsigned int i = 0; i < caloRegions->size(); i++) {
1687     if (caloRegions->at(i).bx() == RCTTrigBx_) {
1688       if (caloRegions->at(i).et() > 0)
1689         errorRegionEtEtaPhi_->Fill(caloRegions->at(i).gctEta(), caloRegions->at(i).gctPhi(), caloRegions->at(i).et());
1690       if (caloRegions->at(i).overFlow())
1691         errorRegionOfEtaPhi_->Fill(caloRegions->at(i).gctEta(), caloRegions->at(i).gctPhi());
1692     }
1693   }
1694   //now plot the error ET
1695   for (unsigned int i = 0; i < totalEtD->size(); i++) {
1696     if (totalEtD->at(i).bx() == GCTTrigBx_) {
1697       errorTotalEtD_Of_->Fill(totalEtD->at(i).overFlow());
1698       if (!totalEtD->at(i).overFlow())
1699         errorTotalEtD_->Fill(totalEtD->at(i).et());
1700     }
1701   }
1702   for (unsigned int i = 0; i < totalEtE->size(); i++) {
1703     if (totalEtE->at(i).bx() == EmuTrigBx_) {
1704       errorTotalEtE_Of_->Fill(totalEtE->at(i).overFlow());
1705       if (!totalEtE->at(i).overFlow())
1706         errorTotalEtE_->Fill(totalEtE->at(i).et());
1707     }
1708   }
1709   //and now the error HT
1710   for (unsigned int i = 0; i < totalHtD->size(); i++) {
1711     if (totalHtD->at(i).bx() == GCTTrigBx_) {
1712       errorTotalHtD_Of_->Fill(totalHtD->at(i).overFlow());
1713       if (!totalHtD->at(i).overFlow())
1714         errorTotalHtD_->Fill(totalHtD->at(i).et());
1715     }
1716   }
1717   for (unsigned int i = 0; i < totalHtE->size(); i++) {
1718     if (totalHtE->at(i).bx() == EmuTrigBx_) {
1719       errorTotalHtE_Of_->Fill(totalHtE->at(i).overFlow());
1720       if (!totalHtE->at(i).overFlow())
1721         errorTotalHtE_->Fill(totalHtE->at(i).et());
1722     }
1723   }
1724 }
1725 
1726 void GctErrorAnalyzer::plotMissingEt(const edm::Handle<L1GctEtMissCollection> &missingEtD,
1727                                      const edm::Handle<L1GctEtMissCollection> &missingEtE) {
1728   for (unsigned int i = 0; i < missingEtD->size(); i++) {
1729     if (doGCTMBx_ || missingEtD->at(i).bx() == GCTTrigBx_) {
1730       missingEtD_Of_->Fill(missingEtD->at(i).overFlow());
1731       if (!missingEtD->at(i).overFlow() && missingEtD->at(i).et() > 0) {
1732         missingEtD_->Fill(missingEtD->at(i).et());
1733         missingEtD_Phi_->Fill(missingEtD->at(i).phi());
1734       }
1735     }
1736   }
1737 
1738   for (unsigned int i = 0; i < missingEtE->size(); i++) {
1739     if (doEmuMBx_ || missingEtE->at(i).bx() == EmuTrigBx_) {
1740       missingEtE_Of_->Fill(missingEtE->at(i).overFlow());
1741       if (!missingEtE->at(i).overFlow() && missingEtE->at(i).et()) {
1742         missingEtE_->Fill(missingEtE->at(i).et());
1743         missingEtE_Phi_->Fill(missingEtE->at(i).phi());
1744       }
1745     }
1746   }
1747 }
1748 
1749 void GctErrorAnalyzer::plotMissingHt(const edm::Handle<L1GctHtMissCollection> &missingHtD,
1750                                      const edm::Handle<L1GctHtMissCollection> &missingHtE) {
1751   for (unsigned int i = 0; i < missingHtD->size(); i++) {
1752     if (doGCTMBx_ || missingHtD->at(i).bx() == GCTTrigBx_) {
1753       missingHtD_Of_->Fill(missingHtD->at(i).overFlow());
1754       if (!missingHtD->at(i).overFlow() && missingHtD->at(i).et() > 0) {
1755         missingHtD_->Fill(missingHtD->at(i).et());
1756         missingHtD_Phi_->Fill(missingHtD->at(i).phi());
1757       }
1758     }
1759   }
1760 
1761   for (unsigned int i = 0; i < missingHtE->size(); i++) {
1762     if (doEmuMBx_ || missingHtE->at(i).bx() == EmuTrigBx_) {
1763       missingHtE_Of_->Fill(missingHtE->at(i).overFlow());
1764       if (!missingHtE->at(i).overFlow() && missingHtE->at(i).et() > 0) {
1765         missingHtE_->Fill(missingHtE->at(i).et());
1766         missingHtE_Phi_->Fill(missingHtE->at(i).phi());
1767       }
1768     }
1769   }
1770 }
1771 
1772 void GctErrorAnalyzer::plotMissingEErrors(const edm::Handle<L1GctEtMissCollection> &missingEtD,
1773                                           const edm::Handle<L1GctEtMissCollection> &missingEtE,
1774                                           const edm::Handle<L1GctHtMissCollection> &missingHtD,
1775                                           const edm::Handle<L1GctHtMissCollection> &missingHtE,
1776                                           edm::Handle<L1CaloRegionCollection> &caloRegions,
1777                                           const edm::Handle<L1GctInternJetDataCollection> &intJetsE,
1778                                           const edm::Handle<L1GctInternHtMissCollection> intMissingHtD) {
1779   std::string errorDirName = "err_";
1780   if (isMissingEError)
1781     errorDirName.append("E");
1782   if (isMissingHError)
1783     errorDirName.append("H");
1784 
1785   //added 05.03.2010 to highlight overflow errors in the missing energy sum calculation
1786   for (unsigned int i = 0; i < missingHtE->size(); i++) {
1787     if (missingHtE->at(i).bx() == EmuTrigBx_) {
1788       for (unsigned int j = 0; j < missingHtD->size(); j++) {
1789         if (missingHtD->at(j).bx() == GCTTrigBx_) {
1790           if (missingHtD->at(j).overFlow() != missingHtE->at(i).overFlow())
1791             errorDirName.append("O");
1792         }
1793       }
1794     }
1795   }
1796 
1797   std::stringstream caseNumber;
1798   caseNumber << eventNumber;
1799   errorDirName.append(caseNumber.str());
1800   TFileDirectory errorDir = errorHistCat.at(4).mkdir(errorDirName);
1801 
1802   TH2I *errorRegionEtEtaPhi_ = errorDir.make<TH2I>(
1803       "errorRegionEtEtaPhi", "errorRegionEtEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
1804   TH2I *errorRegionOfEtaPhi_ = errorDir.make<TH2I>(
1805       "errorRegionOfEtaPhi", "errorRegionOfEtaPhi;#eta (GCT Units);#phi (GCT Units)", 22, -0.5, 21.5, 18, -0.5, 17.5);
1806   TH1I *errorMissingEtD_ =
1807       errorDir.make<TH1I>("errorMissingEtD", "errorMissingEtD;E_{T};Number of Events", 1024, -0.5, 1023.5);
1808   TH1I *errorMissingEtD_Of_ = errorDir.make<TH1I>(
1809       "errorMissingEtD_Of", "errorMissingEtD_Of;Overflow Bit Status;Number of Events", 2, -0.5, 1.5);
1810   TH1I *errorMissingEtD_Phi_ = errorDir.make<TH1I>(
1811       "errorMissingEtD_Phi", "errorMissingEtD_Phi;Missing E_{T} #phi;Number of Events", 72, -0.5, 71.5);
1812   TH1I *errorMissingEtE_ =
1813       errorDir.make<TH1I>("errorMissingEtE", "errorMissingEtE;E_{T};Number of Events", 1024, -0.5, 1023.5);
1814   TH1I *errorMissingEtE_Of_ = errorDir.make<TH1I>(
1815       "errorMissingEtE_Of", "errorMissingEtE_Of;Overflow Bit Status;Number of Events", 2, -0.5, 1.5);
1816   TH1I *errorMissingEtE_Phi_ = errorDir.make<TH1I>(
1817       "errorMissingEtE_Phi", "errorMissingEtE_Phi;Missing E_{T} #phi;Number of Events", 72, -0.5, 71.5);
1818   TH1I *errorMissingHtD_ =
1819       errorDir.make<TH1I>("errorMissingHtD", "errorMissingHtD;H_{T};Number of Events", 1024, -0.5, 1023.5);
1820   TH1I *errorMissingHtD_Of_ = errorDir.make<TH1I>(
1821       "errorMissingHtD_Of", "errorMissingHtD_Of;Overflow Bit Status;Number of Events", 2, -0.5, 1.5);
1822   TH1I *errorMissingHtD_Phi_ = errorDir.make<TH1I>(
1823       "errorMissingHtD_Phi", "errorMissingHtD_Phi;Missing H_{T} #phi;Number of Events", 72, -0.5, 71.5);
1824   TH1I *errorMissingHtE_ =
1825       errorDir.make<TH1I>("errorMissingHtE", "errorMissingHtE;H_{T};Number of Events", 1024, -0.5, 1023.5);
1826   TH1I *errorMissingHtE_Of_ = errorDir.make<TH1I>(
1827       "errorMissingHtE_Of", "errorMissingHtE_Of;Overflow Bit Status;Number of Events", 2, -0.5, 1.5);
1828   TH1I *errorMissingHtE_Phi_ = errorDir.make<TH1I>(
1829       "errorMissingHtE_Phi", "errorMissingHtE_Phi;Missing H_{T} #phi;Number of Events", 72, -0.5, 71.5);
1830 
1831   //Added 19.03.2010 to provide additional information in the case of missingHt failures
1832   //1. The MHT from both wheel inputs (i.e. the leaf cards)
1833   //2. The emulator jet et,eta,phi for all jets found in an event
1834   if (doExtraMissingHTDebug_) {
1835     if (checkCollections(intMissingHtD, GCT_INT_HTMISS_QUANTA, "Internal Missing Ht Data")) {
1836       TH1I *errorMissingHtD_HtXPosLeaf1 = errorDir.make<TH1I>(
1837           "errorMissingHtD_HtXPosLeaf1", "missingHtD;Missing H_{T} X PosLeaf1;Number of Events", 4096, -2048.5, 2047.5);
1838       TH1I *errorMissingHtD_HtXPosLeaf2 = errorDir.make<TH1I>(
1839           "errorMissingHtD_HtXPosLeaf2", "missingHtD;Missing H_{T} X PosLeaf2;Number of Events", 4096, -2048.5, 2047.5);
1840       TH1I *errorMissingHtD_HtXPosLeaf3 = errorDir.make<TH1I>(
1841           "errorMissingHtD_HtXPosLeaf3", "missingHtD;Missing H_{T} X PosLeaf3;Number of Events", 4096, -2048.5, 2047.5);
1842       TH1I *errorMissingHtD_HtXNegLeaf1 = errorDir.make<TH1I>(
1843           "errorMissingHtD_HtXNegLeaf1", "missingHtD;Missing H_{T} X NegLeaf1;Number of Events", 4096, -2048.5, 2047.5);
1844       TH1I *errorMissingHtD_HtXNegLeaf2 = errorDir.make<TH1I>(
1845           "errorMissingHtD_HtXNegLeaf2", "missingHtD;Missing H_{T} X NegLeaf2;Number of Events", 4096, -2048.5, 2047.5);
1846       TH1I *errorMissingHtD_HtXNegLeaf3 = errorDir.make<TH1I>(
1847           "errorMissingHtD_HtXNegLeaf3", "missingHtD;Missing H_{T} X NegLeaf3;Number of Events", 4096, -2048.5, 2047.5);
1848 
1849       TH1I *errorMissingHtD_HtYPosLeaf1 = errorDir.make<TH1I>(
1850           "errorMissingHtD_HtYPosLeaf1", "missingHtD;Missing H_{T} Y PosLeaf1;Number of Events", 4096, -2048.5, 2047.5);
1851       TH1I *errorMissingHtD_HtYPosLeaf2 = errorDir.make<TH1I>(
1852           "errorMissingHtD_HtYPosLeaf2", "missingHtD;Missing H_{T} Y PosLeaf2;Number of Events", 4096, -2048.5, 2047.5);
1853       TH1I *errorMissingHtD_HtYPosLeaf3 = errorDir.make<TH1I>(
1854           "errorMissingHtD_HtYPosLeaf3", "missingHtD;Missing H_{T} Y PosLeaf3;Number of Events", 4096, -2048.5, 2047.5);
1855       TH1I *errorMissingHtD_HtYNegLeaf1 = errorDir.make<TH1I>(
1856           "errorMissingHtD_HtYNegLeaf1", "missingHtD;Missing H_{T} Y NegLeaf1;Number of Events", 4096, -2048.5, 2047.5);
1857       TH1I *errorMissingHtD_HtYNegLeaf2 = errorDir.make<TH1I>(
1858           "errorMissingHtD_HtYNegLeaf2", "missingHtD;Missing H_{T} Y NegLeaf2;Number of Events", 4096, -2048.5, 2047.5);
1859       TH1I *errorMissingHtD_HtYNegLeaf3 = errorDir.make<TH1I>(
1860           "errorMissingHtD_HtYNegLeaf3", "missingHtD;Missing H_{T} Y NegLeaf3;Number of Events", 4096, -2048.5, 2047.5);
1861 
1862       for (unsigned int i = 0; i < intMissingHtD->size(); i++) {
1863         if (intMissingHtD->at(i).bx() == GCTTrigBx_) {
1864           if (!intMissingHtD->at(i).overflow()) {
1865             if (intMissingHtD->at(i).capBlock() == 0x301 && intMissingHtD->at(i).capIndex() == 0 &&
1866                 intMissingHtD->at(i).isThereHtx())
1867               errorMissingHtD_HtXPosLeaf1->Fill(intMissingHtD->at(i).htx());
1868             if (intMissingHtD->at(i).capBlock() == 0x301 && intMissingHtD->at(i).capIndex() == 1 &&
1869                 intMissingHtD->at(i).isThereHtx())
1870               errorMissingHtD_HtXPosLeaf2->Fill(intMissingHtD->at(i).htx());
1871             if (intMissingHtD->at(i).capBlock() == 0x301 && intMissingHtD->at(i).capIndex() == 2 &&
1872                 intMissingHtD->at(i).isThereHtx())
1873               errorMissingHtD_HtXPosLeaf3->Fill(intMissingHtD->at(i).htx());
1874             if (intMissingHtD->at(i).capBlock() == 0x701 && intMissingHtD->at(i).capIndex() == 0 &&
1875                 intMissingHtD->at(i).isThereHtx())
1876               errorMissingHtD_HtXNegLeaf1->Fill(intMissingHtD->at(i).htx());
1877             if (intMissingHtD->at(i).capBlock() == 0x701 && intMissingHtD->at(i).capIndex() == 1 &&
1878                 intMissingHtD->at(i).isThereHtx())
1879               errorMissingHtD_HtXNegLeaf2->Fill(intMissingHtD->at(i).htx());
1880             if (intMissingHtD->at(i).capBlock() == 0x701 && intMissingHtD->at(i).capIndex() == 2 &&
1881                 intMissingHtD->at(i).isThereHtx())
1882               errorMissingHtD_HtXNegLeaf3->Fill(intMissingHtD->at(i).htx());
1883 
1884             if (intMissingHtD->at(i).capBlock() == 0x301 && intMissingHtD->at(i).capIndex() == 0 &&
1885                 intMissingHtD->at(i).isThereHty())
1886               errorMissingHtD_HtYPosLeaf1->Fill(intMissingHtD->at(i).hty());
1887             if (intMissingHtD->at(i).capBlock() == 0x301 && intMissingHtD->at(i).capIndex() == 1 &&
1888                 intMissingHtD->at(i).isThereHty())
1889               errorMissingHtD_HtYPosLeaf2->Fill(intMissingHtD->at(i).hty());
1890             if (intMissingHtD->at(i).capBlock() == 0x301 && intMissingHtD->at(i).capIndex() == 2 &&
1891                 intMissingHtD->at(i).isThereHty())
1892               errorMissingHtD_HtYPosLeaf3->Fill(intMissingHtD->at(i).hty());
1893             if (intMissingHtD->at(i).capBlock() == 0x701 && intMissingHtD->at(i).capIndex() == 0 &&
1894                 intMissingHtD->at(i).isThereHty())
1895               errorMissingHtD_HtYNegLeaf1->Fill(intMissingHtD->at(i).hty());
1896             if (intMissingHtD->at(i).capBlock() == 0x701 && intMissingHtD->at(i).capIndex() == 1 &&
1897                 intMissingHtD->at(i).isThereHty())
1898               errorMissingHtD_HtYNegLeaf2->Fill(intMissingHtD->at(i).hty());
1899             if (intMissingHtD->at(i).capBlock() == 0x701 && intMissingHtD->at(i).capIndex() == 2 &&
1900                 intMissingHtD->at(i).isThereHty())
1901               errorMissingHtD_HtYNegLeaf3->Fill(intMissingHtD->at(i).hty());
1902           }
1903         }
1904       }
1905     }
1906   }
1907 
1908   if (checkCollections(intJetsE, NUM_INT_JETS, "Intermediate Jets Emulator")) {
1909     TH2I *errorIntJetsE_EtEtaPhi = errorDir.make<TH2I>("errorIntJetsE_EtEtaPhi",
1910                                                        "errorIntJetsE_EtEtaPhi;#eta (GCT Units);#phi (GCT Units)",
1911                                                        22,
1912                                                        -0.5,
1913                                                        21.5,
1914                                                        18,
1915                                                        -0.5,
1916                                                        17.5);
1917 
1918     for (unsigned int i = 0; i < intJetsE->size(); i++) {
1919       if (intJetsE->at(i).bx() == EmuTrigBx_) {
1920         if (!intJetsE->at(i).oflow() && intJetsE->at(i).et())
1921           errorIntJetsE_EtEtaPhi->Fill(
1922               intJetsE->at(i).regionId().ieta(), intJetsE->at(i).regionId().iphi(), intJetsE->at(i).et());
1923       }
1924     }
1925   }
1926 
1927   for (unsigned int i = 0; i < caloRegions->size(); i++) {
1928     if (caloRegions->at(i).bx() == RCTTrigBx_) {
1929       if (caloRegions->at(i).et() > 0)
1930         errorRegionEtEtaPhi_->Fill(caloRegions->at(i).gctEta(), caloRegions->at(i).gctPhi(), caloRegions->at(i).et());
1931       if (caloRegions->at(i).overFlow())
1932         errorRegionOfEtaPhi_->Fill(caloRegions->at(i).gctEta(), caloRegions->at(i).gctPhi());
1933     }
1934   }
1935 
1936   //plot the data candidates relating to this mismatch
1937   for (unsigned int i = 0; i < missingEtD->size(); i++) {
1938     if (missingEtD->at(i).bx() == GCTTrigBx_) {
1939       errorMissingEtD_Of_->Fill(missingEtD->at(i).overFlow());
1940       if (!missingEtD->at(i).overFlow()) {
1941         errorMissingEtD_->Fill(missingEtD->at(i).et());
1942         errorMissingEtD_Phi_->Fill(missingEtD->at(i).phi());
1943       }
1944     }
1945   }
1946   for (unsigned int i = 0; i < missingHtD->size(); i++) {
1947     if (missingHtD->at(i).bx() == GCTTrigBx_) {
1948       errorMissingHtD_Of_->Fill(missingHtD->at(i).overFlow());
1949       if (!missingHtD->at(i).overFlow()) {
1950         errorMissingHtD_->Fill(missingHtD->at(i).et());
1951         errorMissingHtD_Phi_->Fill(missingHtD->at(i).phi());
1952       }
1953     }
1954   }
1955   //and now for the emulator candidates
1956   for (unsigned int i = 0; i < missingEtE->size(); i++) {
1957     if (missingEtE->at(i).bx() == EmuTrigBx_) {
1958       errorMissingEtE_Of_->Fill(missingEtE->at(i).overFlow());
1959       if (!missingEtE->at(i).overFlow()) {
1960         errorMissingEtE_->Fill(missingEtE->at(i).et());
1961         errorMissingEtE_Phi_->Fill(missingEtE->at(i).phi());
1962       }
1963     }
1964   }
1965   for (unsigned int i = 0; i < missingHtE->size(); i++) {
1966     if (missingHtE->at(i).bx() == EmuTrigBx_) {
1967       errorMissingHtE_Of_->Fill(missingHtE->at(i).overFlow());
1968       if (!missingHtE->at(i)
1969                .overFlow()) {  //to see what values the emulator outputs despite the o/f bit being set comment out this statement
1970         errorMissingHtE_->Fill(missingHtE->at(i).et());
1971         errorMissingHtE_Phi_->Fill(missingHtE->at(i).phi());
1972       }
1973     }
1974   }
1975 }
1976 
1977 //define this as a plug-in
1978 DEFINE_FWK_MODULE(GctErrorAnalyzer);