CMS 3D CMS Logo

HGCalHitCalibration.cc
Go to the documentation of this file.
1 // user include files
2 
6 
10 
12 
14 
17 
24 
29 
35 
37 
38 #include <map>
39 #include <array>
40 #include <string>
41 #include <numeric>
42 
44 public:
45  explicit HGCalHitCalibration(const edm::ParameterSet&);
46  ~HGCalHitCalibration() override;
47 
48  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
49 
50 private:
51  void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
52  void analyze(const edm::Event&, const edm::EventSetup&) override;
53 
54  void fillWithRecHits(std::map<DetId, const HGCRecHit*>&, DetId, unsigned int, float, int&, float&);
55 
63 
67  int debug_;
69  static constexpr int depletion1_ = 200;
70  static constexpr int depletion2_ = 300;
71  static constexpr int scint_ = 400;
72 
73  std::map<int, MonitorElement*> h_EoP_CPene_calib_fraction_;
74  std::map<int, MonitorElement*> hgcal_EoP_CPene_calib_fraction_;
75  std::map<int, MonitorElement*> hgcal_ele_EoP_CPene_calib_fraction_;
76  std::map<int, MonitorElement*> hgcal_photon_EoP_CPene_calib_fraction_;
78 
79  static constexpr int layers_ = 60;
80  std::array<float, layers_> Energy_layer_calib_;
81  std::array<float, layers_> Energy_layer_calib_fraction_;
82 };
83 
85  : caloGeomToken_(esConsumes<CaloGeometry, CaloGeometryRecord>()),
86  rawRecHits_(iConfig.getParameter<bool>("rawRecHits")),
87  debug_(iConfig.getParameter<int>("debug")) {
88  auto detector = iConfig.getParameter<std::string>("detector");
89  auto recHitsEE = iConfig.getParameter<edm::InputTag>("recHitsEE");
90  auto recHitsFH = iConfig.getParameter<edm::InputTag>("recHitsFH");
91  auto recHitsBH = iConfig.getParameter<edm::InputTag>("recHitsBH");
92  auto caloParticles = iConfig.getParameter<edm::InputTag>("caloParticles");
93  auto hgcalMultiClusters = iConfig.getParameter<edm::InputTag>("hgcalMultiClusters");
94  auto electrons = iConfig.getParameter<edm::InputTag>("electrons");
95  auto photons = iConfig.getParameter<edm::InputTag>("photons");
96  if (detector == "all") {
97  recHitsEE_ = consumes<HGCRecHitCollection>(recHitsEE);
98  recHitsFH_ = consumes<HGCRecHitCollection>(recHitsFH);
99  recHitsBH_ = consumes<HGCRecHitCollection>(recHitsBH);
100  algo_ = 1;
101  } else if (detector == "EM") {
102  recHitsEE_ = consumes<HGCRecHitCollection>(recHitsEE);
103  algo_ = 2;
104  } else if (detector == "HAD") {
105  recHitsFH_ = consumes<HGCRecHitCollection>(recHitsFH);
106  recHitsBH_ = consumes<HGCRecHitCollection>(recHitsBH);
107  algo_ = 3;
108  }
109  depletion0_ = iConfig.getParameter<int>("depletionFine");
110  caloParticles_ = consumes<std::vector<CaloParticle> >(caloParticles);
111  hgcalMultiClusters_ = consumes<std::vector<reco::PFCluster> >(hgcalMultiClusters);
112  electrons_ = consumes<std::vector<reco::GsfElectron> >(electrons);
113  photons_ = consumes<std::vector<reco::Photon> >(photons);
114 
115  // size should be HGC layers 52 is enough
116  Energy_layer_calib_.fill(0.);
118 }
119 
121  // do anything here that needs to be done at desctruction time
122  // (e.g. close files, deallocate resources etc.)
123 }
124 
126  edm::Run const& iRun,
127  edm::EventSetup const& /* iSetup */) {
128  ibooker.cd();
129  ibooker.setCurrentFolder("HGCalHitCalibration");
130  h_EoP_CPene_calib_fraction_[depletion0_] = ibooker.book1D("h_EoP_CPene_100_calib_fraction", "", 1000, -0.5, 2.5);
131  h_EoP_CPene_calib_fraction_[depletion1_] = ibooker.book1D("h_EoP_CPene_200_calib_fraction", "", 1000, -0.5, 2.5);
132  h_EoP_CPene_calib_fraction_[depletion2_] = ibooker.book1D("h_EoP_CPene_300_calib_fraction", "", 1000, -0.5, 2.5);
133  h_EoP_CPene_calib_fraction_[scint_] = ibooker.book1D("h_EoP_CPene_scint_calib_fraction", "", 1000, -0.5, 2.5);
135  ibooker.book1D("hgcal_EoP_CPene_100_calib_fraction", "", 1000, -0.5, 2.5);
137  ibooker.book1D("hgcal_EoP_CPene_200_calib_fraction", "", 1000, -0.5, 2.5);
139  ibooker.book1D("hgcal_EoP_CPene_300_calib_fraction", "", 1000, -0.5, 2.5);
140  hgcal_EoP_CPene_calib_fraction_[scint_] = ibooker.book1D("hgcal_EoP_CPene_scint_calib_fraction", "", 1000, -0.5, 2.5);
142  ibooker.book1D("hgcal_ele_EoP_CPene_100_calib_fraction", "", 1000, -0.5, 2.5);
144  ibooker.book1D("hgcal_ele_EoP_CPene_200_calib_fraction", "", 1000, -0.5, 2.5);
146  ibooker.book1D("hgcal_ele_EoP_CPene_300_calib_fraction", "", 1000, -0.5, 2.5);
148  ibooker.book1D("hgcal_ele_EoP_CPene_scint_calib_fraction", "", 1000, -0.5, 2.5);
150  ibooker.book1D("hgcal_photon_EoP_CPene_100_calib_fraction", "", 1000, -0.5, 2.5);
152  ibooker.book1D("hgcal_photon_EoP_CPene_200_calib_fraction", "", 1000, -0.5, 2.5);
154  ibooker.book1D("hgcal_photon_EoP_CPene_300_calib_fraction", "", 1000, -0.5, 2.5);
156  ibooker.book1D("hgcal_photon_EoP_CPene_scint_calib_fraction", "", 1000, -0.5, 2.5);
157  LayerOccupancy_ = ibooker.book1D("LayerOccupancy", "", layers_, 0., (float)layers_);
158 }
159 
160 void HGCalHitCalibration::fillWithRecHits(std::map<DetId, const HGCRecHit*>& hitmap,
161  const DetId hitid,
162  const unsigned int hitlayer,
163  const float fraction,
164  int& seedDet,
165  float& seedEnergy) {
166  if (hitmap.find(hitid) == hitmap.end()) {
167  // Hit was not reconstructed
168  IfLogTrace(debug_ > 0, "HGCalHitCalibration")
169  << ">>> Failed to find detid " << std::hex << hitid.rawId() << std::dec << " Det " << hitid.det() << " Subdet "
170  << hitid.subdetId() << std::endl;
171  return;
172  }
173  if ((hitid.det() != DetId::Forward) && (hitid.det() != DetId::HGCalEE) && (hitid.det() != DetId::HGCalHSi) &&
174  (hitid.det() != DetId::HGCalHSc)) {
175  IfLogTrace(debug_ > 0, "HGCalHitCalibration")
176  << ">>> Wrong type of detid " << std::hex << hitid.rawId() << std::dec << " Det " << hitid.det() << " Subdet "
177  << hitid.subdetId() << std::endl;
178  return;
179  }
180 
181  unsigned int layer = recHitTools_.getLayerWithOffset(hitid);
182  assert(hitlayer == layer);
183  Energy_layer_calib_fraction_[layer] += hitmap[hitid]->energy() * fraction;
185  if (seedEnergy < hitmap[hitid]->energy()) {
186  seedEnergy = hitmap[hitid]->energy();
187  seedDet = recHitTools_.getSiThickness(hitid);
188  if (hitid.det() == DetId::HGCalHSc) {
189  seedDet = scint_;
190  }
191  }
192 }
193 
195  using namespace edm;
196 
199 
200  Handle<HGCRecHitCollection> recHitHandleEE;
201  Handle<HGCRecHitCollection> recHitHandleFH;
202  Handle<HGCRecHitCollection> recHitHandleBH;
203 
204  Handle<std::vector<CaloParticle> > caloParticleHandle;
205  iEvent.getByToken(caloParticles_, caloParticleHandle);
206  const std::vector<CaloParticle>& caloParticles = *caloParticleHandle;
207 
208  Handle<std::vector<reco::PFCluster> > hgcalMultiClustersHandle;
209  iEvent.getByToken(hgcalMultiClusters_, hgcalMultiClustersHandle);
210 
211  Handle<std::vector<reco::GsfElectron> > PFElectronHandle;
212  iEvent.getByToken(electrons_, PFElectronHandle);
213 
214  Handle<std::vector<reco::Photon> > PFPhotonHandle;
215  iEvent.getByToken(photons_, PFPhotonHandle);
216 
217  // make a map detid-rechit
218  std::map<DetId, const HGCRecHit*> hitmap;
219  switch (algo_) {
220  case 1: {
221  iEvent.getByToken(recHitsEE_, recHitHandleEE);
222  iEvent.getByToken(recHitsFH_, recHitHandleFH);
223  iEvent.getByToken(recHitsBH_, recHitHandleBH);
224  const auto& rechitsEE = *recHitHandleEE;
225  const auto& rechitsFH = *recHitHandleFH;
226  const auto& rechitsBH = *recHitHandleBH;
227  for (unsigned int i = 0; i < rechitsEE.size(); ++i) {
228  hitmap[rechitsEE[i].detid()] = &rechitsEE[i];
229  }
230  for (unsigned int i = 0; i < rechitsFH.size(); ++i) {
231  hitmap[rechitsFH[i].detid()] = &rechitsFH[i];
232  }
233  for (unsigned int i = 0; i < rechitsBH.size(); ++i) {
234  hitmap[rechitsBH[i].detid()] = &rechitsBH[i];
235  }
236  break;
237  }
238  case 2: {
239  iEvent.getByToken(recHitsEE_, recHitHandleEE);
240  const HGCRecHitCollection& rechitsEE = *recHitHandleEE;
241  for (unsigned int i = 0; i < rechitsEE.size(); i++) {
242  hitmap[rechitsEE[i].detid()] = &rechitsEE[i];
243  }
244  break;
245  }
246  case 3: {
247  iEvent.getByToken(recHitsFH_, recHitHandleFH);
248  iEvent.getByToken(recHitsBH_, recHitHandleBH);
249  const auto& rechitsFH = *recHitHandleFH;
250  const auto& rechitsBH = *recHitHandleBH;
251  for (unsigned int i = 0; i < rechitsFH.size(); i++) {
252  hitmap[rechitsFH[i].detid()] = &rechitsFH[i];
253  }
254  for (unsigned int i = 0; i < rechitsBH.size(); i++) {
255  hitmap[rechitsBH[i].detid()] = &rechitsBH[i];
256  }
257  break;
258  }
259  default:
260  assert(false);
261  break;
262  }
263 
264  // loop over caloParticles
265  int seedDet = 0;
266  float seedEnergy = 0.;
267  IfLogTrace(debug_ > 0, "HGCalHitCalibration") << "Number of caloParticles: " << caloParticles.size() << std::endl;
268  for (const auto& it_caloPart : caloParticles) {
269  if (it_caloPart.g4Tracks()[0].eventId().event() != 0 or it_caloPart.g4Tracks()[0].eventId().bunchCrossing() != 0) {
270  LogDebug("HGCalHitCalibration") << "Excluding CaloParticles from event: "
271  << it_caloPart.g4Tracks()[0].eventId().event()
272  << " with BX: " << it_caloPart.g4Tracks()[0].eventId().bunchCrossing()
273  << std::endl;
274  continue;
275  }
276 
277  const SimClusterRefVector& simClusterRefVector = it_caloPart.simClusters();
278  Energy_layer_calib_.fill(0.);
280 
281  seedDet = 0;
282  seedEnergy = 0.;
283  for (const auto& it_sc : simClusterRefVector) {
284  const SimCluster& simCluster = (*(it_sc));
285  IfLogTrace(debug_ > 1, "HGCalHitCalibration")
286  << ">>> SC.energy(): " << simCluster.energy() << " SC.simEnergy(): " << simCluster.simEnergy() << std::endl;
287  const std::vector<std::pair<uint32_t, float> >& hits_and_fractions = simCluster.hits_and_fractions();
288 
289  // loop over hits
290  for (const auto& it_haf : hits_and_fractions) {
291  unsigned int hitlayer = recHitTools_.getLayerWithOffset(it_haf.first);
292  DetId hitid = (it_haf.first);
293  // dump raw RecHits and match
294  if (rawRecHits_) {
295  if ((hitid.det() == DetId::Forward &&
296  (hitid.subdetId() == HGCEE || hitid.subdetId() == HGCHEF || hitid.subdetId() == HGCHEB)) ||
297  (hitid.det() == DetId::HGCalEE) || (hitid.det() == DetId::HGCalHSi) || (hitid.det() == DetId::HGCalHSc) ||
298  (hitid.det() == DetId::Hcal && hitid.subdetId() == HcalEndcap))
299  fillWithRecHits(hitmap, hitid, hitlayer, it_haf.second, seedDet, seedEnergy);
300  }
301  } // end simHit
302  } // end simCluster
303 
304  auto sumCalibRecHitCalib_fraction =
305  std::accumulate(Energy_layer_calib_fraction_.begin(), Energy_layer_calib_fraction_.end(), 0.);
306  IfLogTrace(debug_ > 0, "HGCalHitCalibration")
307  << ">>> MC Energy: " << it_caloPart.energy() << " reco energy: " << sumCalibRecHitCalib_fraction << std::endl;
308  if (h_EoP_CPene_calib_fraction_.count(seedDet))
309  h_EoP_CPene_calib_fraction_[seedDet]->Fill(sumCalibRecHitCalib_fraction / it_caloPart.energy());
310 
311  // Loop on reconstructed SC.
312  const auto& clusters = *hgcalMultiClustersHandle;
313  float total_energy = 0.;
314  float max_dR2 = 0.0025;
315  auto closest =
316  std::min_element(clusters.begin(), clusters.end(), [&](const reco::PFCluster& a, const reco::PFCluster& b) {
317  auto dR2_a = reco::deltaR2(it_caloPart, a);
318  auto dR2_b = reco::deltaR2(it_caloPart, b);
319  auto ERatio_a = a.correctedEnergy() / it_caloPart.energy();
320  auto ERatio_b = b.correctedEnergy() / it_caloPart.energy();
321  // If both clusters are within 0.0025, mark as first (min) the
322  // element with the highest ratio against the SimCluster
323  if (dR2_a < max_dR2 && dR2_b < max_dR2)
324  return ERatio_a > ERatio_b;
325  return dR2_a < dR2_b;
326  });
327  if (closest != clusters.end() && reco::deltaR2(*closest, it_caloPart) < 0.01) {
328  total_energy = closest->correctedEnergy();
329  seedDet = recHitTools_.getSiThickness(closest->seed());
330  if (closest->seed().det() == DetId::HGCalHSc) {
331  seedDet = scint_;
332  }
333  if (hgcal_EoP_CPene_calib_fraction_.count(seedDet)) {
334  hgcal_EoP_CPene_calib_fraction_[seedDet]->Fill(total_energy / it_caloPart.energy());
335  }
336  }
337 
338  auto closest_fcn = [&](auto const& a, auto const& b) {
339  auto dR2_a = reco::deltaR2(it_caloPart, a);
340  auto dR2_b = reco::deltaR2(it_caloPart, b);
341  auto ERatio_a = a.energy() / it_caloPart.energy();
342  auto ERatio_b = b.energy() / it_caloPart.energy();
343  // If both clusters are within 0.0025, mark as first (min) the
344  // element with the highest ratio against the SimCluster
345  if (dR2_a < max_dR2 && dR2_b < max_dR2)
346  return ERatio_a > ERatio_b;
347  return dR2_a < dR2_b;
348  };
349  // ELECTRONS in HGCAL
350  if (PFElectronHandle.isValid()) {
351  auto const& ele = (*PFElectronHandle);
352  auto closest = std::min_element(ele.begin(), ele.end(), closest_fcn);
353  if (closest != ele.end() &&
354  (closest->superCluster()->seed()->seed().det() == DetId::Forward ||
355  closest->superCluster()->seed()->seed().det() == DetId::HGCalEE) &&
356  reco::deltaR2(*closest, it_caloPart) < 0.01) {
357  seedDet = recHitTools_.getSiThickness(closest->superCluster()->seed()->seed());
358  if (closest->superCluster()->seed()->seed().det() == DetId::HGCalHSc) {
359  seedDet = scint_;
360  }
361  if (hgcal_ele_EoP_CPene_calib_fraction_.count(seedDet)) {
362  hgcal_ele_EoP_CPene_calib_fraction_[seedDet]->Fill(closest->energy() / it_caloPart.energy());
363  }
364  }
365  }
366 
367  // PHOTONS in HGCAL
368  if (PFPhotonHandle.isValid()) {
369  auto const& photon = (*PFPhotonHandle);
370  auto closest = std::min_element(photon.begin(), photon.end(), closest_fcn);
371  if (closest != photon.end() &&
372  (closest->superCluster()->seed()->seed().det() == DetId::Forward ||
373  closest->superCluster()->seed()->seed().det() == DetId::HGCalEE) &&
374  reco::deltaR2(*closest, it_caloPart) < 0.01) {
375  seedDet = recHitTools_.getSiThickness(closest->superCluster()->seed()->seed());
376  if (hgcal_photon_EoP_CPene_calib_fraction_.count(seedDet)) {
377  hgcal_photon_EoP_CPene_calib_fraction_[seedDet]->Fill(closest->energy() / it_caloPart.energy());
378  }
379  }
380  }
381  } // end caloparticle
382 }
383 
384 // ------------ method fills 'descriptions' with the allowed parameters for the
385 // module ------------
388  desc.add<int>("debug", 0);
389  desc.add<bool>("rawRecHits", true);
390  desc.add<std::string>("detector", "all");
391  desc.add<int>("depletionFine", 120);
392  desc.add<edm::InputTag>("caloParticles", edm::InputTag("mix", "MergedCaloTruth"));
393  desc.add<edm::InputTag>("recHitsEE", edm::InputTag("HGCalRecHit", "HGCEERecHits"));
394  desc.add<edm::InputTag>("recHitsFH", edm::InputTag("HGCalRecHit", "HGCHEFRecHits"));
395  desc.add<edm::InputTag>("recHitsBH", edm::InputTag("HGCalRecHit", "HGCHEBRecHits"));
396  desc.add<edm::InputTag>("hgcalMultiClusters", edm::InputTag("particleFlowClusterHGCalFromMultiCl"));
397  desc.add<edm::InputTag>("electrons", edm::InputTag("ecalDrivenGsfElectronsFromMultiCl"));
398  desc.add<edm::InputTag>("photons", edm::InputTag("photonsFromMultiCl"));
399  descriptions.add("hgcalHitCalibrationDefault", desc);
400 }
401 
402 // define this as a plug-in
HGCalHitCalibration::recHitsBH_
edm::EDGetTokenT< HGCRecHitCollection > recHitsBH_
Definition: HGCalHitCalibration.cc:58
muons2muons_cfi.photon
photon
Definition: muons2muons_cfi.py:28
dqm::impl::MonitorElement
Definition: MonitorElement.h:98
hgcal::RecHitTools
Definition: RecHitTools.h:23
electrons_cff.bool
bool
Definition: electrons_cff.py:366
EDAnalyzer.h
mps_fire.i
i
Definition: mps_fire.py:428
HGCalHitCalibration::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: HGCalHitCalibration.cc:386
hgcal::RecHitTools::getSiThickness
std::float_t getSiThickness(const DetId &) const
Definition: RecHitTools.cc:176
ESHandle.h
HGCalHitCalibration::algo_
int algo_
Definition: HGCalHitCalibration.cc:64
edm::Run
Definition: Run.h:45
edm::EDGetTokenT
Definition: EDGetToken.h:33
CaloGeometryRecord
Definition: CaloGeometryRecord.h:30
edm
HLT enums.
Definition: AlignableModifier.h:19
DetId::det
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:46
DetId::Hcal
Definition: DetId.h:28
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89301
RecHitTools.h
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
dqm::implementation::NavigatorBase::setCurrentFolder
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
cms::cuda::assert
assert(be >=bs)
HGCalHitCalibration::fillWithRecHits
void fillWithRecHits(std::map< DetId, const HGCRecHit * > &, DetId, unsigned int, float, int &, float &)
Definition: HGCalHitCalibration.cc:160
edm::SortedCollection
Definition: SortedCollection.h:49
DQMStore.h
edm::SortedCollection::size
size_type size() const
Definition: SortedCollection.h:215
TrackingVertex.h
edm::RefVector< SimClusterCollection >
SimCluster::energy
float energy() const
Energy. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:104
SimCluster
Monte Carlo truth information used for tracking validation.
Definition: SimCluster.h:29
HGCalHitCalibration::debug_
int debug_
Definition: HGCalHitCalibration.cc:67
HGCalHitCalibration::scint_
static constexpr int scint_
Definition: HGCalHitCalibration.cc:71
HGCalHitCalibration::hgcal_photon_EoP_CPene_calib_fraction_
std::map< int, MonitorElement * > hgcal_photon_EoP_CPene_calib_fraction_
Definition: HGCalHitCalibration.cc:76
edm::Handle
Definition: AssociativeIterator.h:50
HGCalHitCalibration::h_EoP_CPene_calib_fraction_
std::map< int, MonitorElement * > h_EoP_CPene_calib_fraction_
Definition: HGCalHitCalibration.cc:73
DetId
Definition: DetId.h:17
HGCalHitCalibration::recHitsEE_
edm::EDGetTokenT< HGCRecHitCollection > recHitsEE_
Definition: HGCalHitCalibration.cc:56
HGCalHitCalibration::hgcal_EoP_CPene_calib_fraction_
std::map< int, MonitorElement * > hgcal_EoP_CPene_calib_fraction_
Definition: HGCalHitCalibration.cc:74
DetId::HGCalHSi
Definition: DetId.h:33
DetId::HGCalEE
Definition: DetId.h:32
MakerMacros.h
Photon.h
CaloGeometry
Definition: CaloGeometry.h:21
HGCalHitCalibration::LayerOccupancy_
MonitorElement * LayerOccupancy_
Definition: HGCalHitCalibration.cc:77
SimCluster.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
HGCalHitCalibration::depletion1_
static constexpr int depletion1_
Definition: HGCalHitCalibration.cc:69
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
HLT_FULL_cff.fraction
fraction
Definition: HLT_FULL_cff.py:52823
TruncatedPyramid.h
caloTruthCellsProducer_cfi.caloParticles
caloParticles
Definition: caloTruthCellsProducer_cfi.py:6
dqm::impl::MonitorElement::Fill
void Fill(long long x)
Definition: MonitorElement.h:290
HGCalHitCalibration::Energy_layer_calib_
std::array< float, layers_ > Energy_layer_calib_
Definition: HGCalHitCalibration.cc:80
edm::ESHandle< CaloGeometry >
dqm::implementation::NavigatorBase::cd
virtual void cd()
Definition: DQMStore.cc:29
PFCluster.h
relativeConstraints.geom
geom
Definition: relativeConstraints.py:72
CaloClusterFwd.h
HCALHighEnergyHPDFilter_cfi.energy
energy
Definition: HCALHighEnergyHPDFilter_cfi.py:5
IfLogTrace
#define IfLogTrace(cond, cat)
Definition: MessageLogger.h:260
HGCalHitCalibration::photons_
edm::EDGetTokenT< std::vector< reco::Photon > > photons_
Definition: HGCalHitCalibration.cc:62
GsfElectron.h
b
double b
Definition: hdecay.h:118
HGCalHitCalibration::depletion0_
int depletion0_
Definition: HGCalHitCalibration.cc:64
DetId::subdetId
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector's numbering enum)
Definition: DetId.h:48
DQMEDAnalyzer.h
CaloGeometryRecord.h
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
phase1PixelTopology::layer
constexpr std::array< uint8_t, layerIndexSize > layer
Definition: phase1PixelTopology.h:99
HGCEE
Definition: ForwardSubdetector.h:8
bsc_activity_cfg.clusters
clusters
Definition: bsc_activity_cfg.py:36
DQMEDAnalyzer
Definition: DQMEDAnalyzer.py:1
HGCalHitCalibration::rawRecHits_
bool rawRecHits_
Definition: HGCalHitCalibration.cc:66
CaloSubdetectorGeometry.h
HGCRecHitCollections.h
HGCalGeometry.h
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:233
edm::ParameterSet
Definition: ParameterSet.h:47
a
double a
Definition: hdecay.h:119
Event.h
HGCalHitCalibration::hgcal_ele_EoP_CPene_calib_fraction_
std::map< int, MonitorElement * > hgcal_ele_EoP_CPene_calib_fraction_
Definition: HGCalHitCalibration.cc:75
hgcal::RecHitTools::getLayerWithOffset
unsigned int getLayerWithOffset(const DetId &) const
Definition: RecHitTools.cc:362
HcalDetId.h
reco::deltaR2
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
CaloParticle.h
createfilelist.int
int
Definition: createfilelist.py:10
iEvent
int iEvent
Definition: GenABIO.cc:224
HGCalHitCalibration::HGCalHitCalibration
HGCalHitCalibration(const edm::ParameterSet &)
Definition: HGCalHitCalibration.cc:84
edm::EventSetup::getHandle
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:155
HGCalHitCalibration::Energy_layer_calib_fraction_
std::array< float, layers_ > Energy_layer_calib_fraction_
Definition: HGCalHitCalibration.cc:81
BPHMonitor_cfi.photons
photons
Definition: BPHMonitor_cfi.py:91
edm::EventSetup
Definition: EventSetup.h:58
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
edm::ESGetToken< CaloGeometry, CaloGeometryRecord >
HcalObjRepresent::Fill
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
Definition: HcalObjRepresent.h:1053
HGCalHitCalibration::hgcalMultiClusters_
edm::EDGetTokenT< std::vector< reco::PFCluster > > hgcalMultiClusters_
Definition: HGCalHitCalibration.cc:60
HGCalHitCalibration::analyze
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: HGCalHitCalibration.cc:194
TrackingParticle.h
CaloCellGeometry.h
submitPVResolutionJobs.desc
string desc
Definition: submitPVResolutionJobs.py:251
EcalDeadCellBoundaryEnergyFilter_cfi.recHitsEE
recHitsEE
Definition: EcalDeadCellBoundaryEnergyFilter_cfi.py:18
HGCalDetId.h
DetId::rawId
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
SimCluster::simEnergy
float simEnergy() const
returns the accumulated sim energy in the cluster
Definition: SimCluster.h:213
pwdgSkimBPark_cfi.electrons
electrons
Definition: pwdgSkimBPark_cfi.py:6
HcalEndcap
Definition: HcalAssistant.h:34
Frameworkfwd.h
hgcal::RecHitTools::setGeometry
void setGeometry(CaloGeometry const &)
Definition: RecHitTools.cc:68
HGCalHitCalibration::~HGCalHitCalibration
~HGCalHitCalibration() override
Definition: HGCalHitCalibration.cc:120
DetId::HGCalHSc
Definition: DetId.h:34
HGCalHitCalibration::layers_
static constexpr int layers_
Definition: HGCalHitCalibration.cc:79
CaloGeometry.h
HGCalHitCalibration::depletion2_
static constexpr int depletion2_
Definition: HGCalHitCalibration.cc:70
reco::PFCluster
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:42
HGCalHitCalibration::caloParticles_
edm::EDGetTokenT< std::vector< CaloParticle > > caloParticles_
Definition: HGCalHitCalibration.cc:59
or
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
HGCalHitCalibration::caloGeomToken_
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > caloGeomToken_
Definition: HGCalHitCalibration.cc:65
dqm::implementation::IBooker
Definition: DQMStore.h:43
HGCalHitCalibration
Definition: HGCalHitCalibration.cc:43
hgcalTestNeighbor_cfi.detector
detector
Definition: hgcalTestNeighbor_cfi.py:6
SimCluster::hits_and_fractions
std::vector< std::pair< uint32_t, float > > hits_and_fractions() const
Returns list of rechit IDs and fractions for this SimCluster.
Definition: SimCluster.h:184
HGCalHitCalibration::recHitsFH_
edm::EDGetTokenT< HGCRecHitCollection > recHitsFH_
Definition: HGCalHitCalibration.cc:57
HGCHEF
Definition: ForwardSubdetector.h:9
HGCalHitCalibration::electrons_
edm::EDGetTokenT< std::vector< reco::GsfElectron > > electrons_
Definition: HGCalHitCalibration.cc:61
ParameterSet.h
DetId::Forward
Definition: DetId.h:30
edm::HandleBase::isValid
bool isValid() const
Definition: HandleBase.h:70
DeDxTools::esConsumes
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
edm::Event
Definition: Event.h:73
TauDecayModes.dec
dec
Definition: TauDecayModes.py:142
edm::InputTag
Definition: InputTag.h:15
HGCalHitCalibration::bookHistograms
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Definition: HGCalHitCalibration.cc:125
dqm::implementation::IBooker::book1D
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
HGCHEB
Definition: ForwardSubdetector.h:10
HGCalHitCalibration::recHitTools_
hgcal::RecHitTools recHitTools_
Definition: HGCalHitCalibration.cc:68