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 #include <cmath>
38 
40 public:
41  explicit HGCalHitCalibration(const edm::ParameterSet&);
42  ~HGCalHitCalibration() override;
43 
44  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
45 
46 private:
47  void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
48  void analyze(const edm::Event&, const edm::EventSetup&) override;
49 
50  void fillWithRecHits(std::map<DetId, const HGCRecHit*>&, DetId, unsigned int, float, int&, float&);
51 
59 
63  int debug_;
66  static constexpr int depletion1_ = 200;
67  static constexpr int depletion2_ = 300;
68  static constexpr int scint_ = 400;
69 
70  std::map<int, MonitorElement*> h_EoP_CPene_calib_fraction_;
71  std::map<int, MonitorElement*> hgcal_EoP_CPene_calib_fraction_;
72  std::map<int, MonitorElement*> hgcal_ele_EoP_CPene_calib_fraction_;
73  std::map<int, MonitorElement*> hgcal_photon_EoP_CPene_calib_fraction_;
75 
76  static constexpr int layers_ = 60;
77  std::array<float, layers_> Energy_layer_calib_;
78  std::array<float, layers_> Energy_layer_calib_fraction_;
79 };
80 
82  : caloGeomToken_(esConsumes<CaloGeometry, CaloGeometryRecord>()),
83  rawRecHits_(iConfig.getParameter<bool>("rawRecHits")),
84  debug_(iConfig.getParameter<int>("debug")),
85  folder_(iConfig.getParameter<std::string>("folder")) {
86  auto detector = iConfig.getParameter<std::string>("detector");
87  auto recHitsEE = iConfig.getParameter<edm::InputTag>("recHitsEE");
88  auto recHitsFH = iConfig.getParameter<edm::InputTag>("recHitsFH");
89  auto recHitsBH = iConfig.getParameter<edm::InputTag>("recHitsBH");
90  auto caloParticles = iConfig.getParameter<edm::InputTag>("caloParticles");
91  auto hgcalMultiClusters = iConfig.getParameter<edm::InputTag>("hgcalMultiClusters");
92  auto electrons = iConfig.getParameter<edm::InputTag>("electrons");
93  auto photons = iConfig.getParameter<edm::InputTag>("photons");
94  if (detector == "all") {
95  recHitsEE_ = consumes<HGCRecHitCollection>(recHitsEE);
96  recHitsFH_ = consumes<HGCRecHitCollection>(recHitsFH);
97  recHitsBH_ = consumes<HGCRecHitCollection>(recHitsBH);
98  algo_ = 1;
99  } else if (detector == "EM") {
100  recHitsEE_ = consumes<HGCRecHitCollection>(recHitsEE);
101  algo_ = 2;
102  } else if (detector == "HAD") {
103  recHitsFH_ = consumes<HGCRecHitCollection>(recHitsFH);
104  recHitsBH_ = consumes<HGCRecHitCollection>(recHitsBH);
105  algo_ = 3;
106  }
107  depletion0_ = iConfig.getParameter<int>("depletionFine");
108  caloParticles_ = consumes<std::vector<CaloParticle> >(caloParticles);
109  hgcalMultiClusters_ = consumes<std::vector<reco::PFCluster> >(hgcalMultiClusters);
110  electrons_ = consumes<std::vector<reco::GsfElectron> >(electrons);
111  photons_ = consumes<std::vector<reco::Photon> >(photons);
112 
113  // size should be HGC layers 52 is enough
114  Energy_layer_calib_.fill(0.);
116 }
117 
119  // do anything here that needs to be done at desctruction time
120  // (e.g. close files, deallocate resources etc.)
121 }
122 
124  edm::Run const& iRun,
125  edm::EventSetup const& /* iSetup */) {
126  ibooker.cd();
127  ibooker.setCurrentFolder(folder_);
128  h_EoP_CPene_calib_fraction_[depletion0_] = ibooker.book1D("h_EoP_CPene_100_calib_fraction", "", 1000, -0.5, 2.5);
129  h_EoP_CPene_calib_fraction_[depletion1_] = ibooker.book1D("h_EoP_CPene_200_calib_fraction", "", 1000, -0.5, 2.5);
130  h_EoP_CPene_calib_fraction_[depletion2_] = ibooker.book1D("h_EoP_CPene_300_calib_fraction", "", 1000, -0.5, 2.5);
131  h_EoP_CPene_calib_fraction_[scint_] = ibooker.book1D("h_EoP_CPene_scint_calib_fraction", "", 1000, -0.5, 2.5);
133  ibooker.book1D("hgcal_EoP_CPene_100_calib_fraction", "", 1000, -0.5, 2.5);
135  ibooker.book1D("hgcal_EoP_CPene_200_calib_fraction", "", 1000, -0.5, 2.5);
137  ibooker.book1D("hgcal_EoP_CPene_300_calib_fraction", "", 1000, -0.5, 2.5);
138  hgcal_EoP_CPene_calib_fraction_[scint_] = ibooker.book1D("hgcal_EoP_CPene_scint_calib_fraction", "", 1000, -0.5, 2.5);
140  ibooker.book1D("hgcal_ele_EoP_CPene_100_calib_fraction", "", 1000, -0.5, 2.5);
142  ibooker.book1D("hgcal_ele_EoP_CPene_200_calib_fraction", "", 1000, -0.5, 2.5);
144  ibooker.book1D("hgcal_ele_EoP_CPene_300_calib_fraction", "", 1000, -0.5, 2.5);
146  ibooker.book1D("hgcal_ele_EoP_CPene_scint_calib_fraction", "", 1000, -0.5, 2.5);
148  ibooker.book1D("hgcal_photon_EoP_CPene_100_calib_fraction", "", 1000, -0.5, 2.5);
150  ibooker.book1D("hgcal_photon_EoP_CPene_200_calib_fraction", "", 1000, -0.5, 2.5);
152  ibooker.book1D("hgcal_photon_EoP_CPene_300_calib_fraction", "", 1000, -0.5, 2.5);
154  ibooker.book1D("hgcal_photon_EoP_CPene_scint_calib_fraction", "", 1000, -0.5, 2.5);
155  LayerOccupancy_ = ibooker.book1D("LayerOccupancy", "", layers_, 0., (float)layers_);
156 }
157 
158 void HGCalHitCalibration::fillWithRecHits(std::map<DetId, const HGCRecHit*>& hitmap,
159  const DetId hitid,
160  const unsigned int hitlayer,
161  const float fraction,
162  int& seedDet,
163  float& seedEnergy) {
164  if (hitmap.find(hitid) == hitmap.end()) {
165  // Hit was not reconstructed
166  IfLogTrace(debug_ > 0, "HGCalHitCalibration")
167  << ">>> Failed to find detid " << std::hex << hitid.rawId() << std::dec << " Det " << hitid.det() << " Subdet "
168  << hitid.subdetId() << std::endl;
169  return;
170  }
171  if ((hitid.det() != DetId::Forward) && (hitid.det() != DetId::HGCalEE) && (hitid.det() != DetId::HGCalHSi) &&
172  (hitid.det() != DetId::HGCalHSc)) {
173  IfLogTrace(debug_ > 0, "HGCalHitCalibration")
174  << ">>> Wrong type of detid " << std::hex << hitid.rawId() << std::dec << " Det " << hitid.det() << " Subdet "
175  << hitid.subdetId() << std::endl;
176  return;
177  }
178 
179  unsigned int layer = recHitTools_.getLayerWithOffset(hitid);
180  assert(hitlayer == layer);
181  Energy_layer_calib_fraction_[layer] += hitmap[hitid]->energy() * fraction;
183  if (seedEnergy < hitmap[hitid]->energy()) {
184  seedEnergy = hitmap[hitid]->energy();
185  seedDet = std::rint(recHitTools_.getSiThickness(hitid));
186  if (hitid.det() == DetId::HGCalHSc) {
187  seedDet = scint_;
188  }
189  }
190 }
191 
193  float constexpr max_dR2 = 0.0025;
194 
197 
198  const edm::Handle<std::vector<CaloParticle> >& caloParticleHandle = iEvent.getHandle(caloParticles_);
199  const std::vector<CaloParticle>& caloParticles = *caloParticleHandle;
200 
201  const edm::Handle<std::vector<reco::PFCluster> >& hgcalMultiClustersHandle = iEvent.getHandle(hgcalMultiClusters_);
202 
203  const edm::Handle<std::vector<reco::GsfElectron> >& PFElectronHandle = iEvent.getHandle(electrons_);
204 
205  const edm::Handle<std::vector<reco::Photon> >& PFPhotonHandle = iEvent.getHandle(photons_);
206 
207  // make a map detid-rechit
208  std::map<DetId, const HGCRecHit*> hitmap;
209  switch (algo_) {
210  case 1: {
211  const auto& rechitsEE_handle = iEvent.getHandle(recHitsEE_);
212  if (!rechitsEE_handle.isValid())
213  return;
214 
215  const auto& rechitsFH_handle = iEvent.getHandle(recHitsFH_);
216  if (!rechitsFH_handle.isValid())
217  return;
218 
219  const auto& rechitsBH_handle = iEvent.getHandle(recHitsBH_);
220  if (!rechitsBH_handle.isValid())
221  return;
222 
223  auto const& rechitsEE = *rechitsEE_handle;
224  auto const& rechitsFH = *rechitsFH_handle;
225  auto const& rechitsBH = *rechitsBH_handle;
226  for (unsigned int i = 0; i < rechitsEE.size(); ++i) {
227  hitmap[rechitsEE[i].detid()] = &rechitsEE[i];
228  }
229  for (unsigned int i = 0; i < rechitsFH.size(); ++i) {
230  hitmap[rechitsFH[i].detid()] = &rechitsFH[i];
231  }
232  for (unsigned int i = 0; i < rechitsBH.size(); ++i) {
233  hitmap[rechitsBH[i].detid()] = &rechitsBH[i];
234  }
235  break;
236  }
237  case 2: {
238  const auto& rechitsEE_handle = iEvent.getHandle(recHitsEE_);
239  if (!rechitsEE_handle.isValid())
240  return;
241 
242  auto const& rechitsEE = *rechitsEE_handle;
243  for (unsigned int i = 0; i < rechitsEE.size(); i++) {
244  hitmap[rechitsEE[i].detid()] = &rechitsEE[i];
245  }
246  break;
247  }
248  case 3: {
249  const auto& rechitsFH_handle = iEvent.getHandle(recHitsFH_);
250  if (!rechitsFH_handle.isValid())
251  return;
252 
253  const auto& rechitsBH_handle = iEvent.getHandle(recHitsBH_);
254  if (!rechitsBH_handle.isValid())
255  return;
256 
257  auto const& rechitsFH = *rechitsFH_handle;
258  auto const& rechitsBH = *rechitsBH_handle;
259  for (unsigned int i = 0; i < rechitsFH.size(); i++) {
260  hitmap[rechitsFH[i].detid()] = &rechitsFH[i];
261  }
262  for (unsigned int i = 0; i < rechitsBH.size(); i++) {
263  hitmap[rechitsBH[i].detid()] = &rechitsBH[i];
264  }
265  break;
266  }
267  default:
268  assert(false);
269  break;
270  }
271 
272  // loop over caloParticles
273  int seedDet = 0;
274  float seedEnergy = 0.;
275  IfLogTrace(debug_ > 0, "HGCalHitCalibration") << "Number of caloParticles: " << caloParticles.size() << std::endl;
276  for (const auto& it_caloPart : caloParticles) {
277  if (it_caloPart.g4Tracks()[0].eventId().event() != 0 or it_caloPart.g4Tracks()[0].eventId().bunchCrossing() != 0) {
278  LogDebug("HGCalHitCalibration") << "Excluding CaloParticles from event: "
279  << it_caloPart.g4Tracks()[0].eventId().event()
280  << " with BX: " << it_caloPart.g4Tracks()[0].eventId().bunchCrossing()
281  << std::endl;
282  continue;
283  }
284 
285  const SimClusterRefVector& simClusterRefVector = it_caloPart.simClusters();
286  Energy_layer_calib_.fill(0.);
288 
289  seedDet = 0;
290  seedEnergy = 0.;
291  for (const auto& it_sc : simClusterRefVector) {
292  const SimCluster& simCluster = (*(it_sc));
293  IfLogTrace(debug_ > 1, "HGCalHitCalibration")
294  << ">>> SC.energy(): " << simCluster.energy() << " SC.simEnergy(): " << simCluster.simEnergy() << std::endl;
295  const std::vector<std::pair<uint32_t, float> >& hits_and_fractions = simCluster.hits_and_fractions();
296 
297  // loop over hits
298  for (const auto& it_haf : hits_and_fractions) {
299  unsigned int hitlayer = recHitTools_.getLayerWithOffset(it_haf.first);
300  DetId hitid = (it_haf.first);
301  // dump raw RecHits and match
302  if (rawRecHits_) {
303  if ((hitid.det() == DetId::Forward) || (hitid.det() == DetId::HGCalEE) || (hitid.det() == DetId::HGCalHSi) ||
304  (hitid.det() == DetId::HGCalHSc))
305  fillWithRecHits(hitmap, hitid, hitlayer, it_haf.second, seedDet, seedEnergy);
306  }
307  } // end simHit
308  } // end simCluster
309 
310  auto sumCalibRecHitCalib_fraction =
311  std::accumulate(Energy_layer_calib_fraction_.begin(), Energy_layer_calib_fraction_.end(), 0.);
312  IfLogTrace(debug_ > 0, "HGCalHitCalibration")
313  << ">>> MC Energy: " << it_caloPart.energy() << " reco energy: " << sumCalibRecHitCalib_fraction << std::endl;
314  if (h_EoP_CPene_calib_fraction_.count(seedDet))
315  h_EoP_CPene_calib_fraction_[seedDet]->Fill(sumCalibRecHitCalib_fraction / it_caloPart.energy());
316 
317  // Loop on reconstructed SC.
318  if (hgcalMultiClustersHandle.isValid()) {
319  const auto& clusters = *hgcalMultiClustersHandle;
320  float total_energy = 0.;
321  auto closest =
322  std::min_element(clusters.begin(), clusters.end(), [&](const reco::PFCluster& a, const reco::PFCluster& b) {
323  auto dR2_a = reco::deltaR2(it_caloPart, a);
324  auto dR2_b = reco::deltaR2(it_caloPart, b);
325  auto ERatio_a = a.correctedEnergy() / it_caloPart.energy();
326  auto ERatio_b = b.correctedEnergy() / it_caloPart.energy();
327  // If both clusters are within 0.0025, mark as first (min) the
328  // element with the highest ratio against the SimCluster
329  if (dR2_a < max_dR2 && dR2_b < max_dR2)
330  return ERatio_a > ERatio_b;
331  return dR2_a < dR2_b;
332  });
333  if (closest != clusters.end() && reco::deltaR2(*closest, it_caloPart) < 0.01) {
334  total_energy = closest->correctedEnergy();
335  seedDet = recHitTools_.getSiThickness(closest->seed());
336  if (closest->seed().det() == DetId::HGCalHSc) {
337  seedDet = scint_;
338  }
339  if (hgcal_EoP_CPene_calib_fraction_.count(seedDet)) {
340  hgcal_EoP_CPene_calib_fraction_[seedDet]->Fill(total_energy / it_caloPart.energy());
341  }
342  }
343  }
344 
345  auto closest_fcn = [&](auto const& a, auto const& b) {
346  auto dR2_a = reco::deltaR2(it_caloPart, a);
347  auto dR2_b = reco::deltaR2(it_caloPart, b);
348  auto ERatio_a = a.energy() / it_caloPart.energy();
349  auto ERatio_b = b.energy() / it_caloPart.energy();
350  // If both clusters are within 0.0025, mark as first (min) the
351  // element with the highest ratio against the SimCluster
352  if (dR2_a < max_dR2 && dR2_b < max_dR2)
353  return ERatio_a > ERatio_b;
354  return dR2_a < dR2_b;
355  };
356 
357  // ELECTRONS in HGCAL
358  if (PFElectronHandle.isValid()) {
359  auto const& ele = (*PFElectronHandle);
360  auto closest = std::min_element(ele.begin(), ele.end(), closest_fcn);
361  if (closest != ele.end() &&
362  ((closest->superCluster()->seed()->seed().det() == DetId::Forward) ||
363  (closest->superCluster()->seed()->seed().det() == DetId::HGCalEE) ||
364  (closest->superCluster()->seed()->seed().det() == DetId::HGCalHSi)) &&
365  reco::deltaR2(*closest, it_caloPart) < 0.01) {
366  seedDet = recHitTools_.getSiThickness(closest->superCluster()->seed()->seed());
367  if (closest->superCluster()->seed()->seed().det() == DetId::HGCalHSc) {
368  seedDet = scint_;
369  }
370  if (hgcal_ele_EoP_CPene_calib_fraction_.count(seedDet)) {
371  hgcal_ele_EoP_CPene_calib_fraction_[seedDet]->Fill(closest->energy() / it_caloPart.energy());
372  }
373  }
374  }
375 
376  // PHOTONS in HGCAL
377  if (PFPhotonHandle.isValid()) {
378  auto const& photon = (*PFPhotonHandle);
379  auto closest = std::min_element(photon.begin(), photon.end(), closest_fcn);
380  if (closest != photon.end() &&
381  ((closest->superCluster()->seed()->seed().det() == DetId::Forward) ||
382  (closest->superCluster()->seed()->seed().det() == DetId::HGCalEE) ||
383  (closest->superCluster()->seed()->seed().det() == DetId::HGCalHSi)) &&
384  reco::deltaR2(*closest, it_caloPart) < 0.01) {
385  seedDet = recHitTools_.getSiThickness(closest->superCluster()->seed()->seed());
386  if (hgcal_photon_EoP_CPene_calib_fraction_.count(seedDet)) {
387  hgcal_photon_EoP_CPene_calib_fraction_[seedDet]->Fill(closest->energy() / it_caloPart.energy());
388  }
389  }
390  }
391  } // end caloparticle
392 }
393 
394 // ------------ method fills 'descriptions' with the allowed parameters for the
395 // module ------------
398  desc.add<int>("debug", 0);
399  desc.add<bool>("rawRecHits", true);
400  desc.add<std::string>("folder", "HGCalHitCalibration");
401  desc.add<std::string>("detector", "all");
402  desc.add<int>("depletionFine", 120);
403  desc.add<edm::InputTag>("caloParticles", edm::InputTag("mix", "MergedCaloTruth"));
404  desc.add<edm::InputTag>("recHitsEE", edm::InputTag("HGCalRecHit", "HGCEERecHits"));
405  desc.add<edm::InputTag>("recHitsFH", edm::InputTag("HGCalRecHit", "HGCHEFRecHits"));
406  desc.add<edm::InputTag>("recHitsBH", edm::InputTag("HGCalRecHit", "HGCHEBRecHits"));
407  desc.add<edm::InputTag>("hgcalMultiClusters", edm::InputTag("particleFlowClusterHGCal"));
408  desc.add<edm::InputTag>("electrons", edm::InputTag("ecalDrivenGsfElectronsHGC"));
409  desc.add<edm::InputTag>("photons", edm::InputTag("photonsHGC"));
410  descriptions.add("hgcalHitCalibrationDefault", desc);
411 }
412 
413 // 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:307
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_
std::string folder_
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_
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:188
float simEnergy() const
returns the accumulated sim energy in the cluster
Definition: SimCluster.h:241
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:33
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:79
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:108
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:190
Definition: Run.h:45
static constexpr int depletion1_
#define LogDebug(id)
unsigned int getLayerWithOffset(const DetId &) const
Definition: RecHitTools.cc:381