CMS 3D CMS Logo

EnergyScaleAnalyzer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: EnergyScaleAnalyzer
4 // Class: EnergyScaleAnalyzer
5 //
14 // Original Author: Keti Kaadze
15 // Created: Thu Jun 21 08:59:42 CDT 2007
16 //
17 
19 
20 // Framework
30 
31 // Geometry
38 
39 #include "TFile.h"
40 // Reconstruction classes
49 
53 
60 
61 #include "TString.h"
63 #include <fstream>
64 #include <iomanip>
65 #include <ios>
66 #include <iostream>
67 #include <map>
68 
69 //========================================================================
71 //========================================================================
72 {
73  hepMCLabel_ = consumes<edm::HepMCProduct>(ps.getParameter<std::string>("hepMCLabel"));
74  hybridSuperClusters_token = consumes<reco::SuperClusterCollection>(
75  ps.getUntrackedParameter<edm::InputTag>("hybridSuperClusters", edm::InputTag("hybridSuperClusters")));
76  dynamicHybridSuperClusters_token = consumes<reco::SuperClusterCollection>(ps.getUntrackedParameter<edm::InputTag>(
77  "dynamicHybridSuperClusters", edm::InputTag("dynamicHybridSuperClusters")));
78  correctedHybridSuperClusters_token = consumes<reco::SuperClusterCollection>(ps.getUntrackedParameter<edm::InputTag>(
79  "correctedHybridSuperClusters", edm::InputTag("correctedHybridSuperClusters")));
81  consumes<reco::SuperClusterCollection>(ps.getUntrackedParameter<edm::InputTag>(
82  "correctedDynamicHybridSuperClusters", edm::InputTag("correctedDynamicHybridSuperClusters")));
83  correctedFixedMatrixSuperClustersWithPreshower_token = consumes<reco::SuperClusterCollection>(
84  ps.getUntrackedParameter<edm::InputTag>("correctedFixedMatrixSuperClustersWithPreshower",
85  edm::InputTag("correctedFixedMatrixSuperClustersWithPreshower")));
87  consumes<reco::SuperClusterCollection>(ps.getUntrackedParameter<edm::InputTag>(
88  "fixedMatrixSuperClustersWithPreshower", edm::InputTag("fixedMatrixSuperClustersWithPreshower")));
89 
90  outputFile_ = ps.getParameter<std::string>("outputFile");
91  rootFile_ = TFile::Open(outputFile_.c_str(),
92  "RECREATE"); // open output file to store histograms
93 
94  evtN = 0;
95 }
96 
97 //========================================================================
99 //========================================================================
100 {
101  delete rootFile_;
102 }
103 
104 //========================================================================
106  //========================================================================
107 
108  mytree_ = new TTree("energyScale", "");
109  TString treeVariables =
110  "mc_npar/I:parID:mc_sep/"
111  "F:mc_e:mc_et:mc_phi:mc_eta:mc_theta:"; // MC information
112  treeVariables += "em_dR/F:"; // MC <-> EM matching information
113  treeVariables +=
114  "em_isInCrack/I:em_scType:em_e/F:em_et:em_phi:em_eta:em_theta:em_nCell/"
115  "I:em_nBC:"; // EM SC info
116  treeVariables += "em_pet/F:em_pe:em_peta:em_ptheta:"; // EM SC physics (eta corrected
117  // information)
118 
119  treeVariables += "emCorr_e/F:emCorr_et:emCorr_eta:emCorr_phi:emCorr_theta:"; // CMSSW
120  // standard
121  // corrections
122  treeVariables += "emCorr_pet/F:emCorr_peta:emCorr_ptheta:"; // CMSSW standard physics
123 
124  treeVariables += "em_pw/F:em_ew:em_br"; // EM widths pw -- phiWidth, ew --
125  // etaWidth, ratios of pw/ew
126 
127  mytree_->Branch("energyScale", &(tree_.mc_npar), treeVariables);
128 }
129 
130 //========================================================================
132  using namespace edm; // needed for all fwk related classes
133  using namespace std;
134 
135  // std::cout << "Proceccing event # " << ++evtN << std::endl;
136 
137  // Get containers for MC truth, SC etc.
138  // ===================================================
139  // =======================================================================================
140  // =======================================================================================
143 
144  Labels l;
146 
147  const HepMC::GenEvent *genEvent = hepMC->GetEvent();
148  if (!(hepMC.isValid())) {
149  LogInfo("EnergyScaleAnalyzer") << "Could not get MC Product!";
150  return;
151  }
152 
153  //=======================For Vertex correction
154  std::vector<Handle<HepMCProduct>> evtHandles;
155 
156  //evt.getManyByType(evtHandles);
157  throw cms::Exception("UnsupportedFunction") << "EnergyScaleAnalyzer::analyze: "
158  << "getManyByType has not been supported by the Framework since 2015. "
159  << "This module has been broken since then. Maybe it should be deleted. "
160  << "Another possibility is to upgrade to use GetterOfProducts instead.";
161 
162  for (unsigned int i = 0; i < evtHandles.size(); ++i) {
163  if (evtHandles[i].isValid()) {
164  const HepMC::GenEvent *evt = evtHandles[i]->GetEvent();
165 
166  // take only 1st vertex for now - it's been tested only of PGuns...
167  //
168  HepMC::GenEvent::vertex_const_iterator vtx = evt->vertices_begin();
169  if (evtHandles[i].provenance()->moduleLabel() == std::string(l.module)) {
170  // Corrdinates of Vertex w.r.o. the point (0,0,0)
171  xVtx_ = 0.1 * (*vtx)->position().x();
172  yVtx_ = 0.1 * (*vtx)->position().y();
173  zVtx_ = 0.1 * (*vtx)->position().z();
174  }
175  }
176  }
177  //==============================================================================
178  // Get handle to SC collections
179 
181  try {
183  } catch (cms::Exception &ex) {
184  edm::LogError("EnergyScaleAnalyzer") << "Can't get collection with producer hybridSuperClusters.";
185  }
186 
188  try {
190  } catch (cms::Exception &ex) {
191  edm::LogError("EnergyScaleAnalyzer") << "Can't get collection with producer dynamicHybridSuperClusters.";
192  }
193 
194  Handle<reco::SuperClusterCollection> fixedMatrixSuperClustersWithPS;
195  try {
196  evt.getByToken(fixedMatrixSuperClustersWithPreshower_token, fixedMatrixSuperClustersWithPS);
197  } catch (cms::Exception &ex) {
198  edm::LogError("EnergyScaleAnalyzer") << "Can't get collection with producer "
199  "fixedMatrixSuperClustersWithPreshower.";
200  }
201 
202  // Corrected collections
203  Handle<reco::SuperClusterCollection> correctedHybridSC;
204  try {
205  evt.getByToken(correctedHybridSuperClusters_token, correctedHybridSC);
206  } catch (cms::Exception &ex) {
207  edm::LogError("EnergyScaleAnalyzer") << "Can't get collection with producer correctedHybridSuperClusters.";
208  }
209 
210  Handle<reco::SuperClusterCollection> correctedDynamicHybridSC;
211  try {
212  evt.getByToken(correctedDynamicHybridSuperClusters_token, correctedDynamicHybridSC);
213  } catch (cms::Exception &ex) {
214  edm::LogError("EnergyScaleAnalyzer") << "Can't get collection with producer "
215  "correctedDynamicHybridSuperClusters.";
216  }
217 
218  Handle<reco::SuperClusterCollection> correctedFixedMatrixSCWithPS;
219  try {
220  evt.getByToken(correctedFixedMatrixSuperClustersWithPreshower_token, correctedFixedMatrixSCWithPS);
221  } catch (cms::Exception &ex) {
222  edm::LogError("EnergyScaleAnalyzer") << "Can't get collection with producer "
223  "correctedFixedMatrixSuperClustersWithPreshower.";
224  }
225 
227  const reco::SuperClusterCollection *hSC = hybridSuperClusters.product();
228  const reco::SuperClusterCollection *fmSC = fixedMatrixSuperClustersWithPS.product();
229  const reco::SuperClusterCollection *chSC = correctedHybridSC.product();
230  const reco::SuperClusterCollection *cdSC = correctedDynamicHybridSC.product();
231  const reco::SuperClusterCollection *cfmSC = correctedFixedMatrixSCWithPS.product();
232 
233  // ----------------------- Print outs for debugging
234  /*
235  std::cout << "MC truth" << std::endl;
236  int counterI = 0;
237  for( HepMC::GenEvent::particle_const_iterator p =
238  genEvent->particles_begin(); p != genEvent->particles_end(); ++p ) { if (
239  fabs((*p)->momentum().eta()) < 1.5 ) { std::cout << ++counterI << " " <<
240  (*p)->momentum().e() << " "
241  << (*p)->momentum().phi() << " " << (*p)->momentum().eta() <<
242  std::endl;
243  }
244  }
245 
246  std::cout << "Standard clusters" << std::endl;
247  counterI = 0;
248  for(reco::SuperClusterCollection::const_iterator em = hSC->begin();
249  em != hSC->end(); ++em )
250  std::cout << ++counterI << " " << em->energy() << " " <<
251  em->position().phi() << " " << em->position().eta() << std::endl;
252 
253  std::cout << "Dynamic clusters" << std::endl;
254  counterI = 0;
255  for(reco::SuperClusterCollection::const_iterator em = dSC->begin();
256  em != dSC->end(); ++em )
257  std::cout << ++counterI << " " << em->energy() << " " <<
258  em->position().phi() << " " << em->position().eta() << std::endl;
259 
260  std::cout << "FixedMatrix clusters with PS" << std::endl;
261  counterI = 0;
262  for(reco::SuperClusterCollection::const_iterator em = fmSC->begin();
263  em != fmSC->end(); ++em )
264  std::cout << ++counterI << " " << em->energy() << " " <<
265  em->position().phi() << " " << em->position().eta() << std::endl;
266  */
267  // -----------------------------
268  //=====================================================================
269  // All containers are loaded, perform the analysis
270  //====================================================================
271 
272  // --------------------------- Store MC particles
273  HepMC::GenEvent::particle_const_iterator p = genEvent->particles_begin();
274 
275  // Search for MC electrons or photons that satisfy the criteria
276  float min_eT = 5;
277  float max_eta = 2.5;
278 
279  std::vector<HepMC::GenParticle *> mcParticles;
280  // int counter = 0;
281  for (HepMC::GenEvent::particle_const_iterator p = genEvent->particles_begin(); p != genEvent->particles_end(); ++p) {
282  // LogInfo("EnergyScaleAnalyzer") << "Particle " << ++counter
283  //<< " PDG ID = " << (*p)->pdg_id() << " pT = " << (*p)->momentum().perp();
284  // require photon or electron
285  if ((*p)->pdg_id() != 22 && abs((*p)->pdg_id()) != 11)
286  continue;
287 
288  // require selection criteria
289  bool satisfySelectionCriteria = (*p)->momentum().perp() > min_eT && fabs((*p)->momentum().eta()) < max_eta;
290 
291  if (!satisfySelectionCriteria)
292  continue;
293 
294  // EM MC particle is found, save it in the vector
295  mcParticles.push_back(*p);
296  }
297  // separation in dR between 2 first MC particles
298  // should not be used for MC samples with > 2 em objects generated!
299  if (mcParticles.size() == 2) {
300  HepMC::GenParticle *mc1 = mcParticles[0];
301  HepMC::GenParticle *mc2 = mcParticles[1];
302  tree_.mc_sep =
303  kinem::delta_R(mc1->momentum().eta(), mc1->momentum().phi(), mc2->momentum().eta(), mc2->momentum().phi());
304  } else
305  tree_.mc_sep = -100;
306 
307  // now loop over MC particles, find the match with SC and do everything we
308  // need then save info in the tree for every MC particle
309  for (std::vector<HepMC::GenParticle *>::const_iterator p = mcParticles.begin(); p != mcParticles.end(); ++p) {
310  HepMC::GenParticle *mc = *p;
311 
312  // Fill MC information
313  tree_.mc_npar = mcParticles.size();
314  tree_.parID = mc->pdg_id();
315  tree_.mc_e = mc->momentum().e();
316  tree_.mc_et = mc->momentum().e() * sin(mc->momentum().theta());
317  tree_.mc_phi = mc->momentum().phi();
318  tree_.mc_eta = mc->momentum().eta();
319  tree_.mc_theta = mc->momentum().theta();
320 
321  // Call function to fill tree
322  // scType coprreponds:
323  // HybridSuperCluster -- 1
324  // DynamicHybridSuperCluster -- 2
325  // FixedMatrixSuperClustersWithPreshower -- 3
326 
327  fillTree(hSC, chSC, mc, tree_, xVtx_, yVtx_, zVtx_, 1);
328  // std::cout << " TYPE " << 1 << " : " << tree_.em_e << " : " <<
329  // tree_.em_phi << " : " << tree_.em_eta << std::endl;
330 
331  fillTree(dSC, cdSC, mc, tree_, xVtx_, yVtx_, zVtx_, 2);
332  // std::cout << " TYPE " << 2 << " : " << tree_.em_e << " : " <<
333  // tree_.em_phi << " : " << tree_.em_eta << std::endl;
334 
335  fillTree(fmSC, cfmSC, mc, tree_, xVtx_, yVtx_, zVtx_, 3);
336  // std::cout << " TYPE " << 3 << " : " << tree_.em_e << " : " <<
337  // tree_.em_phi << " : " << tree_.em_eta << std::endl;
338 
339  // mytree_->Fill();
340  } // loop over particles
341 }
342 
344  const reco::SuperClusterCollection *corrSCColl,
347  float xV,
348  float yV,
349  float zV,
350  int scType) {
351  // ----------------------------- SuperClusters before energy correction
352  reco::SuperClusterCollection::const_iterator em = scColl->end();
353  float energyMax = -100.0; // dummy energy of the matched SC
354  for (reco::SuperClusterCollection::const_iterator aClus = scColl->begin(); aClus != scColl->end(); ++aClus) {
355  // check the matching
356  float dR =
357  kinem::delta_R(mc->momentum().eta(), mc->momentum().phi(), aClus->position().eta(), aClus->position().phi());
358  if (dR < 0.4) { // a rather loose matching cut
359  float energy = aClus->energy();
360  if (energy < energyMax)
361  continue;
362  energyMax = energy;
363  em = aClus;
364  }
365  }
366 
367  if (em == scColl->end()) {
368  // std::cout << "No matching SC with type " << scType << " was found for
369  // MC particle! " << std::endl; std::cout << "Going to next type of SC.
370  // " << std::endl;
371  return;
372  }
373  // ------------
374 
375  tree_.em_scType = scType;
376 
377  tree_.em_isInCrack = 0;
378  double emAbsEta = fabs(em->position().eta());
379  // copied from
380  // RecoEgama/EgammaElectronAlgos/src/EgammaElectronClassification.cc
381  if (emAbsEta < 0.018 || (emAbsEta > 0.423 && emAbsEta < 0.461) || (emAbsEta > 0.770 && emAbsEta < 0.806) ||
382  (emAbsEta > 1.127 && emAbsEta < 1.163) || (emAbsEta > 1.460 && emAbsEta < 1.558))
383  tree_.em_isInCrack = 1;
384 
385  tree_.em_dR = kinem::delta_R(mc->momentum().eta(), mc->momentum().phi(), em->position().eta(), em->position().phi());
386  tree_.em_e = em->energy();
387  tree_.em_et = em->energy() * sin(em->position().theta());
388  tree_.em_phi = em->position().phi();
389  tree_.em_eta = em->position().eta();
390  tree_.em_theta = em->position().theta();
391  tree_.em_nCell = em->size();
392  tree_.em_nBC = em->clustersSize();
393 
394  // Get physics e, et etc:
395  // Coordinates of EM object with respect of the point (0,0,0)
396  xClust_zero_ = em->position().x();
397  yClust_zero_ = em->position().y();
398  zClust_zero_ = em->position().z();
399  // Coordinates of EM object w.r.o. the Vertex position
400  xClust_vtx_ = xClust_zero_ - xV;
401  yClust_vtx_ = yClust_zero_ - yV;
402  zClust_vtx_ = zClust_zero_ - zV;
403 
404  energyMax_ = em->energy();
405  thetaMax_ = em->position().theta();
406  etaMax_ = em->position().eta();
407  phiMax_ = em->position().phi();
409  if (phiMax_ < 0)
410  phiMax_ += 2 * 3.14159;
411 
414  etaMaxVtx_ = -log(tan(thetaMaxVtx_ / 2));
417  if (phiMaxVtx_ < 0)
418  phiMaxVtx_ += 2 * 3.14159;
419  //=============================
420  // parametres of EM object after vertex correction
425 
426  //------------------------------- Get SC after energy correction
427  em = corrSCColl->end();
428  energyMax = -100.0; // dummy energy of the matched SC
429  for (reco::SuperClusterCollection::const_iterator aClus = corrSCColl->begin(); aClus != corrSCColl->end(); ++aClus) {
430  // check the matching
431  float dR =
432  kinem::delta_R(mc->momentum().eta(), mc->momentum().phi(), aClus->position().eta(), aClus->position().phi());
433  if (dR < 0.4) {
434  float energy = aClus->energy();
435  if (energy < energyMax)
436  continue;
437  energyMax = energy;
438  em = aClus;
439  }
440  }
441 
442  if (em == corrSCColl->end()) {
443  // std::cout << "No matching corrected SC with type " << scType << " was
444  // found for MC particle! " << std::endl; std::cout << "Going to next
445  // type of SC. " << std::endl;
446  return;
447  }
448  // ------------
449 
451  tree_.emCorr_e = em->energy();
452  tree_.emCorr_et = em->energy() * sin(em->position().theta());
453  tree_.emCorr_phi = em->position().phi();
454  tree_.emCorr_eta = em->position().eta();
455  tree_.emCorr_theta = em->position().theta();
456 
457  // =========== Eta and Theta wrt Vertex does not change after energy
458  // corrections are applied
459  // =========== So, no need to calculate them again
460 
464 
465  tree_.em_pw = em->phiWidth();
466  tree_.em_ew = em->etaWidth();
468 
469  mytree_->Fill();
470 }
471 
472 //==========================================================================
474  //========================================================================
475  // Fill ROOT tree
476  rootFile_->Write();
477 }
478 
edm::EDGetTokenT< reco::SuperClusterCollection > dynamicHybridSuperClusters_token
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
edm::EDGetTokenT< reco::SuperClusterCollection > fixedMatrixSuperClustersWithPreshower_token
const bool isValid(const Frame &aFrame, const FrameQuality &aQuality, const uint16_t aExpectedPos)
void analyze(const edm::Event &, const edm::EventSetup &) override
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T const * product() const
Definition: Handle.h:70
edm::EDGetTokenT< edm::HepMCProduct > hepMCLabel_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:528
edm::EDGetTokenT< reco::SuperClusterCollection > correctedFixedMatrixSuperClustersWithPreshower_token
Log< level::Error, false > LogError
edm::EDGetTokenT< reco::SuperClusterCollection > hybridSuperClusters_token
T getUntrackedParameter(std::string const &, T const &) const
void fillTree(const reco::SuperClusterCollection *scColl, const reco::SuperClusterCollection *corrSCColl, HepMC::GenParticle *mc, tree_structure_ &tree_, float xV, float yV, float zV, int scType)
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
T sqrt(T t)
Definition: SSEVec.h:19
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
double delta_R(double eta1, double phi1, double eta2, double phi2)
Definition: AnglesUtil.h:106
edm::EDGetTokenT< reco::SuperClusterCollection > correctedDynamicHybridSuperClusters_token
edm::EDGetTokenT< reco::SuperClusterCollection > correctedHybridSuperClusters_token
Log< level::Info, false > LogInfo
EnergyScaleAnalyzer(const edm::ParameterSet &)
HLT enums.
void labelsForToken(EDGetToken iToken, Labels &oLabels) const