CMS 3D CMS Logo

HcalCollapseAnalyzer.cc
Go to the documentation of this file.
1 // -*- C++ -*-//
2 // Package: Hcal
3 // Class: HcalCollapseAnalyzer
4 //
13 //
14 // Original Author: Sunanda Banerjee
15 // Created: Thu Dec 26 18:52:02 CST 2017
16 //
17 //
18 
19 // system include files
20 #include <string>
21 #include <vector>
22 
23 // Root objects
24 #include "TH1.h"
25 
26 // user include files
29 
32 
40 
43 
45 public:
46  explicit HcalCollapseAnalyzer(const edm::ParameterSet &);
47  ~HcalCollapseAnalyzer() override {}
48 
49  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
50 
51 private:
52  void analyze(edm::Event const &, edm::EventSetup const &) override;
53  void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
54 
55  // ----------member data ---------------------------
57  const int verbosity_;
59  const bool doHE_, doHB_;
60  const HcalTopology *theHBHETopology = nullptr;
61 
65 
68 };
69 
71  : topFolderName_(iConfig.getParameter<std::string>("topFolderName")),
72  verbosity_(iConfig.getUntrackedParameter<int>("verbosity", 0)),
73  recHitHBHE_(iConfig.getUntrackedParameter<edm::InputTag>("recHitHBHE", edm::InputTag("hbhereco"))),
74  preRecHitHBHE_(iConfig.getUntrackedParameter<edm::InputTag>("preRecHitHBHE", edm::InputTag("hbheprereco"))),
75  doHE_(iConfig.getUntrackedParameter<bool>("doHE", true)),
76  doHB_(iConfig.getUntrackedParameter<bool>("doHB", false)),
77  hcalTopologyToken_{esConsumes<HcalTopology, HcalRecNumberingRecord>()} {
78  // define tokens for access
79  tok_hbhe_ = consumes<HBHERecHitCollection>(recHitHBHE_);
80  tok_prehbhe_ = consumes<HBHERecHitCollection>(preRecHitHBHE_);
81 
82  edm::LogVerbatim("Collapse") << "Verbosity " << verbosity_ << " with tags " << recHitHBHE_ << " and "
83  << preRecHitHBHE_ << " and Do " << doHB_ << ":" << doHE_;
84 }
85 
88  desc.add<std::string>("topFolderName", "HcalCollapse");
89  desc.addUntracked<int>("verbosity", 0);
90  desc.addUntracked<edm::InputTag>("recHitHBHE", edm::InputTag("hbhereco"));
91  desc.addUntracked<edm::InputTag>("preRecHitHBHE", edm::InputTag("hbheprereco"));
92  desc.addUntracked<bool>("doHE", true);
93  desc.addUntracked<bool>("doHB", false);
94  descriptions.add("hcalCollapseAnalyzer", desc);
95 }
96 
98  if (verbosity_ > 0)
99  edm::LogVerbatim("Collapse") << "Run " << iEvent.id().run() << " Event " << iEvent.id().event()
100  << " starts ==============";
101 
103 
105  iEvent.getByToken(tok_hbhe_, hbhereco);
107  iEvent.getByToken(tok_prehbhe_, hbheprereco);
108  if (verbosity_ > 0) {
109  edm::LogVerbatim("Collapse") << "Handle Reco " << hbhereco << " Size " << hbhereco->size();
110  edm::LogVerbatim("Collapse") << "Handle PreReco " << hbheprereco << " Size " << hbheprereco->size();
111  }
112  if (hbhereco.isValid() && hbheprereco.isValid()) {
113  const HBHERecHitCollection *recohbhe = hbhereco.product();
114  const HBHERecHitCollection *prerecohbhe = hbheprereco.product();
115  if (verbosity_ > 0)
116  edm::LogVerbatim("Collapse") << "Size of hbhereco: " << recohbhe->size()
117  << " and hbheprereco: " << prerecohbhe->size();
118  double sfrac = (prerecohbhe->empty()) ? 1 : ((double)(recohbhe->size())) / ((double)(prerecohbhe->size()));
119  h_sfrac->Fill(sfrac);
120  h_size->Fill(recohbhe->size());
121  for (const auto &hit : (*recohbhe)) {
122  HcalDetId id = hit.id();
123  if (((id.subdet() == HcalEndcap) && doHE_) || ((id.subdet() == HcalBarrel) && doHB_)) {
124  h_depth->Fill(id.depth());
125  std::vector<HcalDetId> ids;
127  if (verbosity_ > 0) {
128  edm::LogVerbatim("Collapse") << id << " is derived from " << ids.size() << " components";
129  for (unsigned int k = 0; k < ids.size(); ++k)
130  edm::LogVerbatim("Collapse") << "[" << k << "] " << ids[k];
131  }
132  h_merge->Fill(ids.size());
133  double energy = hit.energy();
134  double etot(0);
135  unsigned int k(0);
136  for (const auto phit : (*prerecohbhe)) {
137  if (std::find(ids.begin(), ids.end(), phit.id()) != ids.end()) {
138  etot += phit.energy();
139  double frac = (energy > 0) ? phit.energy() / energy : 1.0;
140  h_frac->Fill(frac);
141  if (verbosity_ > 0)
142  edm::LogVerbatim("Collapse") << "Frac[" << k << "] " << frac << " for " << phit.id();
143  ++k;
144  }
145  }
146  double frac = (energy > 0) ? etot / energy : 1.0;
147  h_balance->Fill(frac);
148  if (verbosity_ > 0)
149  edm::LogVerbatim("Collapse") << "Total energy " << energy << ":" << etot << ":" << frac;
150  }
151  }
152  }
153 }
154 
156  // Book histograms
158  h_merge = ibooker.book1D("h_merge", "Number of hits merged", 10, 0.0, 10.0);
159  h_size = ibooker.book1D("h_size", "Size of the RecHit collection", 100, 500.0, 1500.0);
160  h_depth = ibooker.book1D("h_depth", "Depth of the Id's used", 10, 0.0, 10.0);
161  h_sfrac = ibooker.book1D("h_sfrac", "Ratio of sizes of preRecHit and RecHit collections", 150, 0.0, 1.5);
162  h_frac = ibooker.book1D("h_frac", "Fraction of energy before collapse", 150, 0.0, 1.5);
163  h_balance = ibooker.book1D("h_balance", "Balance of energy between pre- and post-collapse", 100, 0.5, 1.5);
164 }
165 
Log< level::Info, true > LogVerbatim
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
size_type size() const
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:224
edm::ESGetToken< HcalTopology, HcalRecNumberingRecord > hcalTopologyToken_
const HcalTopology * theHBHETopology
const edm::InputTag preRecHitHBHE_
void analyze(edm::Event const &, edm::EventSetup const &) override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
unsigned int id
HcalCollapseAnalyzer(const edm::ParameterSet &)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::EDGetTokenT< HBHERecHitCollection > tok_prehbhe_
const edm::InputTag recHitHBHE_
HLT enums.
void unmergeDepthDetId(const HcalDetId &id, std::vector< HcalDetId > &ids) const
Definition: HcalTopology.h:165
const std::string topFolderName_
MonitorElement * h_balance
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
Definition: Run.h:45
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::EDGetTokenT< HBHERecHitCollection > tok_hbhe_