Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:25:50

0001 // -*- C++ -*-
0002 //
0003 // Package:    HcalSevLvlAnalyzer
0004 // Class:      HcalSevLvlAnalyzer
0005 //
0006 /**\class HcalSevLvlAnalyzer HcalSevLvlAnalyzer.cc 
0007 
0008  Description: <one line class summary>
0009 
0010  Implementation:
0011      <Notes on implementation>
0012 */
0013 //
0014 // Original Author:  Radek Ofierzynski
0015 //         Created:  Wed Jan 21 13:46:27 CET 2009
0016 //
0017 //
0018 
0019 // system include files
0020 #include <iostream>
0021 #include <memory>
0022 
0023 // user include files
0024 #include "FWCore/Framework/interface/Frameworkfwd.h"
0025 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0026 
0027 #include "FWCore/Framework/interface/Event.h"
0028 #include "FWCore/Framework/interface/MakerMacros.h"
0029 
0030 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0031 
0032 #include "DataFormats/HcalDetId/interface/HcalGenericDetId.h"
0033 #include "RecoLocalCalo/HcalRecAlgos/interface/HcalSeverityLevelComputer.h"
0034 #include "RecoLocalCalo/HcalRecAlgos/interface/HcalSeverityLevelComputerRcd.h"
0035 //
0036 // class decleration
0037 //
0038 
0039 class HcalSevLvlAnalyzer : public edm::one::EDAnalyzer<> {
0040 public:
0041   explicit HcalSevLvlAnalyzer(const edm::ParameterSet&);
0042   ~HcalSevLvlAnalyzer();
0043 
0044 private:
0045   void beginJob() override;
0046   void analyze(const edm::Event&, const edm::EventSetup&) override;
0047   void endJob() override;
0048 
0049   // ----------member data ---------------------------
0050   const edm::ESGetToken<HcalSeverityLevelComputer, HcalSeverityLevelComputerRcd> tokSev_;
0051 };
0052 
0053 //
0054 // constants, enums and typedefs
0055 //
0056 
0057 //
0058 // static data member definitions
0059 //
0060 
0061 //
0062 // constructors and destructor
0063 //
0064 HcalSevLvlAnalyzer::HcalSevLvlAnalyzer(const edm::ParameterSet& iConfig) : tokSev_(esConsumes()) {
0065   //now do what ever initialization is needed
0066 
0067   // initialize the severity level code
0068 }
0069 
0070 HcalSevLvlAnalyzer::~HcalSevLvlAnalyzer() {
0071   // do anything here that needs to be done at desctruction time
0072   // (e.g. close files, deallocate resources etc.)
0073 }
0074 
0075 //
0076 // member functions
0077 //
0078 
0079 // ------------ method called to for each event  ------------
0080 void HcalSevLvlAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0081   // here's how to access the severity level computer:
0082   const HcalSeverityLevelComputer* myProd = &iSetup.getData(tokSev_);
0083 
0084   // create some example cases:
0085   std::cout << std::showbase;
0086   //HB
0087   HcalGenericDetId myIdHB(1107313548);
0088   uint32_t sampleHBChSt[5] = {0,     //nothing
0089                               4,     //nothing significant
0090                               64,    //hot
0091                               1,     //off
0092                               32};   //dead
0093   uint32_t sampleHBRHFlag[4] = {0,   //nothing
0094                                 32,  //nothing significant
0095                                 2,   //Pulseshape
0096                                 1};  //Multiplicity
0097 
0098   for (unsigned i = 0; i < 4; i++)    //loop over rechitflag
0099     for (unsigned k = 0; k < 5; k++)  //loop over channelstatus
0100     {
0101       int theLevel = myProd->getSeverityLevel(myIdHB, sampleHBRHFlag[i], sampleHBChSt[k]);
0102 
0103       std::cout << "Status for " << myIdHB << " with RHFlag " << sampleHBRHFlag[i] << " (" << std::hex
0104                 << sampleHBRHFlag[i] << ") and ChStFlag " << std::dec << sampleHBChSt[k] << " (" << std::hex
0105                 << sampleHBChSt[k] << std::dec << ") is: " << theLevel;
0106       std::cout << std::endl;
0107 
0108       bool dropchannel = myProd->dropChannel(sampleHBChSt[k]);
0109       bool recovered = myProd->recoveredRecHit(myIdHB, sampleHBRHFlag[i]);
0110 
0111       std::cout << "DropChannel status for " << myIdHB << " with RHFlag " << sampleHBRHFlag[i] << " (" << std::hex
0112                 << sampleHBRHFlag[i] << ") and ChStFlag " << std::dec << sampleHBChSt[k] << " (" << std::hex
0113                 << sampleHBChSt[k] << std::dec << ") is: " << dropchannel << ", recovered status is: " << recovered
0114                 << std::endl;
0115       std::cout << std::endl;
0116     }
0117   std::cout << std::endl;
0118 
0119   //   sampleChSt = 64; //hot
0120   //   sampleRHFlag = 2; //HBHEPulseShape
0121   //   theLevel = myProd->getSeverityLevel(myId, sampleRHFlag, sampleChSt);
0122   //   std::cout << "Status for " << myId << " with RHFlag " <<  std::hex << sampleRHFlag
0123   //         << " and ChStFlag " << sampleChSt << std::dec << " is: " << theLevel << std::endl;
0124   //   std::cout << std::endl;
0125 
0126   //HE
0127   HcalGenericDetId myIdHE(1140869515);
0128   uint32_t sampleHEChSt[5] = {0,     //nothing
0129                               4,     //nothing significant
0130                               64,    //hot
0131                               1,     //off
0132                               32};   //dead
0133   uint32_t sampleHERHFlag[4] = {0,   //nothing
0134                                 32,  //nothing significant
0135                                 2,   //Pulseshape
0136                                 1};  //Multiplicity
0137 
0138   for (unsigned i = 0; i < 4; i++)    //loop over rechitflag
0139     for (unsigned k = 0; k < 5; k++)  //loop over channelstatus
0140     {
0141       int theLevel = myProd->getSeverityLevel(myIdHE, sampleHERHFlag[i], sampleHEChSt[k]);
0142 
0143       std::cout << "Status for " << myIdHE << " with RHFlag " << sampleHERHFlag[i] << " (" << std::hex
0144                 << sampleHERHFlag[i] << ") and ChStFlag " << std::dec << sampleHEChSt[k] << " (" << std::hex
0145                 << sampleHEChSt[k] << std::dec << ") is: " << theLevel;
0146       std::cout << std::endl;
0147     }
0148   std::cout << std::endl;
0149 
0150   //HF
0151   HcalGenericDetId myIdHF(1207979653);
0152   uint32_t sampleHFChSt[5] = {0,     //nothing
0153                               4,     //nothing significant
0154                               64,    //hot
0155                               1,     //off
0156                               32};   //dead
0157   uint32_t sampleHFRHFlag[4] = {0,   //nothing
0158                                 32,  //nothing significant
0159                                 1,   //HFDigiTime
0160                                 2};  //HFLongShort
0161 
0162   for (unsigned i = 0; i < 4; i++)    //loop over rechitflag
0163     for (unsigned k = 0; k < 5; k++)  //loop over channelstatus
0164     {
0165       int theLevel = myProd->getSeverityLevel(myIdHF, sampleHFRHFlag[i], sampleHFChSt[k]);
0166 
0167       std::cout << "Status for " << myIdHF << " with RHFlag " << sampleHFRHFlag[i] << " (" << std::hex
0168                 << sampleHFRHFlag[i] << ") and ChStFlag " << std::dec << sampleHFChSt[k] << " (" << std::hex
0169                 << sampleHFChSt[k] << std::dec << ") is: " << theLevel;
0170       std::cout << std::endl;
0171     }
0172   std::cout << std::endl;
0173 
0174   //HO
0175   HcalGenericDetId myIdHO(1174471682);
0176   uint32_t sampleHOChSt[5] = {0,     //nothing
0177                               4,     //nothing significant
0178                               64,    //hot
0179                               1,     //off
0180                               32};   //dead
0181   uint32_t sampleHORHFlag[3] = {0,   //nothing
0182                                 32,  //nothing significant
0183                                 1};  //HOBit
0184 
0185   for (unsigned i = 0; i < 3; i++)    //loop over rechitflag
0186     for (unsigned k = 0; k < 5; k++)  //loop over channelstatus
0187     {
0188       int theLevel = myProd->getSeverityLevel(myIdHO, sampleHORHFlag[i], sampleHOChSt[k]);
0189 
0190       std::cout << "Status for " << myIdHO << " with RHFlag " << sampleHORHFlag[i] << " (" << std::hex
0191                 << sampleHORHFlag[i] << ") and ChStFlag " << std::dec << sampleHOChSt[k] << " (" << std::hex
0192                 << sampleHOChSt[k] << std::dec << ") is: " << theLevel;
0193       std::cout << std::endl;
0194     }
0195   std::cout << std::endl;
0196 
0197   //ZDC
0198 
0199   //Calib
0200 
0201   std::cout << std::noshowbase;
0202 }
0203 
0204 // ------------ method called once each job just before starting event loop  ------------
0205 void HcalSevLvlAnalyzer::beginJob() {}
0206 
0207 // ------------ method called once each job just after ending the event loop  ------------
0208 void HcalSevLvlAnalyzer::endJob() {}
0209 
0210 //define this as a plug-in
0211 DEFINE_FWK_MODULE(HcalSevLvlAnalyzer);