Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:58:09

0001 // system include files
0002 #include <memory>
0003 
0004 // user include files
0005 #include "FWCore/Framework/interface/Frameworkfwd.h"
0006 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0007 #include "FWCore/Framework/interface/Event.h"
0008 #include "FWCore/Framework/interface/EventSetup.h"
0009 #include "FWCore/Framework/interface/MakerMacros.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 #include "FWCore/ParameterSet/interface/FileInPath.h"
0012 #include "FWCore/PluginManager/interface/ModuleDef.h"
0013 #include "FWCore/ServiceRegistry/interface/Service.h"
0014 
0015 //calo headers
0016 #include "CondFormats/HcalObjects/interface/HBHEDarkening.h"
0017 #include "CondFormats/DataRecord/interface/HBHEDarkeningRecord.h"
0018 #include "CalibCalorimetry/HcalAlgos/interface/HBHERecalibration.h"
0019 #include "Geometry/CaloTopology/interface/HcalTopology.h"
0020 #include "Geometry/Records/interface/HcalRecNumberingRecord.h"
0021 
0022 //STL headers
0023 #include <iostream>
0024 #include <iomanip>
0025 #include <vector>
0026 
0027 //
0028 // class declaration
0029 //
0030 
0031 class HBHEDarkeningAnalyzer : public edm::one::EDAnalyzer<edm::one::WatchRuns> {
0032 public:
0033   explicit HBHEDarkeningAnalyzer(const edm::ParameterSet&);
0034   ~HBHEDarkeningAnalyzer();
0035 
0036   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0037 
0038 private:
0039   void beginJob() override;
0040   void beginRun(const edm::Run&, const edm::EventSetup&) override;
0041   void analyze(const edm::Event&, const edm::EventSetup&) override;
0042   void endRun(const edm::Run&, const edm::EventSetup&) override {}
0043   void endJob() override {}
0044   void print(int ieta_min,
0045              int ieta_max,
0046              int lay_min,
0047              int lay_max,
0048              const HBHEDarkening* darkening,
0049              const HBHERecalibration& recalibration);
0050 
0051   // ----------member data ---------------------------
0052   double intlumi;
0053   edm::ESGetToken<HcalTopology, HcalRecNumberingRecord> tok_topo_;
0054   edm::ESGetToken<HBHEDarkening, HBHEDarkeningRecord> tok_hbdark_;
0055   edm::ESGetToken<HBHEDarkening, HBHEDarkeningRecord> tok_hedark_;
0056 
0057   const HBHEDarkening* hb_darkening;
0058   const HBHEDarkening* he_darkening;
0059   HBHERecalibration hb_recalibration, he_recalibration;
0060   const HcalTopology* theTopology;
0061 };
0062 
0063 //
0064 // constructors and destructor
0065 //
0066 HBHEDarkeningAnalyzer::HBHEDarkeningAnalyzer(const edm::ParameterSet& iConfig)
0067     : intlumi(iConfig.getParameter<double>("deliveredLumi")),
0068       hb_darkening(NULL),
0069       he_darkening(NULL),
0070       hb_recalibration(intlumi, 0, iConfig.getParameter<edm::FileInPath>("HBmeanenergies").fullPath()),
0071       he_recalibration(intlumi, 0, iConfig.getParameter<edm::FileInPath>("HEmeanenergies").fullPath()),
0072       theTopology(NULL) {
0073   tok_topo_ = esConsumes<HcalTopology, HcalRecNumberingRecord, edm::Transition::BeginRun>();
0074   tok_hbdark_ = esConsumes<HBHEDarkening, HBHEDarkeningRecord, edm::Transition::BeginRun>(edm::ESInputTag("", "HB"));
0075   tok_hedark_ = esConsumes<HBHEDarkening, HBHEDarkeningRecord, edm::Transition::BeginRun>(edm::ESInputTag("", "HE"));
0076 }
0077 
0078 HBHEDarkeningAnalyzer::~HBHEDarkeningAnalyzer() {}
0079 
0080 //
0081 // member functions
0082 //
0083 
0084 void HBHEDarkeningAnalyzer::beginJob() {}
0085 
0086 void HBHEDarkeningAnalyzer::beginRun(const edm::Run& iRun, const edm::EventSetup& iSetup) {
0087   theTopology = &iSetup.getData(tok_topo_);
0088   hb_darkening = &iSetup.getData(tok_hbdark_);
0089   he_darkening = &iSetup.getData(tok_hedark_);
0090 
0091   //initialize recalibration classes
0092   std::vector<std::vector<int>> m_segmentation;
0093   int maxEta = theTopology->lastHBHERing();
0094   m_segmentation.resize(maxEta);
0095   for (int i = 0; i < maxEta; i++) {
0096     theTopology->getDepthSegmentation(i + 1, m_segmentation[i]);
0097   }
0098   std::cout << "HB: Eta " << theTopology->firstHBRing() << ":" << theTopology->lastHBRing() << " HE: Eta "
0099             << theTopology->firstHERing() << ":" << theTopology->lastHERing() << std::endl;
0100   hb_recalibration.setup(m_segmentation, hb_darkening);
0101   he_recalibration.setup(m_segmentation, he_darkening);
0102 }
0103 
0104 // ------------ method called on each new Event  ------------
0105 void HBHEDarkeningAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0106   std::cout << std::setprecision(2);
0107 
0108   //HB tests
0109   int hb_ieta_min = theTopology->firstHBRing();
0110   int hb_ieta_max = theTopology->lastHBRing();
0111   int hb_lay_min = 1;
0112   int hb_lay_max = 17;
0113   std::cout << "HB" << std::endl;
0114   print(hb_ieta_min, hb_ieta_max, hb_lay_min, hb_lay_max, hb_darkening, hb_recalibration);
0115 
0116   std::cout << std::endl;
0117 
0118   //HE tests
0119   int he_ieta_min = theTopology->firstHERing();
0120   int he_ieta_max = theTopology->lastHERing();
0121   int he_lay_min = 1;
0122   int he_lay_max = 19;
0123   std::cout << "HE" << std::endl;
0124   print(he_ieta_min, he_ieta_max, he_lay_min, he_lay_max, he_darkening, he_recalibration);
0125 }
0126 
0127 void HBHEDarkeningAnalyzer::print(int ieta_min,
0128                                   int ieta_max,
0129                                   int lay_min,
0130                                   int lay_max,
0131                                   const HBHEDarkening* darkening,
0132                                   const HBHERecalibration& recalibration) {
0133   std::cout << "Darkening: ieta " << ieta_min << ":" << ieta_max << " layer " << lay_min << ":" << lay_max << std::endl;
0134   for (int ieta = ieta_min; ieta <= ieta_max; ++ieta) {
0135     std::cout << "Tower " << ieta << ": ";
0136     for (int lay = lay_min; lay <= lay_max; ++lay) {
0137       std::cout << darkening->degradation(intlumi, ieta, lay) << " ";
0138     }
0139     std::cout << std::endl;
0140   }
0141 
0142   std::cout << "Recalibration: ieta " << ieta_min << ":" << ieta_max << " layer " << lay_min << ":" << lay_max
0143             << std::endl;
0144   for (int ieta = ieta_min; ieta <= ieta_max; ++ieta) {
0145     std::cout << "Tower " << ieta << ": ";
0146     for (int depth = 1; depth <= recalibration.maxDepth(); ++depth) {
0147       std::cout << recalibration.getCorr(ieta, depth) << " ";
0148     }
0149     std::cout << std::endl;
0150   }
0151 }
0152 
0153 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
0154 void HBHEDarkeningAnalyzer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0155   //The following says we do not know what parameters are allowed so do no validation
0156   // Please change this to state exactly what you do use, even if it is no parameters
0157   edm::ParameterSetDescription desc;
0158   desc.setUnknown();
0159   descriptions.addDefault(desc);
0160 }
0161 //define this as a plug-in
0162 DEFINE_FWK_MODULE(HBHEDarkeningAnalyzer);