HBHEDarkeningAnalyzer

Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162
// system include files
#include <memory>

// user include files
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/one/EDAnalyzer.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/ParameterSet/interface/FileInPath.h"
#include "FWCore/PluginManager/interface/ModuleDef.h"
#include "FWCore/ServiceRegistry/interface/Service.h"

//calo headers
#include "CondFormats/HcalObjects/interface/HBHEDarkening.h"
#include "CondFormats/DataRecord/interface/HBHEDarkeningRecord.h"
#include "CalibCalorimetry/HcalAlgos/interface/HBHERecalibration.h"
#include "Geometry/CaloTopology/interface/HcalTopology.h"
#include "Geometry/Records/interface/HcalRecNumberingRecord.h"

//STL headers
#include <iostream>
#include <iomanip>
#include <vector>

//
// class declaration
//

class HBHEDarkeningAnalyzer : public edm::one::EDAnalyzer<edm::one::WatchRuns> {
public:
  explicit HBHEDarkeningAnalyzer(const edm::ParameterSet&);
  ~HBHEDarkeningAnalyzer();

  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);

private:
  void beginJob() override;
  void beginRun(const edm::Run&, const edm::EventSetup&) override;
  void analyze(const edm::Event&, const edm::EventSetup&) override;
  void endRun(const edm::Run&, const edm::EventSetup&) override {}
  void endJob() override {}
  void print(int ieta_min,
             int ieta_max,
             int lay_min,
             int lay_max,
             const HBHEDarkening* darkening,
             const HBHERecalibration& recalibration);

  // ----------member data ---------------------------
  double intlumi;
  edm::ESGetToken<HcalTopology, HcalRecNumberingRecord> tok_topo_;
  edm::ESGetToken<HBHEDarkening, HBHEDarkeningRecord> tok_hbdark_;
  edm::ESGetToken<HBHEDarkening, HBHEDarkeningRecord> tok_hedark_;

  const HBHEDarkening* hb_darkening;
  const HBHEDarkening* he_darkening;
  HBHERecalibration hb_recalibration, he_recalibration;
  const HcalTopology* theTopology;
};

//
// constructors and destructor
//
HBHEDarkeningAnalyzer::HBHEDarkeningAnalyzer(const edm::ParameterSet& iConfig)
    : intlumi(iConfig.getParameter<double>("deliveredLumi")),
      hb_darkening(NULL),
      he_darkening(NULL),
      hb_recalibration(intlumi, 0, iConfig.getParameter<edm::FileInPath>("HBmeanenergies").fullPath()),
      he_recalibration(intlumi, 0, iConfig.getParameter<edm::FileInPath>("HEmeanenergies").fullPath()),
      theTopology(NULL) {
  tok_topo_ = esConsumes<HcalTopology, HcalRecNumberingRecord, edm::Transition::BeginRun>();
  tok_hbdark_ = esConsumes<HBHEDarkening, HBHEDarkeningRecord, edm::Transition::BeginRun>(edm::ESInputTag("", "HB"));
  tok_hedark_ = esConsumes<HBHEDarkening, HBHEDarkeningRecord, edm::Transition::BeginRun>(edm::ESInputTag("", "HE"));
}

HBHEDarkeningAnalyzer::~HBHEDarkeningAnalyzer() {}

//
// member functions
//

void HBHEDarkeningAnalyzer::beginJob() {}

void HBHEDarkeningAnalyzer::beginRun(const edm::Run& iRun, const edm::EventSetup& iSetup) {
  theTopology = &iSetup.getData(tok_topo_);
  hb_darkening = &iSetup.getData(tok_hbdark_);
  he_darkening = &iSetup.getData(tok_hedark_);

  //initialize recalibration classes
  std::vector<std::vector<int>> m_segmentation;
  int maxEta = theTopology->lastHBHERing();
  m_segmentation.resize(maxEta);
  for (int i = 0; i < maxEta; i++) {
    theTopology->getDepthSegmentation(i + 1, m_segmentation[i]);
  }
  std::cout << "HB: Eta " << theTopology->firstHBRing() << ":" << theTopology->lastHBRing() << " HE: Eta "
            << theTopology->firstHERing() << ":" << theTopology->lastHERing() << std::endl;
  hb_recalibration.setup(m_segmentation, hb_darkening);
  he_recalibration.setup(m_segmentation, he_darkening);
}

// ------------ method called on each new Event  ------------
void HBHEDarkeningAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
  std::cout << std::setprecision(2);

  //HB tests
  int hb_ieta_min = theTopology->firstHBRing();
  int hb_ieta_max = theTopology->lastHBRing();
  int hb_lay_min = 1;
  int hb_lay_max = 17;
  std::cout << "HB" << std::endl;
  print(hb_ieta_min, hb_ieta_max, hb_lay_min, hb_lay_max, hb_darkening, hb_recalibration);

  std::cout << std::endl;

  //HE tests
  int he_ieta_min = theTopology->firstHERing();
  int he_ieta_max = theTopology->lastHERing();
  int he_lay_min = 1;
  int he_lay_max = 19;
  std::cout << "HE" << std::endl;
  print(he_ieta_min, he_ieta_max, he_lay_min, he_lay_max, he_darkening, he_recalibration);
}

void HBHEDarkeningAnalyzer::print(int ieta_min,
                                  int ieta_max,
                                  int lay_min,
                                  int lay_max,
                                  const HBHEDarkening* darkening,
                                  const HBHERecalibration& recalibration) {
  std::cout << "Darkening: ieta " << ieta_min << ":" << ieta_max << " layer " << lay_min << ":" << lay_max << std::endl;
  for (int ieta = ieta_min; ieta <= ieta_max; ++ieta) {
    std::cout << "Tower " << ieta << ": ";
    for (int lay = lay_min; lay <= lay_max; ++lay) {
      std::cout << darkening->degradation(intlumi, ieta, lay) << " ";
    }
    std::cout << std::endl;
  }

  std::cout << "Recalibration: ieta " << ieta_min << ":" << ieta_max << " layer " << lay_min << ":" << lay_max
            << std::endl;
  for (int ieta = ieta_min; ieta <= ieta_max; ++ieta) {
    std::cout << "Tower " << ieta << ": ";
    for (int depth = 1; depth <= recalibration.maxDepth(); ++depth) {
      std::cout << recalibration.getCorr(ieta, depth) << " ";
    }
    std::cout << std::endl;
  }
}

// ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
void HBHEDarkeningAnalyzer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
  //The following says we do not know what parameters are allowed so do no validation
  // Please change this to state exactly what you do use, even if it is no parameters
  edm::ParameterSetDescription desc;
  desc.setUnknown();
  descriptions.addDefault(desc);
}
//define this as a plug-in
DEFINE_FWK_MODULE(HBHEDarkeningAnalyzer);