CMS 3D CMS Logo

HGCalHitCalibration.cc
Go to the documentation of this file.
1 // user include files
2 
6 
10 
12 
14 
22 
27 
30 
32 
33 #include <map>
34 #include <array>
35 #include <string>
36 #include <numeric>
37 
39 public:
40  explicit HGCalHitCalibration(const edm::ParameterSet&);
41  ~HGCalHitCalibration() override;
42 
43  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
44 
45 private:
46  void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
47  void analyze(const edm::Event&, const edm::EventSetup&) override;
48 
49  void fillWithRecHits(std::map<DetId, const HGCRecHit*>&, DetId, unsigned int, float, int&, float&);
50 
58 
62  int debug_;
64  static constexpr int depletion1_ = 200;
65  static constexpr int depletion2_ = 300;
66  static constexpr int scint_ = 400;
67 
68  std::map<int, MonitorElement*> h_EoP_CPene_calib_fraction_;
69  std::map<int, MonitorElement*> hgcal_EoP_CPene_calib_fraction_;
70  std::map<int, MonitorElement*> hgcal_ele_EoP_CPene_calib_fraction_;
71  std::map<int, MonitorElement*> hgcal_photon_EoP_CPene_calib_fraction_;
73 
74  static constexpr int layers_ = 60;
75  std::array<float, layers_> Energy_layer_calib_;
76  std::array<float, layers_> Energy_layer_calib_fraction_;
77 };
78 
80  : caloGeomToken_(esConsumes<CaloGeometry, CaloGeometryRecord>()),
81  rawRecHits_(iConfig.getParameter<bool>("rawRecHits")),
82  debug_(iConfig.getParameter<int>("debug")) {
83  auto detector = iConfig.getParameter<std::string>("detector");
84  auto recHitsEE = iConfig.getParameter<edm::InputTag>("recHitsEE");
85  auto recHitsFH = iConfig.getParameter<edm::InputTag>("recHitsFH");
86  auto recHitsBH = iConfig.getParameter<edm::InputTag>("recHitsBH");
87  auto caloParticles = iConfig.getParameter<edm::InputTag>("caloParticles");
88  auto hgcalMultiClusters = iConfig.getParameter<edm::InputTag>("hgcalMultiClusters");
89  auto electrons = iConfig.getParameter<edm::InputTag>("electrons");
90  auto photons = iConfig.getParameter<edm::InputTag>("photons");
91  if (detector == "all") {
92  recHitsEE_ = consumes<HGCRecHitCollection>(recHitsEE);
93  recHitsFH_ = consumes<HGCRecHitCollection>(recHitsFH);
94  recHitsBH_ = consumes<HGCRecHitCollection>(recHitsBH);
95  algo_ = 1;
96  } else if (detector == "EM") {
97  recHitsEE_ = consumes<HGCRecHitCollection>(recHitsEE);
98  algo_ = 2;
99  } else if (detector == "HAD") {
100  recHitsFH_ = consumes<HGCRecHitCollection>(recHitsFH);
101  recHitsBH_ = consumes<HGCRecHitCollection>(recHitsBH);
102  algo_ = 3;
103  }
104  depletion0_ = iConfig.getParameter<int>("depletionFine");
105  caloParticles_ = consumes<std::vector<CaloParticle> >(caloParticles);
106  hgcalMultiClusters_ = consumes<std::vector<reco::PFCluster> >(hgcalMultiClusters);
107  electrons_ = consumes<std::vector<reco::GsfElectron> >(electrons);
108  photons_ = consumes<std::vector<reco::Photon> >(photons);
109 
110  // size should be HGC layers 52 is enough
111  Energy_layer_calib_.fill(0.);
113 }
114 
116  // do anything here that needs to be done at desctruction time
117  // (e.g. close files, deallocate resources etc.)
118 }
119 
121  edm::Run const& iRun,
122  edm::EventSetup const& /* iSetup */) {
123  ibooker.cd();
124  ibooker.setCurrentFolder("HGCalHitCalibration");
125  h_EoP_CPene_calib_fraction_[depletion0_] = ibooker.book1D("h_EoP_CPene_100_calib_fraction", "", 1000, -0.5, 2.5);
126  h_EoP_CPene_calib_fraction_[depletion1_] = ibooker.book1D("h_EoP_CPene_200_calib_fraction", "", 1000, -0.5, 2.5);
127  h_EoP_CPene_calib_fraction_[depletion2_] = ibooker.book1D("h_EoP_CPene_300_calib_fraction", "", 1000, -0.5, 2.5);
128  h_EoP_CPene_calib_fraction_[scint_] = ibooker.book1D("h_EoP_CPene_scint_calib_fraction", "", 1000, -0.5, 2.5);
130  ibooker.book1D("hgcal_EoP_CPene_100_calib_fraction", "", 1000, -0.5, 2.5);
132  ibooker.book1D("hgcal_EoP_CPene_200_calib_fraction", "", 1000, -0.5, 2.5);
134  ibooker.book1D("hgcal_EoP_CPene_300_calib_fraction", "", 1000, -0.5, 2.5);
135  hgcal_EoP_CPene_calib_fraction_[scint_] = ibooker.book1D("hgcal_EoP_CPene_scint_calib_fraction", "", 1000, -0.5, 2.5);
137  ibooker.book1D("hgcal_ele_EoP_CPene_100_calib_fraction", "", 1000, -0.5, 2.5);
139  ibooker.book1D("hgcal_ele_EoP_CPene_200_calib_fraction", "", 1000, -0.5, 2.5);
141  ibooker.book1D("hgcal_ele_EoP_CPene_300_calib_fraction", "", 1000, -0.5, 2.5);
143  ibooker.book1D("hgcal_ele_EoP_CPene_scint_calib_fraction", "", 1000, -0.5, 2.5);
145  ibooker.book1D("hgcal_photon_EoP_CPene_100_calib_fraction", "", 1000, -0.5, 2.5);
147  ibooker.book1D("hgcal_photon_EoP_CPene_200_calib_fraction", "", 1000, -0.5, 2.5);
149  ibooker.book1D("hgcal_photon_EoP_CPene_300_calib_fraction", "", 1000, -0.5, 2.5);
151  ibooker.book1D("hgcal_photon_EoP_CPene_scint_calib_fraction", "", 1000, -0.5, 2.5);
152  LayerOccupancy_ = ibooker.book1D("LayerOccupancy", "", layers_, 0., (float)layers_);
153 }
154 
155 void HGCalHitCalibration::fillWithRecHits(std::map<DetId, const HGCRecHit*>& hitmap,
156  const DetId hitid,
157  const unsigned int hitlayer,
158  const float fraction,
159  int& seedDet,
160  float& seedEnergy) {
161  if (hitmap.find(hitid) == hitmap.end()) {
162  // Hit was not reconstructed
163  IfLogTrace(debug_ > 0, "HGCalHitCalibration")
164  << ">>> Failed to find detid " << std::hex << hitid.rawId() << std::dec << " Det " << hitid.det() << " Subdet "
165  << hitid.subdetId() << std::endl;
166  return;
167  }
168  if ((hitid.det() != DetId::Forward) && (hitid.det() != DetId::HGCalEE) && (hitid.det() != DetId::HGCalHSi) &&
169  (hitid.det() != DetId::HGCalHSc)) {
170  IfLogTrace(debug_ > 0, "HGCalHitCalibration")
171  << ">>> Wrong type of detid " << std::hex << hitid.rawId() << std::dec << " Det " << hitid.det() << " Subdet "
172  << hitid.subdetId() << std::endl;
173  return;
174  }
175 
176  unsigned int layer = recHitTools_.getLayerWithOffset(hitid);
177  assert(hitlayer == layer);
178  Energy_layer_calib_fraction_[layer] += hitmap[hitid]->energy() * fraction;
180  if (seedEnergy < hitmap[hitid]->energy()) {
181  seedEnergy = hitmap[hitid]->energy();
182  seedDet = recHitTools_.getSiThickness(hitid);
183  if (hitid.det() == DetId::HGCalHSc) {
184  seedDet = scint_;
185  }
186  }
187 }
188 
192 
193  const edm::Handle<std::vector<CaloParticle> >& caloParticleHandle = iEvent.getHandle(caloParticles_);
194  const std::vector<CaloParticle>& caloParticles = *caloParticleHandle;
195 
196  const edm::Handle<std::vector<reco::PFCluster> >& hgcalMultiClustersHandle = iEvent.getHandle(hgcalMultiClusters_);
197 
198  const edm::Handle<std::vector<reco::GsfElectron> >& PFElectronHandle = iEvent.getHandle(electrons_);
199 
200  const edm::Handle<std::vector<reco::Photon> >& PFPhotonHandle = iEvent.getHandle(photons_);
201 
202  // make a map detid-rechit
203  std::map<DetId, const HGCRecHit*> hitmap;
204  switch (algo_) {
205  case 1: {
206  const auto& rechitsEE = iEvent.get(recHitsEE_);
207  const auto& rechitsFH = iEvent.get(recHitsFH_);
208  const auto& rechitsBH = iEvent.get(recHitsBH_);
209  for (unsigned int i = 0; i < rechitsEE.size(); ++i) {
210  hitmap[rechitsEE[i].detid()] = &rechitsEE[i];
211  }
212  for (unsigned int i = 0; i < rechitsFH.size(); ++i) {
213  hitmap[rechitsFH[i].detid()] = &rechitsFH[i];
214  }
215  for (unsigned int i = 0; i < rechitsBH.size(); ++i) {
216  hitmap[rechitsBH[i].detid()] = &rechitsBH[i];
217  }
218  break;
219  }
220  case 2: {
221  const HGCRecHitCollection& rechitsEE = iEvent.get(recHitsEE_);
222  for (unsigned int i = 0; i < rechitsEE.size(); i++) {
223  hitmap[rechitsEE[i].detid()] = &rechitsEE[i];
224  }
225  break;
226  }
227  case 3: {
228  const auto& rechitsFH = iEvent.get(recHitsFH_);
229  const auto& rechitsBH = iEvent.get(recHitsBH_);
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  default:
239  assert(false);
240  break;
241  }
242 
243  // loop over caloParticles
244  int seedDet = 0;
245  float seedEnergy = 0.;
246  IfLogTrace(debug_ > 0, "HGCalHitCalibration") << "Number of caloParticles: " << caloParticles.size() << std::endl;
247  for (const auto& it_caloPart : caloParticles) {
248  if (it_caloPart.g4Tracks()[0].eventId().event() != 0 or it_caloPart.g4Tracks()[0].eventId().bunchCrossing() != 0) {
249  LogDebug("HGCalHitCalibration") << "Excluding CaloParticles from event: "
250  << it_caloPart.g4Tracks()[0].eventId().event()
251  << " with BX: " << it_caloPart.g4Tracks()[0].eventId().bunchCrossing()
252  << std::endl;
253  continue;
254  }
255 
256  const SimClusterRefVector& simClusterRefVector = it_caloPart.simClusters();
257  Energy_layer_calib_.fill(0.);
259 
260  seedDet = 0;
261  seedEnergy = 0.;
262  for (const auto& it_sc : simClusterRefVector) {
263  const SimCluster& simCluster = (*(it_sc));
264  IfLogTrace(debug_ > 1, "HGCalHitCalibration")
265  << ">>> SC.energy(): " << simCluster.energy() << " SC.simEnergy(): " << simCluster.simEnergy() << std::endl;
266  const std::vector<std::pair<uint32_t, float> >& hits_and_fractions = simCluster.hits_and_fractions();
267 
268  // loop over hits
269  for (const auto& it_haf : hits_and_fractions) {
270  unsigned int hitlayer = recHitTools_.getLayerWithOffset(it_haf.first);
271  DetId hitid = (it_haf.first);
272  // dump raw RecHits and match
273  if (rawRecHits_) {
274  if ((hitid.det() == DetId::Forward) || (hitid.det() == DetId::HGCalEE) || (hitid.det() == DetId::HGCalHSi) ||
275  (hitid.det() == DetId::HGCalHSc))
276  fillWithRecHits(hitmap, hitid, hitlayer, it_haf.second, seedDet, seedEnergy);
277  }
278  } // end simHit
279  } // end simCluster
280 
281  auto sumCalibRecHitCalib_fraction =
282  std::accumulate(Energy_layer_calib_fraction_.begin(), Energy_layer_calib_fraction_.end(), 0.);
283  IfLogTrace(debug_ > 0, "HGCalHitCalibration")
284  << ">>> MC Energy: " << it_caloPart.energy() << " reco energy: " << sumCalibRecHitCalib_fraction << std::endl;
285  if (h_EoP_CPene_calib_fraction_.count(seedDet))
286  h_EoP_CPene_calib_fraction_[seedDet]->Fill(sumCalibRecHitCalib_fraction / it_caloPart.energy());
287 
288  // Loop on reconstructed SC.
289  const auto& clusters = *hgcalMultiClustersHandle;
290  float total_energy = 0.;
291  float max_dR2 = 0.0025;
292  auto closest =
293  std::min_element(clusters.begin(), clusters.end(), [&](const reco::PFCluster& a, const reco::PFCluster& b) {
294  auto dR2_a = reco::deltaR2(it_caloPart, a);
295  auto dR2_b = reco::deltaR2(it_caloPart, b);
296  auto ERatio_a = a.correctedEnergy() / it_caloPart.energy();
297  auto ERatio_b = b.correctedEnergy() / it_caloPart.energy();
298  // If both clusters are within 0.0025, mark as first (min) the
299  // element with the highest ratio against the SimCluster
300  if (dR2_a < max_dR2 && dR2_b < max_dR2)
301  return ERatio_a > ERatio_b;
302  return dR2_a < dR2_b;
303  });
304  if (closest != clusters.end() && reco::deltaR2(*closest, it_caloPart) < 0.01) {
305  total_energy = closest->correctedEnergy();
306  seedDet = recHitTools_.getSiThickness(closest->seed());
307  if (closest->seed().det() == DetId::HGCalHSc) {
308  seedDet = scint_;
309  }
310  if (hgcal_EoP_CPene_calib_fraction_.count(seedDet)) {
311  hgcal_EoP_CPene_calib_fraction_[seedDet]->Fill(total_energy / it_caloPart.energy());
312  }
313  }
314 
315  auto closest_fcn = [&](auto const& a, auto const& b) {
316  auto dR2_a = reco::deltaR2(it_caloPart, a);
317  auto dR2_b = reco::deltaR2(it_caloPart, b);
318  auto ERatio_a = a.energy() / it_caloPart.energy();
319  auto ERatio_b = b.energy() / it_caloPart.energy();
320  // If both clusters are within 0.0025, mark as first (min) the
321  // element with the highest ratio against the SimCluster
322  if (dR2_a < max_dR2 && dR2_b < max_dR2)
323  return ERatio_a > ERatio_b;
324  return dR2_a < dR2_b;
325  };
326  // ELECTRONS in HGCAL
327  if (PFElectronHandle.isValid()) {
328  auto const& ele = (*PFElectronHandle);
329  auto closest = std::min_element(ele.begin(), ele.end(), closest_fcn);
330  if (closest != ele.end() &&
331  ((closest->superCluster()->seed()->seed().det() == DetId::Forward) ||
332  (closest->superCluster()->seed()->seed().det() == DetId::HGCalEE) ||
333  (closest->superCluster()->seed()->seed().det() == DetId::HGCalHSi)) &&
334  reco::deltaR2(*closest, it_caloPart) < 0.01) {
335  seedDet = recHitTools_.getSiThickness(closest->superCluster()->seed()->seed());
336  if (closest->superCluster()->seed()->seed().det() == DetId::HGCalHSc) {
337  seedDet = scint_;
338  }
339  if (hgcal_ele_EoP_CPene_calib_fraction_.count(seedDet)) {
340  hgcal_ele_EoP_CPene_calib_fraction_[seedDet]->Fill(closest->energy() / it_caloPart.energy());
341  }
342  }
343  }
344 
345  // PHOTONS in HGCAL
346  if (PFPhotonHandle.isValid()) {
347  auto const& photon = (*PFPhotonHandle);
348  auto closest = std::min_element(photon.begin(), photon.end(), closest_fcn);
349  if (closest != photon.end() &&
350  ((closest->superCluster()->seed()->seed().det() == DetId::Forward) ||
351  (closest->superCluster()->seed()->seed().det() == DetId::HGCalEE) ||
352  (closest->superCluster()->seed()->seed().det() == DetId::HGCalHSi)) &&
353  reco::deltaR2(*closest, it_caloPart) < 0.01) {
354  seedDet = recHitTools_.getSiThickness(closest->superCluster()->seed()->seed());
355  if (hgcal_photon_EoP_CPene_calib_fraction_.count(seedDet)) {
356  hgcal_photon_EoP_CPene_calib_fraction_[seedDet]->Fill(closest->energy() / it_caloPart.energy());
357  }
358  }
359  }
360  } // end caloparticle
361 }
362 
363 // ------------ method fills 'descriptions' with the allowed parameters for the
364 // module ------------
367  desc.add<int>("debug", 0);
368  desc.add<bool>("rawRecHits", true);
369  desc.add<std::string>("detector", "all");
370  desc.add<int>("depletionFine", 120);
371  desc.add<edm::InputTag>("caloParticles", edm::InputTag("mix", "MergedCaloTruth"));
372  desc.add<edm::InputTag>("recHitsEE", edm::InputTag("HGCalRecHit", "HGCEERecHits"));
373  desc.add<edm::InputTag>("recHitsFH", edm::InputTag("HGCalRecHit", "HGCHEFRecHits"));
374  desc.add<edm::InputTag>("recHitsBH", edm::InputTag("HGCalRecHit", "HGCHEBRecHits"));
375  desc.add<edm::InputTag>("hgcalMultiClusters", edm::InputTag("particleFlowClusterHGCal"));
376  desc.add<edm::InputTag>("electrons", edm::InputTag("ecalDrivenGsfElectronsHGC"));
377  desc.add<edm::InputTag>("photons", edm::InputTag("photonsHGC"));
378  descriptions.add("hgcalHitCalibrationDefault", desc);
379 }
380 
381 // define this as a plug-in
static constexpr int layers_
edm::EDGetTokenT< std::vector< CaloParticle > > caloParticles_
static constexpr int depletion2_
static constexpr int scint_
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
edm::EDGetTokenT< std::vector< reco::Photon > > photons_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::array< float, layers_ > Energy_layer_calib_fraction_
edm::EDGetTokenT< std::vector< reco::PFCluster > > hgcalMultiClusters_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::map< int, MonitorElement * > hgcal_EoP_CPene_calib_fraction_
edm::EDGetTokenT< HGCRecHitCollection > recHitsFH_
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:42
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
edm::EDGetTokenT< HGCRecHitCollection > recHitsEE_
#define IfLogTrace(cond, cat)
std::map< int, MonitorElement * > hgcal_ele_EoP_CPene_calib_fraction_
size_type size() const
int closest(std::vector< int > const &vec, int value)
assert(be >=bs)
HGCalHitCalibration(const edm::ParameterSet &)
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:46
std::map< int, MonitorElement * > hgcal_photon_EoP_CPene_calib_fraction_
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
float simEnergy() const
returns the accumulated sim energy in the cluster
Definition: SimCluster.h:213
void Fill(long long x)
MonitorElement * LayerOccupancy_
int iEvent
Definition: GenABIO.cc:224
edm::EDGetTokenT< HGCRecHitCollection > recHitsBH_
Monte Carlo truth information used for tracking validation.
Definition: SimCluster.h:29
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
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
void fillWithRecHits(std::map< DetId, const HGCRecHit *> &, DetId, unsigned int, float, int &, float &)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
hgcal::RecHitTools recHitTools_
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Definition: DetId.h:17
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
edm::EDGetTokenT< std::vector< reco::GsfElectron > > electrons_
std::map< int, MonitorElement * > h_EoP_CPene_calib_fraction_
double b
Definition: hdecay.h:120
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool isValid() const
Definition: HandleBase.h:70
void setGeometry(CaloGeometry const &)
Definition: RecHitTools.cc:68
void analyze(const edm::Event &, const edm::EventSetup &) override
float energy() const
Energy. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:104
double a
Definition: hdecay.h:121
std::array< float, layers_ > Energy_layer_calib_
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > caloGeomToken_
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
std::float_t getSiThickness(const DetId &) const
Definition: RecHitTools.cc:179
Definition: Run.h:45
static constexpr int depletion1_
#define LogDebug(id)
unsigned int getLayerWithOffset(const DetId &) const
Definition: RecHitTools.cc:365