CMS 3D CMS Logo

HGCalHitCalibration.cc
Go to the documentation of this file.
1 // user include files
2 
6 
10 
12 
14 
19 
26 
31 
37 
39 
40 #include <map>
41 #include <array>
42 #include <string>
43 #include <numeric>
44 
46  public:
47  explicit HGCalHitCalibration(const edm::ParameterSet&);
48  ~HGCalHitCalibration() override;
49 
50  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
51 
52  private:
54  edm::EventSetup const&) override;
55  void analyze(const edm::Event&, const edm::EventSetup&) override;
56 
57  void fillWithRecHits(std::map<DetId, const HGCRecHit*>&, DetId, unsigned int,
58  float, int&, float&);
59 
67 
68  int algo_;
70  int debug_;
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 const int layers_ = 60;
80  std::array<float, layers_> Energy_layer_calib_;
81  std::array<float, layers_> Energy_layer_calib_fraction_;
82 };
83 
85  : rawRecHits_(iConfig.getParameter<bool>("rawRecHits")),
86  debug_(iConfig.getParameter<int>("debug")) {
87  auto detector = iConfig.getParameter<std::string>("detector");
88  auto recHitsEE = iConfig.getParameter<edm::InputTag>("recHitsEE");
89  auto recHitsFH = iConfig.getParameter<edm::InputTag>("recHitsFH");
90  auto recHitsBH = iConfig.getParameter<edm::InputTag>("recHitsBH");
91  auto caloParticles = iConfig.getParameter<edm::InputTag>("caloParticles");
92  auto hgcalMultiClusters = iConfig.getParameter<edm::InputTag>("hgcalMultiClusters");
93  auto electrons = iConfig.getParameter<edm::InputTag>("electrons");
94  auto photons = iConfig.getParameter<edm::InputTag>("photons");
95  if (detector == "all") {
96  recHitsEE_ = consumes<HGCRecHitCollection>(recHitsEE);
97  recHitsFH_ = consumes<HGCRecHitCollection>(recHitsFH);
98  recHitsBH_ = consumes<HGCRecHitCollection>(recHitsBH);
99  algo_ = 1;
100  } else if (detector == "EM") {
101  recHitsEE_ = consumes<HGCRecHitCollection>(recHitsEE);
102  algo_ = 2;
103  } else if (detector == "HAD") {
104  recHitsFH_ = consumes<HGCRecHitCollection>(recHitsFH);
105  recHitsBH_ = consumes<HGCRecHitCollection>(recHitsBH);
106  algo_ = 3;
107  }
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("HGCalHitCalibration");
129  ibooker.book1D("h_EoP_CPene_100_calib_fraction", "", 1000, -0.5, 2.5);
131  ibooker.book1D("h_EoP_CPene_200_calib_fraction", "", 1000, -0.5, 2.5);
133  ibooker.book1D("h_EoP_CPene_300_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);
141  ibooker.book1D("hgcal_ele_EoP_CPene_100_calib_fraction", "", 1000, -0.5, 2.5);
143  ibooker.book1D("hgcal_ele_EoP_CPene_200_calib_fraction", "", 1000, -0.5, 2.5);
145  ibooker.book1D("hgcal_ele_EoP_CPene_300_calib_fraction", "", 1000, -0.5, 2.5);
147  ibooker.book1D("hgcal_photon_EoP_CPene_100_calib_fraction", "", 1000, -0.5, 2.5);
149  ibooker.book1D("hgcal_photon_EoP_CPene_200_calib_fraction", "", 1000, -0.5, 2.5);
151  ibooker.book1D("hgcal_photon_EoP_CPene_300_calib_fraction", "", 1000, -0.5, 2.5);
152  LayerOccupancy_ = ibooker.book1D("LayerOccupancy", "", layers_, 0., (float)layers_);
153 }
154 
156  std::map<DetId, const HGCRecHit*>& hitmap, const DetId hitid,
157  const unsigned int hitlayer, const float fraction, int& seedDet,
158  float& seedEnergy) {
159  if (hitmap.find(hitid) == hitmap.end()) {
160  // Hit was not reconstructed
161  IfLogTrace(debug_ > 0, "HGCalHitCalibration")
162  << ">>> Failed to find detid " << hitid.rawId()
163  << " Det " << hitid.det()
164  << " Subdet " << hitid.subdetId() << std::endl;
165  return;
166  }
167  unsigned int layer = recHitTools_.getLayerWithOffset(hitid);
168  assert(hitlayer == layer);
169  Energy_layer_calib_fraction_[layer] += hitmap[hitid]->energy() * fraction;
170  LayerOccupancy_->Fill(layer);
171  if (seedEnergy < hitmap[hitid]->energy()) {
172  seedEnergy = hitmap[hitid]->energy();
173  seedDet = recHitTools_.getSiThickness(hitid);
174  }
175 }
176 
178  const edm::EventSetup& iSetup) {
179  using namespace edm;
180 
181  recHitTools_.getEventSetup(iSetup);
182 
183  Handle<HGCRecHitCollection> recHitHandleEE;
184  Handle<HGCRecHitCollection> recHitHandleFH;
185  Handle<HGCRecHitCollection> recHitHandleBH;
186 
187  Handle<std::vector<CaloParticle> > caloParticleHandle;
188  iEvent.getByToken(caloParticles_, caloParticleHandle);
189  const std::vector<CaloParticle>& caloParticles = *caloParticleHandle;
190 
191  Handle<std::vector<reco::PFCluster> > hgcalMultiClustersHandle;
192  iEvent.getByToken(hgcalMultiClusters_, hgcalMultiClustersHandle);
193 
194  Handle<std::vector<reco::GsfElectron> > PFElectronHandle;
195  iEvent.getByToken(electrons_, PFElectronHandle);
196 
197  Handle<std::vector<reco::Photon> > PFPhotonHandle;
198  iEvent.getByToken(photons_, PFPhotonHandle);
199 
200  // make a map detid-rechit
201  std::map<DetId, const HGCRecHit*> hitmap;
202  switch (algo_) {
203  case 1: {
204  iEvent.getByToken(recHitsEE_, recHitHandleEE);
205  iEvent.getByToken(recHitsFH_, recHitHandleFH);
206  iEvent.getByToken(recHitsBH_, recHitHandleBH);
207  const auto& rechitsEE = *recHitHandleEE;
208  const auto& rechitsFH = *recHitHandleFH;
209  const auto& rechitsBH = *recHitHandleBH;
210  for (unsigned int i = 0; i < rechitsEE.size(); ++i) {
211  hitmap[rechitsEE[i].detid()] = &rechitsEE[i];
212  }
213  for (unsigned int i = 0; i < rechitsFH.size(); ++i) {
214  hitmap[rechitsFH[i].detid()] = &rechitsFH[i];
215  }
216  for (unsigned int i = 0; i < rechitsBH.size(); ++i) {
217  hitmap[rechitsBH[i].detid()] = &rechitsBH[i];
218  }
219  break;
220  }
221  case 2: {
222  iEvent.getByToken(recHitsEE_, recHitHandleEE);
223  const HGCRecHitCollection& rechitsEE = *recHitHandleEE;
224  for (unsigned int i = 0; i < rechitsEE.size(); i++) {
225  hitmap[rechitsEE[i].detid()] = &rechitsEE[i];
226  }
227  break;
228  }
229  case 3: {
230  iEvent.getByToken(recHitsFH_, recHitHandleFH);
231  iEvent.getByToken(recHitsBH_, recHitHandleBH);
232  const auto& rechitsFH = *recHitHandleFH;
233  const auto& rechitsBH = *recHitHandleBH;
234  for (unsigned int i = 0; i < rechitsFH.size(); i++) {
235  hitmap[rechitsFH[i].detid()] = &rechitsFH[i];
236  }
237  for (unsigned int i = 0; i < rechitsBH.size(); i++) {
238  hitmap[rechitsBH[i].detid()] = &rechitsBH[i];
239  }
240  break;
241  }
242  default:
243  assert(false);
244  break;
245  }
246 
247  // loop over caloParticles
248  int seedDet = 0;
249  float seedEnergy = 0.;
250  IfLogTrace(debug_ > 0, "HGCalHitCalibration")
251  << "Number of caloParticles: " << caloParticles.size() << std::endl;
252  for (const auto& it_caloPart : caloParticles) {
253  const SimClusterRefVector& simClusterRefVector = it_caloPart.simClusters();
254  Energy_layer_calib_.fill(0.);
256 
257  seedDet = 0;
258  seedEnergy = 0.;
259  for (const auto& it_sc : simClusterRefVector) {
260  const SimCluster& simCluster = (*(it_sc));
261  IfLogTrace(debug_ > 1, "HGCalHitCalibration")
262  << ">>> SC.energy(): " << simCluster.energy()
263  << " SC.simEnergy(): " << simCluster.simEnergy()
264  << std::endl;
265  const std::vector<std::pair<uint32_t, float> >& hits_and_fractions =
266  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 &&
275  (hitid.subdetId() == HGCEE or hitid.subdetId() == HGCHEF or
276  hitid.subdetId() == HGCHEB)) ||
277  (hitid.det() == DetId::Hcal && hitid.subdetId() == HcalEndcap))
278  fillWithRecHits(hitmap, hitid, hitlayer, it_haf.second, seedDet,
279  seedEnergy);
280  }
281  } // end simHit
282  } // end simCluster
283 
284  auto sumCalibRecHitCalib_fraction = std::accumulate(Energy_layer_calib_fraction_.begin(),
285  Energy_layer_calib_fraction_.end(), 0.);
286  IfLogTrace(debug_ > 0, "HGCalHitCalibration")
287  << ">>> MC Energy: " << it_caloPart.energy()
288  << " reco energy: " << sumCalibRecHitCalib_fraction
289  << std::endl;
290  if (h_EoP_CPene_calib_fraction_.count(seedDet))
291  h_EoP_CPene_calib_fraction_[seedDet]->Fill(sumCalibRecHitCalib_fraction /
292  it_caloPart.energy());
293 
294  // Loop on reconstructed SC.
295  const auto& clusters = *hgcalMultiClustersHandle;
296  float total_energy = 0.;
297  float max_dR2 = 0.0025;
298  auto closest = std::min_element(clusters.begin(),
299  clusters.end(),
300  [&](const reco::PFCluster &a,
301  const reco::PFCluster &b) {
302  auto dR2_a = reco::deltaR2(it_caloPart, a);
303  auto dR2_b = reco::deltaR2(it_caloPart, b);
304  auto ERatio_a = a.correctedEnergy()/it_caloPart.energy();
305  auto ERatio_b = b.correctedEnergy()/it_caloPart.energy();
306  // If both clusters are within 0.0025, mark as first (min) the
307  // element with the highest ratio against the SimCluster
308  if (dR2_a < max_dR2 && dR2_b < max_dR2)
309  return ERatio_a > ERatio_b;
310  return dR2_a < dR2_b;
311  });
312  if (closest != clusters.end()
313  && reco::deltaR2(*closest, it_caloPart) < 0.01) {
314  total_energy = closest->correctedEnergy();
315  seedDet = recHitTools_.getSiThickness(closest->seed());
316  if (hgcal_EoP_CPene_calib_fraction_.count(seedDet)) {
317  hgcal_EoP_CPene_calib_fraction_[seedDet]->Fill(total_energy /
318  it_caloPart.energy());
319  }
320  }
321 
322  auto closest_fcn = [&](auto const &a, auto const&b){
323  auto dR2_a = reco::deltaR2(it_caloPart, a);
324  auto dR2_b = reco::deltaR2(it_caloPart, b);
325  auto ERatio_a = a.energy()/it_caloPart.energy();
326  auto ERatio_b = b.energy()/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  // ELECTRONS in HGCAL
334  if (PFElectronHandle.isValid()) {
335  auto const & ele = (*PFElectronHandle);
336  auto closest = std::min_element(ele.begin(),
337  ele.end(),
338  closest_fcn);
339  if (closest != ele.end()
340  && closest->superCluster()->seed()->seed().det() == DetId::Forward
341  && reco::deltaR2(*closest, it_caloPart) < 0.01) {
342  seedDet = recHitTools_.getSiThickness(closest->superCluster()->seed()->seed());
343  if (hgcal_ele_EoP_CPene_calib_fraction_.count(seedDet)) {
344  hgcal_ele_EoP_CPene_calib_fraction_[seedDet]->Fill(closest->energy() /
345  it_caloPart.energy());
346  }
347  }
348  }
349 
350  // PHOTONS in HGCAL
351  if (PFPhotonHandle.isValid()) {
352  auto const & photon = (*PFPhotonHandle);
353  auto closest = std::min_element(photon.begin(),
354  photon.end(),
355  closest_fcn);
356  if (closest != photon.end()
357  && closest->superCluster()->seed()->seed().det() == DetId::Forward
358  && reco::deltaR2(*closest, it_caloPart) < 0.01) {
359  seedDet = recHitTools_.getSiThickness(closest->superCluster()->seed()->seed());
360  if (hgcal_photon_EoP_CPene_calib_fraction_.count(seedDet)) {
361  hgcal_photon_EoP_CPene_calib_fraction_[seedDet]->Fill(closest->energy() /
362  it_caloPart.energy());
363  }
364  }
365  }
366  } // end caloparticle
367 }
368 
369 // ------------ method fills 'descriptions' with the allowed parameters for the
370 // module ------------
372  edm::ConfigurationDescriptions& descriptions) {
374  desc.add<int>("debug", 0);
375  desc.add<bool>("rawRecHits", true);
376  desc.add<std::string>("detector", "all");
377  desc.add<edm::InputTag>("caloParticles", edm::InputTag("mix", "MergedCaloTruth"));
378  desc.add<edm::InputTag>("recHitsEE", edm::InputTag("HGCalRecHit", "HGCEERecHits"));
379  desc.add<edm::InputTag>("recHitsFH", edm::InputTag("HGCalRecHit", "HGCHEFRecHits"));
380  desc.add<edm::InputTag>("recHitsBH", edm::InputTag("HGCalRecHit", "HGCHEBRecHits"));
381  desc.add<edm::InputTag>("hgcalMultiClusters", edm::InputTag("particleFlowClusterHGCalFromMultiCl"));
382  desc.add<edm::InputTag>("electrons", edm::InputTag("ecalDrivenGsfElectronsFromMultiCl"));
383  desc.add<edm::InputTag>("photons", edm::InputTag("photonsFromMultiCl"));
384  descriptions.add("hgcalHitCalibration", desc);
385 }
386 
387 // define this as a plug-in
T getParameter(std::string const &) const
edm::EDGetTokenT< std::vector< CaloParticle > > caloParticles_
edm::EDGetTokenT< std::vector< reco::Photon > > photons_
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_
double correctedEnergy() const
Definition: CaloCluster.h:125
edm::EDGetTokenT< HGCRecHitCollection > recHitsFH_
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:47
edm::EDGetTokenT< HGCRecHitCollection > recHitsEE_
std::vector< std::pair< uint32_t, float > > hits_and_fractions() const
Returns list of rechit IDs and fractions for this SimCluster.
Definition: SimCluster.h:210
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
std::map< int, MonitorElement * > hgcal_ele_EoP_CPene_calib_fraction_
void cd(void)
Definition: DQMStore.cc:269
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
#define IfLogTrace(cond, cat)
HGCalHitCalibration(const edm::ParameterSet &)
std::map< int, MonitorElement * > hgcal_photon_EoP_CPene_calib_fraction_
void getEventSetup(const edm::EventSetup &)
Definition: RecHitTools.cc:61
void Fill(long long x)
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
MonitorElement * LayerOccupancy_
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< HGCRecHitCollection > recHitsBH_
Monte Carlo truth information used for tracking validation.
Definition: SimCluster.h:28
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
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
float energy() const
Energy. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:111
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:74
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
double energy() const
cluster energy
Definition: PFCluster.h:82
hgcal::RecHitTools recHitTools_
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
float simEnergy() const
returns the accumulated sim energy in the cluster
Definition: SimCluster.h:225
std::float_t getSiThickness(const DetId &) const
Definition: RecHitTools.cc:103
Definition: DetId.h:18
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:277
unsigned int getLayerWithOffset(const DetId &) const
Definition: RecHitTools.cc:176
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)
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36
void analyze(const edm::Event &, const edm::EventSetup &) override
HLT enums.
void fillWithRecHits(std::map< DetId, const HGCRecHit * > &, DetId, unsigned int, float, int &, float &)
size_type size() const
double a
Definition: hdecay.h:121
std::array< float, layers_ > Energy_layer_calib_
Detector det() const
get the detector field from this detid
Definition: DetId.h:35
Definition: Run.h:42
static const int layers_