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  evt.getManyByType(evtHandles);
156 
157  for (unsigned int i = 0; i < evtHandles.size(); ++i) {
158  if (evtHandles[i].isValid()) {
159  const HepMC::GenEvent *evt = evtHandles[i]->GetEvent();
160 
161  // take only 1st vertex for now - it's been tested only of PGuns...
162  //
163  HepMC::GenEvent::vertex_const_iterator vtx = evt->vertices_begin();
164  if (evtHandles[i].provenance()->moduleLabel() == std::string(l.module)) {
165  // Corrdinates of Vertex w.r.o. the point (0,0,0)
166  xVtx_ = 0.1 * (*vtx)->position().x();
167  yVtx_ = 0.1 * (*vtx)->position().y();
168  zVtx_ = 0.1 * (*vtx)->position().z();
169  }
170  }
171  }
172  //==============================================================================
173  // Get handle to SC collections
174 
176  try {
178  } catch (cms::Exception &ex) {
179  edm::LogError("EnergyScaleAnalyzer") << "Can't get collection with producer hybridSuperClusters.";
180  }
181 
183  try {
185  } catch (cms::Exception &ex) {
186  edm::LogError("EnergyScaleAnalyzer") << "Can't get collection with producer dynamicHybridSuperClusters.";
187  }
188 
189  Handle<reco::SuperClusterCollection> fixedMatrixSuperClustersWithPS;
190  try {
191  evt.getByToken(fixedMatrixSuperClustersWithPreshower_token, fixedMatrixSuperClustersWithPS);
192  } catch (cms::Exception &ex) {
193  edm::LogError("EnergyScaleAnalyzer") << "Can't get collection with producer "
194  "fixedMatrixSuperClustersWithPreshower.";
195  }
196 
197  // Corrected collections
198  Handle<reco::SuperClusterCollection> correctedHybridSC;
199  try {
200  evt.getByToken(correctedHybridSuperClusters_token, correctedHybridSC);
201  } catch (cms::Exception &ex) {
202  edm::LogError("EnergyScaleAnalyzer") << "Can't get collection with producer correctedHybridSuperClusters.";
203  }
204 
205  Handle<reco::SuperClusterCollection> correctedDynamicHybridSC;
206  try {
207  evt.getByToken(correctedDynamicHybridSuperClusters_token, correctedDynamicHybridSC);
208  } catch (cms::Exception &ex) {
209  edm::LogError("EnergyScaleAnalyzer") << "Can't get collection with producer "
210  "correctedDynamicHybridSuperClusters.";
211  }
212 
213  Handle<reco::SuperClusterCollection> correctedFixedMatrixSCWithPS;
214  try {
215  evt.getByToken(correctedFixedMatrixSuperClustersWithPreshower_token, correctedFixedMatrixSCWithPS);
216  } catch (cms::Exception &ex) {
217  edm::LogError("EnergyScaleAnalyzer") << "Can't get collection with producer "
218  "correctedFixedMatrixSuperClustersWithPreshower.";
219  }
220 
222  const reco::SuperClusterCollection *hSC = hybridSuperClusters.product();
223  const reco::SuperClusterCollection *fmSC = fixedMatrixSuperClustersWithPS.product();
224  const reco::SuperClusterCollection *chSC = correctedHybridSC.product();
225  const reco::SuperClusterCollection *cdSC = correctedDynamicHybridSC.product();
226  const reco::SuperClusterCollection *cfmSC = correctedFixedMatrixSCWithPS.product();
227 
228  // ----------------------- Print outs for debugging
229  /*
230  std::cout << "MC truth" << std::endl;
231  int counterI = 0;
232  for( HepMC::GenEvent::particle_const_iterator p =
233  genEvent->particles_begin(); p != genEvent->particles_end(); ++p ) { if (
234  fabs((*p)->momentum().eta()) < 1.5 ) { std::cout << ++counterI << " " <<
235  (*p)->momentum().e() << " "
236  << (*p)->momentum().phi() << " " << (*p)->momentum().eta() <<
237  std::endl;
238  }
239  }
240 
241  std::cout << "Standard clusters" << std::endl;
242  counterI = 0;
243  for(reco::SuperClusterCollection::const_iterator em = hSC->begin();
244  em != hSC->end(); ++em )
245  std::cout << ++counterI << " " << em->energy() << " " <<
246  em->position().phi() << " " << em->position().eta() << std::endl;
247 
248  std::cout << "Dynamic clusters" << std::endl;
249  counterI = 0;
250  for(reco::SuperClusterCollection::const_iterator em = dSC->begin();
251  em != dSC->end(); ++em )
252  std::cout << ++counterI << " " << em->energy() << " " <<
253  em->position().phi() << " " << em->position().eta() << std::endl;
254 
255  std::cout << "FixedMatrix clusters with PS" << std::endl;
256  counterI = 0;
257  for(reco::SuperClusterCollection::const_iterator em = fmSC->begin();
258  em != fmSC->end(); ++em )
259  std::cout << ++counterI << " " << em->energy() << " " <<
260  em->position().phi() << " " << em->position().eta() << std::endl;
261  */
262  // -----------------------------
263  //=====================================================================
264  // All containers are loaded, perform the analysis
265  //====================================================================
266 
267  // --------------------------- Store MC particles
268  HepMC::GenEvent::particle_const_iterator p = genEvent->particles_begin();
269 
270  // Search for MC electrons or photons that satisfy the criteria
271  float min_eT = 5;
272  float max_eta = 2.5;
273 
274  std::vector<HepMC::GenParticle *> mcParticles;
275  // int counter = 0;
276  for (HepMC::GenEvent::particle_const_iterator p = genEvent->particles_begin(); p != genEvent->particles_end(); ++p) {
277  // LogInfo("EnergyScaleAnalyzer") << "Particle " << ++counter
278  //<< " PDG ID = " << (*p)->pdg_id() << " pT = " << (*p)->momentum().perp();
279  // require photon or electron
280  if ((*p)->pdg_id() != 22 && abs((*p)->pdg_id()) != 11)
281  continue;
282 
283  // require selection criteria
284  bool satisfySelectionCriteria = (*p)->momentum().perp() > min_eT && fabs((*p)->momentum().eta()) < max_eta;
285 
286  if (!satisfySelectionCriteria)
287  continue;
288 
289  // EM MC particle is found, save it in the vector
290  mcParticles.push_back(*p);
291  }
292  // separation in dR between 2 first MC particles
293  // should not be used for MC samples with > 2 em objects generated!
294  if (mcParticles.size() == 2) {
295  HepMC::GenParticle *mc1 = mcParticles[0];
296  HepMC::GenParticle *mc2 = mcParticles[1];
297  tree_.mc_sep =
298  kinem::delta_R(mc1->momentum().eta(), mc1->momentum().phi(), mc2->momentum().eta(), mc2->momentum().phi());
299  } else
300  tree_.mc_sep = -100;
301 
302  // now loop over MC particles, find the match with SC and do everything we
303  // need then save info in the tree for every MC particle
304  for (std::vector<HepMC::GenParticle *>::const_iterator p = mcParticles.begin(); p != mcParticles.end(); ++p) {
305  HepMC::GenParticle *mc = *p;
306 
307  // Fill MC information
308  tree_.mc_npar = mcParticles.size();
309  tree_.parID = mc->pdg_id();
310  tree_.mc_e = mc->momentum().e();
311  tree_.mc_et = mc->momentum().e() * sin(mc->momentum().theta());
312  tree_.mc_phi = mc->momentum().phi();
313  tree_.mc_eta = mc->momentum().eta();
314  tree_.mc_theta = mc->momentum().theta();
315 
316  // Call function to fill tree
317  // scType coprreponds:
318  // HybridSuperCluster -- 1
319  // DynamicHybridSuperCluster -- 2
320  // FixedMatrixSuperClustersWithPreshower -- 3
321 
322  fillTree(hSC, chSC, mc, tree_, xVtx_, yVtx_, zVtx_, 1);
323  // std::cout << " TYPE " << 1 << " : " << tree_.em_e << " : " <<
324  // tree_.em_phi << " : " << tree_.em_eta << std::endl;
325 
326  fillTree(dSC, cdSC, mc, tree_, xVtx_, yVtx_, zVtx_, 2);
327  // std::cout << " TYPE " << 2 << " : " << tree_.em_e << " : " <<
328  // tree_.em_phi << " : " << tree_.em_eta << std::endl;
329 
330  fillTree(fmSC, cfmSC, mc, tree_, xVtx_, yVtx_, zVtx_, 3);
331  // std::cout << " TYPE " << 3 << " : " << tree_.em_e << " : " <<
332  // tree_.em_phi << " : " << tree_.em_eta << std::endl;
333 
334  // mytree_->Fill();
335  } // loop over particles
336 }
337 
339  const reco::SuperClusterCollection *corrSCColl,
342  float xV,
343  float yV,
344  float zV,
345  int scType) {
346  // ----------------------------- SuperClusters before energy correction
347  reco::SuperClusterCollection::const_iterator em = scColl->end();
348  float energyMax = -100.0; // dummy energy of the matched SC
349  for (reco::SuperClusterCollection::const_iterator aClus = scColl->begin(); aClus != scColl->end(); ++aClus) {
350  // check the matching
351  float dR =
352  kinem::delta_R(mc->momentum().eta(), mc->momentum().phi(), aClus->position().eta(), aClus->position().phi());
353  if (dR < 0.4) { // a rather loose matching cut
354  float energy = aClus->energy();
355  if (energy < energyMax)
356  continue;
357  energyMax = energy;
358  em = aClus;
359  }
360  }
361 
362  if (em == scColl->end()) {
363  // std::cout << "No matching SC with type " << scType << " was found for
364  // MC particle! " << std::endl; std::cout << "Going to next type of SC.
365  // " << std::endl;
366  return;
367  }
368  // ------------
369 
370  tree_.em_scType = scType;
371 
372  tree_.em_isInCrack = 0;
373  double emAbsEta = fabs(em->position().eta());
374  // copied from
375  // RecoEgama/EgammaElectronAlgos/src/EgammaElectronClassification.cc
376  if (emAbsEta < 0.018 || (emAbsEta > 0.423 && emAbsEta < 0.461) || (emAbsEta > 0.770 && emAbsEta < 0.806) ||
377  (emAbsEta > 1.127 && emAbsEta < 1.163) || (emAbsEta > 1.460 && emAbsEta < 1.558))
378  tree_.em_isInCrack = 1;
379 
380  tree_.em_dR = kinem::delta_R(mc->momentum().eta(), mc->momentum().phi(), em->position().eta(), em->position().phi());
381  tree_.em_e = em->energy();
382  tree_.em_et = em->energy() * sin(em->position().theta());
383  tree_.em_phi = em->position().phi();
384  tree_.em_eta = em->position().eta();
385  tree_.em_theta = em->position().theta();
386  tree_.em_nCell = em->size();
387  tree_.em_nBC = em->clustersSize();
388 
389  // Get physics e, et etc:
390  // Coordinates of EM object with respect of the point (0,0,0)
391  xClust_zero_ = em->position().x();
392  yClust_zero_ = em->position().y();
393  zClust_zero_ = em->position().z();
394  // Coordinates of EM object w.r.o. the Vertex position
395  xClust_vtx_ = xClust_zero_ - xV;
396  yClust_vtx_ = yClust_zero_ - yV;
397  zClust_vtx_ = zClust_zero_ - zV;
398 
399  energyMax_ = em->energy();
400  thetaMax_ = em->position().theta();
401  etaMax_ = em->position().eta();
402  phiMax_ = em->position().phi();
404  if (phiMax_ < 0)
405  phiMax_ += 2 * 3.14159;
406 
409  etaMaxVtx_ = -log(tan(thetaMaxVtx_ / 2));
412  if (phiMaxVtx_ < 0)
413  phiMaxVtx_ += 2 * 3.14159;
414  //=============================
415  // parametres of EM object after vertex correction
420 
421  //------------------------------- Get SC after energy correction
422  em = corrSCColl->end();
423  energyMax = -100.0; // dummy energy of the matched SC
424  for (reco::SuperClusterCollection::const_iterator aClus = corrSCColl->begin(); aClus != corrSCColl->end(); ++aClus) {
425  // check the matching
426  float dR =
427  kinem::delta_R(mc->momentum().eta(), mc->momentum().phi(), aClus->position().eta(), aClus->position().phi());
428  if (dR < 0.4) {
429  float energy = aClus->energy();
430  if (energy < energyMax)
431  continue;
432  energyMax = energy;
433  em = aClus;
434  }
435  }
436 
437  if (em == corrSCColl->end()) {
438  // std::cout << "No matching corrected SC with type " << scType << " was
439  // found for MC particle! " << std::endl; std::cout << "Going to next
440  // type of SC. " << std::endl;
441  return;
442  }
443  // ------------
444 
446  tree_.emCorr_e = em->energy();
447  tree_.emCorr_et = em->energy() * sin(em->position().theta());
448  tree_.emCorr_phi = em->position().phi();
449  tree_.emCorr_eta = em->position().eta();
450  tree_.emCorr_theta = em->position().theta();
451 
452  // =========== Eta and Theta wrt Vertex does not change after energy
453  // corrections are applied
454  // =========== So, no need to calculate them again
455 
459 
460  tree_.em_pw = em->phiWidth();
461  tree_.em_ew = em->etaWidth();
463 
464  mytree_->Fill();
465 }
466 
467 //==========================================================================
469  //========================================================================
470  // Fill ROOT tree
471  rootFile_->Write();
472 }
473 
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)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
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:539
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
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 &)
void getManyByType(std::vector< Handle< PROD >> &results) const
Definition: Event.h:530
HLT enums.
void labelsForToken(EDGetToken iToken, Labels &oLabels) const