CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
HcalCollapseAnalyzer Class Reference

#include <HcalCollapseAnalyzer.cc>

Inheritance diagram for HcalCollapseAnalyzer:
one::DQMEDAnalyzer< T > one::dqmimplementation::DQMBaseClass< T... >

Public Member Functions

 HcalCollapseAnalyzer (const edm::ParameterSet &)
 
 ~HcalCollapseAnalyzer () override
 
- Public Member Functions inherited from one::DQMEDAnalyzer< T >
 DQMEDAnalyzer ()=default
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > const &)=delete
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > &&)=delete
 
 ~DQMEDAnalyzer () override=default
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

Private Member Functions

void analyze (edm::Event const &, edm::EventSetup const &) override
 
void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 

Private Attributes

const bool doHB_
 
const bool doHE_
 
MonitorElementh_balance
 
MonitorElementh_depth
 
MonitorElementh_frac
 
MonitorElementh_merge
 
MonitorElementh_sfrac
 
MonitorElementh_size
 
const edm::InputTag preRecHitHBHE_
 
const edm::InputTag recHitHBHE_
 
edm::EDGetTokenT< HBHERecHitCollectiontok_hbhe_
 
edm::EDGetTokenT< HBHERecHitCollectiontok_prehbhe_
 
const std::string topFolderName_
 
const int verbosity_
 

Detailed Description

DQMOffline/Hcal/src/HcalCollapseAnalyzer.cc

Description: Studies the collapsing of HcalRecHits

Implementation: <Notes on="" implementation>="">

Definition at line 48 of file HcalCollapseAnalyzer.cc.

Constructor & Destructor Documentation

HcalCollapseAnalyzer::HcalCollapseAnalyzer ( const edm::ParameterSet iConfig)
explicit

Definition at line 72 of file HcalCollapseAnalyzer.cc.

References doHB_, doHE_, preRecHitHBHE_, recHitHBHE_, tok_hbhe_, tok_prehbhe_, and verbosity_.

73  : topFolderName_(iConfig.getParameter<std::string>("topFolderName")),
74  verbosity_(iConfig.getUntrackedParameter<int>("verbosity", 0)),
75  recHitHBHE_(iConfig.getUntrackedParameter<edm::InputTag>("recHitHBHE", edm::InputTag("hbhereco"))),
76  preRecHitHBHE_(iConfig.getUntrackedParameter<edm::InputTag>("preRecHitHBHE", edm::InputTag("hbheprereco"))),
77  doHE_(iConfig.getUntrackedParameter<bool>("doHE", true)),
78  doHB_(iConfig.getUntrackedParameter<bool>("doHB", false)) {
79  // define tokens for access
80  tok_hbhe_ = consumes<HBHERecHitCollection>(recHitHBHE_);
81  tok_prehbhe_ = consumes<HBHERecHitCollection>(preRecHitHBHE_);
82 
83  edm::LogVerbatim("Collapse") << "Verbosity " << verbosity_ << " with tags " << recHitHBHE_ << " and "
84  << preRecHitHBHE_ << " and Do " << doHB_ << ":" << doHE_;
85 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
const edm::InputTag preRecHitHBHE_
edm::EDGetTokenT< HBHERecHitCollection > tok_prehbhe_
const edm::InputTag recHitHBHE_
const std::string topFolderName_
edm::EDGetTokenT< HBHERecHitCollection > tok_hbhe_
HcalCollapseAnalyzer::~HcalCollapseAnalyzer ( )
inlineoverride

Definition at line 51 of file HcalCollapseAnalyzer.cc.

References analyze(), bookHistograms(), and fillDescriptions().

51 {}

Member Function Documentation

void HcalCollapseAnalyzer::analyze ( edm::Event const &  iEvent,
edm::EventSetup const &  iSetup 
)
overrideprivate

Definition at line 98 of file HcalCollapseAnalyzer.cc.

References egammaForCoreTracking_cff::depth, doHB_, doHE_, edm::SortedCollection< T, SORT >::empty(), edm::EventID::event(), MonitorElement::Fill(), spr::find(), DivergingColor::frac, edm::EventSetup::get(), edm::Event::getByToken(), h_balance, h_depth, h_frac, h_merge, h_sfrac, h_size, hcalLocalReco_cff::hbheprereco, HiIsolationCommonParameters_cff::hbhereco, HcalBarrel, HcalEndcap, edm::EventBase::id(), hit::id, edm::HandleBase::isValid(), gen::k, edm::Handle< T >::product(), edm::ESHandle< T >::product(), edm::EventID::run(), edm::SortedCollection< T, SORT >::size(), tok_hbhe_, tok_prehbhe_, HcalTopology::unmergeDepthDetId(), and verbosity_.

Referenced by ~HcalCollapseAnalyzer().

98  {
99  if (verbosity_ > 0)
100  edm::LogVerbatim("Collapse") << "Run " << iEvent.id().run() << " Event " << iEvent.id().event()
101  << " starts ==============";
102 
104  iSetup.get<HcalRecNumberingRecord>().get(htopo);
105  const HcalTopology *theHBHETopology = htopo.product();
106 
108  iEvent.getByToken(tok_hbhe_, hbhereco);
110  iEvent.getByToken(tok_prehbhe_, hbheprereco);
111  if (verbosity_ > 0) {
112  edm::LogVerbatim("Collapse") << "Handle Reco " << hbhereco << " Size " << hbhereco->size();
113  edm::LogVerbatim("Collapse") << "Handle PreReco " << hbheprereco << " Size " << hbheprereco->size();
114  }
115  if (hbhereco.isValid() && hbheprereco.isValid()) {
116  const HBHERecHitCollection *recohbhe = hbhereco.product();
117  const HBHERecHitCollection *prerecohbhe = hbheprereco.product();
118  if (verbosity_ > 0)
119  edm::LogVerbatim("Collapse") << "Size of hbhereco: " << recohbhe->size()
120  << " and hbheprereco: " << prerecohbhe->size();
121  double sfrac = (prerecohbhe->empty()) ? 1 : ((double)(recohbhe->size())) / ((double)(prerecohbhe->size()));
122  h_sfrac->Fill(sfrac);
123  h_size->Fill(recohbhe->size());
124  for (const auto &hit : (*recohbhe)) {
125  HcalDetId id = hit.id();
126  if (((id.subdet() == HcalEndcap) && doHE_) || ((id.subdet() == HcalBarrel) && doHB_)) {
127  h_depth->Fill(id.depth());
128  std::vector<HcalDetId> ids;
129  theHBHETopology->unmergeDepthDetId(id, ids);
130  if (verbosity_ > 0) {
131  edm::LogVerbatim("Collapse") << id << " is derived from " << ids.size() << " components";
132  for (unsigned int k = 0; k < ids.size(); ++k)
133  edm::LogVerbatim("Collapse") << "[" << k << "] " << ids[k];
134  }
135  h_merge->Fill(ids.size());
136  double energy = hit.energy();
137  double etot(0);
138  unsigned int k(0);
139  for (const auto phit : (*prerecohbhe)) {
140  if (std::find(ids.begin(), ids.end(), phit.id()) != ids.end()) {
141  etot += phit.energy();
142  double frac = (energy > 0) ? phit.energy() / energy : 1.0;
143  h_frac->Fill(frac);
144  if (verbosity_ > 0)
145  edm::LogVerbatim("Collapse") << "Frac[" << k << "] " << frac << " for " << phit.id();
146  ++k;
147  }
148  }
149  double frac = (energy > 0) ? etot / energy : 1.0;
150  h_balance->Fill(frac);
151  if (verbosity_ > 0)
152  edm::LogVerbatim("Collapse") << "Total energy " << energy << ":" << etot << ":" << frac;
153  }
154  }
155  }
156 }
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:224
void unmergeDepthDetId(const HcalDetId &id, std::vector< HcalDetId > &ids) const
Definition: HcalTopology.h:172
bool isValid() const
Definition: HandleBase.h:74
int k[5][pyjets_maxn]
unsigned int id
T const * product() const
Definition: Handle.h:74
edm::EDGetTokenT< HBHERecHitCollection > tok_prehbhe_
size_type size() const
MonitorElement * h_balance
T const * product() const
Definition: ESHandle.h:86
edm::EDGetTokenT< HBHERecHitCollection > tok_hbhe_
void HcalCollapseAnalyzer::bookHistograms ( DQMStore::IBooker ibooker,
edm::Run const &  ,
edm::EventSetup const &   
)
overrideprivate

Definition at line 158 of file HcalCollapseAnalyzer.cc.

References DQMStore::IBooker::book1D(), DEFINE_FWK_MODULE, h_balance, h_depth, h_frac, h_merge, h_sfrac, h_size, DQMStore::IBooker::setCurrentFolder(), and topFolderName_.

Referenced by ~HcalCollapseAnalyzer().

158  {
159  // Book histograms
161  h_merge = ibooker.book1D("h_merge", "Number of hits merged", 10, 0.0, 10.0);
162  h_size = ibooker.book1D("h_size", "Size of the RecHit collection", 100, 500.0, 1500.0);
163  h_depth = ibooker.book1D("h_depth", "Depth of the Id's used", 10, 0.0, 10.0);
164  h_sfrac = ibooker.book1D("h_sfrac", "Ratio of sizes of preRecHit and RecHit collections", 150, 0.0, 1.5);
165  h_frac = ibooker.book1D("h_frac", "Fraction of energy before collapse", 150, 0.0, 1.5);
166  h_balance = ibooker.book1D("h_balance", "Balance of energy between pre- and post-collapse", 100, 0.5, 1.5);
167 }
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
const std::string topFolderName_
MonitorElement * h_balance
void HcalCollapseAnalyzer::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 87 of file HcalCollapseAnalyzer.cc.

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::add(), edm::ParameterSetDescription::addUntracked(), and AlCaHLTBitMon_QueryRunRegistry::string.

Referenced by ~HcalCollapseAnalyzer().

87  {
89  desc.add<std::string>("topFolderName", "HcalCollapse");
90  desc.addUntracked<int>("verbosity", 0);
91  desc.addUntracked<edm::InputTag>("recHitHBHE", edm::InputTag("hbhereco"));
92  desc.addUntracked<edm::InputTag>("preRecHitHBHE", edm::InputTag("hbheprereco"));
93  desc.addUntracked<bool>("doHE", true);
94  desc.addUntracked<bool>("doHB", false);
95  descriptions.add("hcalCollapseAnalyzer", desc);
96 }
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)

Member Data Documentation

const bool HcalCollapseAnalyzer::doHB_
private

Definition at line 63 of file HcalCollapseAnalyzer.cc.

Referenced by analyze(), and HcalCollapseAnalyzer().

const bool HcalCollapseAnalyzer::doHE_
private

Definition at line 63 of file HcalCollapseAnalyzer.cc.

Referenced by analyze(), and HcalCollapseAnalyzer().

MonitorElement * HcalCollapseAnalyzer::h_balance
private

Definition at line 69 of file HcalCollapseAnalyzer.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement * HcalCollapseAnalyzer::h_depth
private

Definition at line 68 of file HcalCollapseAnalyzer.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement * HcalCollapseAnalyzer::h_frac
private

Definition at line 69 of file HcalCollapseAnalyzer.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* HcalCollapseAnalyzer::h_merge
private

Definition at line 68 of file HcalCollapseAnalyzer.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* HcalCollapseAnalyzer::h_sfrac
private

Definition at line 69 of file HcalCollapseAnalyzer.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement * HcalCollapseAnalyzer::h_size
private

Definition at line 68 of file HcalCollapseAnalyzer.cc.

Referenced by analyze(), and bookHistograms().

const edm::InputTag HcalCollapseAnalyzer::preRecHitHBHE_
private

Definition at line 62 of file HcalCollapseAnalyzer.cc.

Referenced by HcalCollapseAnalyzer().

const edm::InputTag HcalCollapseAnalyzer::recHitHBHE_
private

Definition at line 62 of file HcalCollapseAnalyzer.cc.

Referenced by HcalCollapseAnalyzer().

edm::EDGetTokenT<HBHERecHitCollection> HcalCollapseAnalyzer::tok_hbhe_
private

Definition at line 65 of file HcalCollapseAnalyzer.cc.

Referenced by analyze(), and HcalCollapseAnalyzer().

edm::EDGetTokenT<HBHERecHitCollection> HcalCollapseAnalyzer::tok_prehbhe_
private

Definition at line 66 of file HcalCollapseAnalyzer.cc.

Referenced by analyze(), and HcalCollapseAnalyzer().

const std::string HcalCollapseAnalyzer::topFolderName_
private

Definition at line 60 of file HcalCollapseAnalyzer.cc.

Referenced by bookHistograms().

const int HcalCollapseAnalyzer::verbosity_
private

Definition at line 61 of file HcalCollapseAnalyzer.cc.

Referenced by analyze(), and HcalCollapseAnalyzer().