CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EgammaSuperClusters.cc
Go to the documentation of this file.
2 
6 
11 
13 
16 
18 {
19  outputFile_ = ps.getUntrackedParameter<std::string>("outputFile", "");
20  //CMSSW_Version_ = ps.getUntrackedParameter<std::string>("CMSSW_Version", "");
21 
22  verboseDBE_ = ps.getUntrackedParameter<bool>("verboseDBE", false);
23 
24  hist_min_Size_ = ps.getParameter<double>("hist_min_Size");
25  hist_max_Size_ = ps.getParameter<double>("hist_max_Size");
26  hist_bins_Size_ = ps.getParameter<int> ("hist_bins_Size");
27 
28  hist_min_NumBC_ = ps.getParameter<double>("hist_min_NumBC");
29  hist_max_NumBC_ = ps.getParameter<double>("hist_max_NumBC");
30  hist_bins_NumBC_ = ps.getParameter<int> ("hist_bins_NumBC");
31 
32  hist_min_ET_ = ps.getParameter<double>("hist_min_ET");
33  hist_max_ET_ = ps.getParameter<double>("hist_max_ET");
34  hist_bins_ET_ = ps.getParameter<int> ("hist_bins_ET");
35 
36  hist_min_Eta_ = ps.getParameter<double>("hist_min_Eta");
37  hist_max_Eta_ = ps.getParameter<double>("hist_max_Eta");
38  hist_bins_Eta_ = ps.getParameter<int> ("hist_bins_Eta");
39 
40  hist_min_Phi_ = ps.getParameter<double>("hist_min_Phi");
41  hist_max_Phi_ = ps.getParameter<double>("hist_max_Phi");
42  hist_bins_Phi_ = ps.getParameter<int> ("hist_bins_Phi");
43 
44  hist_min_S1toS9_ = ps.getParameter<double>("hist_min_S1toS9");
45  hist_max_S1toS9_ = ps.getParameter<double>("hist_max_S1toS9");
46  hist_bins_S1toS9_ = ps.getParameter<int> ("hist_bins_S1toS9");
47 
48  hist_min_S25toE_ = ps.getParameter<double>("hist_min_S25toE");
49  hist_max_S25toE_ = ps.getParameter<double>("hist_max_S25toE");
50  hist_bins_S25toE_ = ps.getParameter<int> ("hist_bins_S25toE");
51 
52  hist_min_EoverTruth_ = ps.getParameter<double>("hist_min_EoverTruth");
53  hist_max_EoverTruth_ = ps.getParameter<double>("hist_max_EoverTruth");
54  hist_bins_EoverTruth_ = ps.getParameter<int> ("hist_bins_EoverTruth");
55 
56  hist_min_deltaR_ = ps.getParameter<double>("hist_min_deltaR");
57  hist_max_deltaR_ = ps.getParameter<double>("hist_max_deltaR");
58  hist_bins_deltaR_ = ps.getParameter<int> ("hist_bins_deltaR");
59 
60  hist_min_phiWidth_ = ps.getParameter<double>("hist_min_phiWidth");
61  hist_max_phiWidth_ = ps.getParameter<double>("hist_max_phiWidth");
62  hist_bins_phiWidth_ = ps.getParameter<int>("hist_bins_phiWidth");
63 
64  hist_min_etaWidth_ = ps.getParameter<double>("hist_min_etaWidth");
65  hist_max_etaWidth_ = ps.getParameter<double>("hist_max_etaWidth");
66  hist_bins_etaWidth_ = ps.getParameter<int>("hist_bins_etaWidth");
67 
68  hist_bins_preshowerE_ = ps.getParameter<int>("hist_bins_preshowerE");
69  hist_min_preshowerE_ = ps.getParameter<double>("hist_min_preshowerE");
70  hist_max_preshowerE_ = ps.getParameter<double>("hist_max_preshowerE");
71 
72  hist_min_R_ = ps.getParameter<double>("hist_min_R");
73  hist_max_R_ = ps.getParameter<double>("hist_max_R");
74  hist_bins_R_ = ps.getParameter<int> ("hist_bins_R");
75 
76  MCTruthCollection_ = ps.getParameter<edm::InputTag>("MCTruthCollection");
77 
78  barrelRawSuperClusterCollection_ = ps.getParameter<edm::InputTag>("barrelRawSuperClusterCollection");
79  barrelCorSuperClusterCollection_ = ps.getParameter<edm::InputTag>("barrelCorSuperClusterCollection");
80 
81  endcapRawSuperClusterCollection_ = ps.getParameter<edm::InputTag>("endcapRawSuperClusterCollection");
82  endcapPreSuperClusterCollection_ = ps.getParameter<edm::InputTag>("endcapPreSuperClusterCollection");
83  endcapCorSuperClusterCollection_ = ps.getParameter<edm::InputTag>("endcapCorSuperClusterCollection");
84 
85  barrelRecHitCollection_ = ps.getParameter<edm::InputTag>("barrelRecHitCollection");
86  endcapRecHitCollection_ = ps.getParameter<edm::InputTag>("endcapRecHitCollection");
87 }
88 
90 
92 {
94 
95  if ( verboseDBE_ )
96  {
97  dbe_->setVerbose(1);
99  }
100  else
101  dbe_->setVerbose(0);
102 
103  //dbe_->setCurrentFolder("Ecal/CMSSW_"+CMSSW_Version_+"/EcalClusters/SuperClusters/");
104  dbe_->setCurrentFolder("EcalClusterV/EcalSuperClusters/");
105 
106  // Number of SuperClusters
107  //
109  = dbe_->book1D("hist_EB_RawSC_Size_","# Raw SuperClusters in Barrel",
112  = dbe_->book1D("hist_EE_RawSC_Size_","# Raw SuperClusters in Endcap",
115  = dbe_->book1D("hist_EB_CorSC_Size_","# Corrected SuperClusters in Barrel",
118  = dbe_->book1D("hist_EE_CorSC_Size_","# Corrected SuperClusters in Endcap",
121  = dbe_->book1D("hist_EE_PreSC_Size_","# SuperClusters with Preshower in Endcap",
123 
124  // Number of BasicClusters in SuperCluster
125  //
127  = dbe_->book1D("hist_EB_RawSC_NumBC_","# of Basic Clusters in Raw Super Clusters in Barrel",
130  = dbe_->book1D("hist_EE_RawSC_NumBC_","# of Basic Clusters in Raw Super Clusters in Endcap",
133  = dbe_->book1D("hist_EB_CorSC_NumBC_","# of Basic Clusters in Corrected SuperClusters in Barrel",
136  = dbe_->book1D("hist_EE_CorSC_NumBC_","# of Basic Clusters in Corrected SuperClusters in Endcap",
139  = dbe_->book1D("hist_EE_PreSC_NumBC_","# of Basic Clusters in SuperClusters with Preshower in Endcap",
141 
142  // ET distribution of SuperClusters
143  //
145  = dbe_->book1D("hist_EB_RawSC_ET_","ET of Raw SuperClusters in Barrel",
148  = dbe_->book1D("hist_EE_RawSC_ET_","ET of Raw SuperClusters in Endcap",
151  = dbe_->book1D("hist_EB_CorSC_ET_","ET of Corrected SuperClusters in Barrel",
154  = dbe_->book1D("hist_EE_CorSC_ET_","ET of Corrected SuperClusters in Endcap",
157  = dbe_->book1D("hist_EE_PreSC_ET_","ET of SuperClusters with Preshower in Endcap",
159 
160  // Eta distribution of SuperClusters
161  //
163  = dbe_->book1D("hist_EB_RawSC_Eta_","Eta of Raw SuperClusters in Barrel",
166  = dbe_->book1D("hist_EE_RawSC_Eta_","Eta of Raw SuperClusters in Endcap",
169  = dbe_->book1D("hist_EB_CorSC_Eta_","Eta of Corrected SuperClusters in Barrel",
172  = dbe_->book1D("hist_EE_CorSC_Eta_","Eta of Corrected SuperClusters in Endcap",
175  = dbe_->book1D("hist_EE_PreSC_Eta_","Eta of SuperClusters with Preshower in Endcap",
177 
178  // Phi distribution of SuperClusters
179  //
181  = dbe_->book1D("hist_EB_RawSC_Phi_","Phi of Raw SuperClusters in Barrel",
184  = dbe_->book1D("hist_EE_RawSC_Phi_","Phi of Raw SuperClusters in Endcap",
187  = dbe_->book1D("hist_EB_CorSC_Phi_","Phi of Corrected SuperClusters in Barrel",
190  = dbe_->book1D("hist_EE_CorSC_Phi_","Phi of Corrected SuperClusters in Endcap",
193  = dbe_->book1D("hist_EE_PreSC_Phi_","Phi of SuperClusters with Preshower in Endcap",
195 
196  // S1/S9 distribution of SuperClusters
197  //
199  = dbe_->book1D("hist_EB_RawSC_S1toS9_","S1/S9 of Raw Super Clusters in Barrel",
202  = dbe_->book1D("hist_EE_RawSC_S1toS9_","S1/S9 of Raw Super Clusters in Endcap",
205  = dbe_->book1D("hist_EB_CorSC_S1toS9_","S1/S9 of Corrected SuperClusters in Barrel",
208  = dbe_->book1D("hist_EE_CorSC_S1toS9_","S1/S9 of Corrected SuperClusters in Endcap",
211  = dbe_->book1D("hist_EE_PreSC_S1toS9_","S1/S9 of SuperClusters with Preshower in Endcap",
213 
214  // S25/E distribution of SuperClusters
215  //
217  = dbe_->book1D("hist_EB_RawSC_S25toE_","S25/E of Raw Super Clusters in Barrel",
220  = dbe_->book1D("hist_EE_RawSC_S25toE_","S25/E of Raw Super Clusters in Endcap",
223  = dbe_->book1D("hist_EB_CorSC_S25toE_","S25/E of Corrected SuperClusters in Barrel",
226  = dbe_->book1D("hist_EE_CorSC_S25toE_","S25/E of Corrected SuperClusters in Endcap",
229  = dbe_->book1D("hist_EE_PreSC_S25toE_","S25/E of SuperClusters with Preshower in Endcap",
231 
232  // E/E(true) distribution of SuperClusters
233  //
235  = dbe_->book1D("hist_EB_RawSC_EoverTruth_","E/True E of Raw SuperClusters in Barrel",
238  = dbe_->book1D("hist_EE_RawSC_EoverTruth_","E/True E of Raw SuperClusters in Endcap",
241  = dbe_->book1D("hist_EB_CorSC_EoverTruth_","E/True E of Corrected SuperClusters in Barrel",
244  = dbe_->book1D("hist_EE_CorSC_EoverTruth_","E/True E of Corrected SuperClusters in Endcap",
247  = dbe_->book1D("hist_EE_PreSC_EoverTruth_","E/True E of SuperClusters with Preshower in Endcap",
249 
250  // dR distribution of SuperClusters from truth
251  //
253  = dbe_->book1D("hist_EB_RawSC_deltaR_","dR to MC truth of Raw Super Clusters in Barrel",
256  = dbe_->book1D("hist_EE_RawSC_deltaR_","dR to MC truth of Raw Super Clusters in Endcap",
259  = dbe_->book1D("hist_EB_CorSC_deltaR_","dR to MC truth of Corrected SuperClusters in Barrel",
262  = dbe_->book1D("hist_EE_CorSC_deltaR_","dR to MC truth of Corrected SuperClusters in Endcap",
265  = dbe_->book1D("hist_EE_PreSC_deltaR_","dR to MC truth of SuperClusters with Preshower in Endcap",
267 
268  // phi width stored in corrected SuperClusters
270  = dbe_->book1D("hist_EB_CorSC_phiWidth_","phiWidth of Corrected Super Clusters in Barrel",
273  = dbe_->book1D("hist_EE_CorSC_phiWidth_","phiWidth of Corrected Super Clusters in Endcap",
275 
276  // eta width stored in corrected SuperClusters
278  = dbe_->book1D("hist_EB_CorSC_etaWidth_","etaWidth of Corrected Super Clusters in Barrel",
281  = dbe_->book1D("hist_EE_CorSC_etaWidth_","etaWidth of Corrected Super Clusters in Endcap",
283 
284 
285  // preshower energy
287  = dbe_->book1D("hist_EE_PreSC_preshowerE_","preshower energy in Super Clusters with Preshower in Endcap",
290  = dbe_->book1D("hist_EE_CorSC_preshowerE_","preshower energy in Corrected Super Clusters with Preshower in Endcap",
292 
293 
294  //
295  hist_EB_CorSC_ET_vs_Eta_ = dbe_->book2D( "hist_EB_CorSC_ET_vs_Eta_", "Corr Super Cluster ET versus Eta in Barrel",
298 
299  hist_EB_CorSC_ET_vs_Phi_ = dbe_->book2D( "hist_EB_CorSC_ET_vs_Phi_", "Corr Super Cluster ET versus Phi in Barrel",
302 
303  hist_EE_CorSC_ET_vs_Eta_ = dbe_->book2D( "hist_EE_CorSC_ET_vs_Eta_", "Corr Super Cluster ET versus Eta in Endcap",
306 
307  hist_EE_CorSC_ET_vs_Phi_ = dbe_->book2D( "hist_EE_CorSC_ET_vs_Phi_", "Corr Super Cluster ET versus Phi in Endcap",
310 
311  hist_EE_CorSC_ET_vs_R_ = dbe_->book2D( "hist_EE_CorSC_ET_vs_R_", "Corr Super Cluster ET versus Radius in Endcap",
314 
315 
316 }
317 
319 {
320 
321  bool skipMC = false;
322  bool skipBarrel = false;
323  bool skipEndcap = false;
324 
325  //
326  // Get MCTRUTH
327  //
329  evt.getByLabel(MCTruthCollection_, pMCTruth);
330  if (!pMCTruth.isValid()) {
331  edm::LogError("EgammaSuperClusters") << "Error! can't get collection with label "
333  skipMC = true;
334  }
335  const HepMC::GenEvent* genEvent = pMCTruth->GetEvent();
336 
337  if( skipMC ) return;
338 
339  //
340  // Get the BARREL products
341  //
342  edm::Handle<reco::SuperClusterCollection> pBarrelRawSuperClusters;
343  evt.getByLabel(barrelRawSuperClusterCollection_, pBarrelRawSuperClusters);
344  if (!pBarrelRawSuperClusters.isValid()) {
345  edm::LogError("EgammaSuperClusters") << "Error! can't get collection with label "
347  skipBarrel = true;
348  }
349 
350  edm::Handle<reco::SuperClusterCollection> pBarrelCorSuperClusters;
351  evt.getByLabel(barrelCorSuperClusterCollection_, pBarrelCorSuperClusters);
352  if (!pBarrelCorSuperClusters.isValid()) {
353  edm::LogError("EgammaSuperClusters") << "Error! can't get collection with label "
355  skipBarrel = true;
356  }
357 
358  edm::Handle< EBRecHitCollection > pBarrelRecHitCollection;
359  evt.getByLabel( barrelRecHitCollection_, pBarrelRecHitCollection );
360  if ( ! pBarrelRecHitCollection.isValid() ) {
361  skipBarrel = true;
362  }
363  edm::Handle< EERecHitCollection > pEndcapRecHitCollection;
364  evt.getByLabel( endcapRecHitCollection_, pEndcapRecHitCollection );
365  if ( ! pEndcapRecHitCollection.isValid() ) {
366  skipEndcap = true;
367  }
368 
369  if( skipBarrel || skipEndcap ) return;
370 
372 
373  // Get the BARREL collections
374  const reco::SuperClusterCollection* barrelRawSuperClusters = pBarrelRawSuperClusters.product();
375  const reco::SuperClusterCollection* barrelCorSuperClusters = pBarrelCorSuperClusters.product();
376 
377  // Number of entries in collections
378  hist_EB_RawSC_Size_->Fill(barrelRawSuperClusters->size());
379  hist_EB_CorSC_Size_->Fill(barrelCorSuperClusters->size());
380 
381  // Do RAW BARREL SuperClusters
382  for(reco::SuperClusterCollection::const_iterator aClus = barrelRawSuperClusters->begin();
383  aClus != barrelRawSuperClusters->end(); aClus++)
384  {
385  // kinematics
386  hist_EB_RawSC_NumBC_->Fill(aClus->clustersSize());
387  hist_EB_RawSC_ET_->Fill(aClus->energy()/std::cosh(aClus->position().eta()));
388  hist_EB_RawSC_Eta_->Fill(aClus->position().eta());
389  hist_EB_RawSC_Phi_->Fill(aClus->position().phi());
390 
391  // cluster shape
392  const reco::CaloClusterPtr seed = aClus->seed();
393  hist_EB_RawSC_S1toS9_->Fill( lazyTool.eMax( *seed ) / lazyTool.e3x3( *seed ) );
394  hist_EB_RawSC_S25toE_->Fill( lazyTool.e5x5( *seed ) / aClus->energy() );
395 
396  // truth
397  double dRClosest = 999.9;
398  double energyClosest = 0;
399  closestMCParticle(genEvent, *aClus, dRClosest, energyClosest);
400 
401  if (dRClosest < 0.1)
402  {
403  hist_EB_RawSC_EoverTruth_->Fill(aClus->energy()/energyClosest);
404  hist_EB_RawSC_deltaR_->Fill(dRClosest);
405 
406  }
407 
408  }
409 
410  // Do CORRECTED BARREL SuperClusters
411  for(reco::SuperClusterCollection::const_iterator aClus = barrelCorSuperClusters->begin();
412  aClus != barrelCorSuperClusters->end(); aClus++)
413  {
414  // kinematics
415  hist_EB_CorSC_NumBC_->Fill(aClus->clustersSize());
416  hist_EB_CorSC_ET_->Fill(aClus->energy()/std::cosh(aClus->position().eta()));
417  hist_EB_CorSC_Eta_->Fill(aClus->position().eta());
418  hist_EB_CorSC_Phi_->Fill(aClus->position().phi());
419 
420  hist_EB_CorSC_ET_vs_Eta_->Fill( aClus->energy()/std::cosh(aClus->position().eta()), aClus->eta() );
421  hist_EB_CorSC_ET_vs_Phi_->Fill( aClus->energy()/std::cosh(aClus->position().eta()), aClus->phi() );
422 
423 
424  // cluster shape
425  const reco::CaloClusterPtr seed = aClus->seed();
426  hist_EB_CorSC_S1toS9_->Fill( lazyTool.eMax( *seed ) / lazyTool.e3x3( *seed ) );
427  hist_EB_CorSC_S25toE_->Fill( lazyTool.e5x5( *seed ) / aClus->energy() );
428 
429  // correction variables
430  hist_EB_CorSC_phiWidth_->Fill(aClus->phiWidth());
431  hist_EB_CorSC_etaWidth_->Fill(aClus->etaWidth());
432 
433  // truth
434  double dRClosest = 999.9;
435  double energyClosest = 0;
436  closestMCParticle(genEvent, *aClus, dRClosest, energyClosest);
437 
438  if (dRClosest < 0.1)
439  {
440  hist_EB_CorSC_EoverTruth_->Fill(aClus->energy()/energyClosest);
441  hist_EB_CorSC_deltaR_->Fill(dRClosest);
442 
443  }
444 
445  }
446 
447  //
448  // Get the ENDCAP products
449  //
450  edm::Handle<reco::SuperClusterCollection> pEndcapRawSuperClusters;
451  evt.getByLabel(endcapRawSuperClusterCollection_, pEndcapRawSuperClusters);
452  if (!pEndcapRawSuperClusters.isValid()) {
453  edm::LogError("EgammaSuperClusters") << "Error! can't get collection with label "
455  }
456 
457  edm::Handle<reco::SuperClusterCollection> pEndcapPreSuperClusters;
458  evt.getByLabel(endcapPreSuperClusterCollection_, pEndcapPreSuperClusters);
459  if (!pEndcapPreSuperClusters.isValid()) {
460  edm::LogError("EgammaSuperClusters") << "Error! can't get collection with label "
462  }
463 
464  edm::Handle<reco::SuperClusterCollection> pEndcapCorSuperClusters;
465  evt.getByLabel(endcapCorSuperClusterCollection_, pEndcapCorSuperClusters);
466  if (!pEndcapCorSuperClusters.isValid()) {
467  edm::LogError("EgammaSuperClusters") << "Error! can't get collection with label "
469  }
470 
471  // Get the ENDCAP collections
472  const reco::SuperClusterCollection* endcapRawSuperClusters = pEndcapRawSuperClusters.product();
473  const reco::SuperClusterCollection* endcapPreSuperClusters = pEndcapPreSuperClusters.product();
474  const reco::SuperClusterCollection* endcapCorSuperClusters = pEndcapCorSuperClusters.product();
475 
476  // Number of entries in collections
477  hist_EE_RawSC_Size_->Fill(endcapRawSuperClusters->size());
478  hist_EE_PreSC_Size_->Fill(endcapPreSuperClusters->size());
479  hist_EE_CorSC_Size_->Fill(endcapCorSuperClusters->size());
480 
481  // Do RAW ENDCAP SuperClusters
482  for(reco::SuperClusterCollection::const_iterator aClus = endcapRawSuperClusters->begin();
483  aClus != endcapRawSuperClusters->end(); aClus++)
484  {
485  hist_EE_RawSC_NumBC_->Fill(aClus->clustersSize());
486  hist_EE_RawSC_ET_->Fill(aClus->energy()/std::cosh(aClus->position().eta()));
487  hist_EE_RawSC_Eta_->Fill(aClus->position().eta());
488  hist_EE_RawSC_Phi_->Fill(aClus->position().phi());
489 
490  const reco::CaloClusterPtr seed = aClus->seed();
491  hist_EE_RawSC_S1toS9_->Fill( lazyTool.eMax( *seed ) / lazyTool.e3x3( *seed ) );
492  hist_EE_RawSC_S25toE_->Fill( lazyTool.e5x5( *seed ) / aClus->energy() );
493 
494  // truth
495  double dRClosest = 999.9;
496  double energyClosest = 0;
497  closestMCParticle(genEvent, *aClus, dRClosest, energyClosest);
498 
499  if (dRClosest < 0.1)
500  {
501  hist_EE_RawSC_EoverTruth_->Fill(aClus->energy()/energyClosest);
502  hist_EE_RawSC_deltaR_->Fill(dRClosest);
503  }
504 
505  }
506 
507  // Do ENDCAP SuperClusters with PRESHOWER
508  for(reco::SuperClusterCollection::const_iterator aClus = endcapPreSuperClusters->begin();
509  aClus != endcapPreSuperClusters->end(); aClus++)
510  {
511  hist_EE_PreSC_NumBC_->Fill(aClus->clustersSize());
512  hist_EE_PreSC_ET_->Fill(aClus->energy()/std::cosh(aClus->position().eta()));
513  hist_EE_PreSC_Eta_->Fill(aClus->position().eta());
514  hist_EE_PreSC_Phi_->Fill(aClus->position().phi());
515  hist_EE_PreSC_preshowerE_->Fill(aClus->preshowerEnergy());
516 
517  const reco::CaloClusterPtr seed = aClus->seed();
518  hist_EE_PreSC_S1toS9_->Fill( lazyTool.eMax( *seed ) / lazyTool.e3x3( *seed ) );
519  hist_EE_PreSC_S25toE_->Fill( lazyTool.e5x5( *seed ) / aClus->energy() );
520 
521  // truth
522  double dRClosest = 999.9;
523  double energyClosest = 0;
524  closestMCParticle(genEvent, *aClus, dRClosest, energyClosest);
525 
526  if (dRClosest < 0.1)
527  {
528  hist_EE_PreSC_EoverTruth_->Fill(aClus->energy()/energyClosest);
529  hist_EE_PreSC_deltaR_->Fill(dRClosest);
530  }
531 
532  }
533 
534  // Do CORRECTED ENDCAP SuperClusters
535  for(reco::SuperClusterCollection::const_iterator aClus = endcapCorSuperClusters->begin();
536  aClus != endcapCorSuperClusters->end(); aClus++)
537  {
538  hist_EE_CorSC_NumBC_->Fill(aClus->clustersSize());
539  hist_EE_CorSC_ET_->Fill(aClus->energy()/std::cosh(aClus->position().eta()));
540  hist_EE_CorSC_Eta_->Fill(aClus->position().eta());
541  hist_EE_CorSC_Phi_->Fill(aClus->position().phi());
542  hist_EE_CorSC_preshowerE_->Fill(aClus->preshowerEnergy());
543 
544  hist_EE_CorSC_ET_vs_Eta_->Fill( aClus->energy()/std::cosh(aClus->position().eta()), aClus->eta() );
545  hist_EE_CorSC_ET_vs_Phi_->Fill( aClus->energy()/std::cosh(aClus->position().eta()), aClus->phi() );
546  hist_EE_CorSC_ET_vs_R_->Fill( aClus->energy()/std::cosh(aClus->position().eta()),
547  std::sqrt( std::pow(aClus->x(),2) + std::pow(aClus->y(),2) ) );
548 
549 
550  // correction variables
551  hist_EE_CorSC_phiWidth_->Fill(aClus->phiWidth());
552  hist_EE_CorSC_etaWidth_->Fill(aClus->etaWidth());
553 
554  const reco::CaloClusterPtr seed = aClus->seed();
555  hist_EE_CorSC_S1toS9_->Fill( lazyTool.eMax( *seed ) / lazyTool.e3x3( *seed ) );
556  hist_EE_CorSC_S25toE_->Fill( lazyTool.e5x5( *seed ) / aClus->energy() );
557 
558  // truth
559  double dRClosest = 999.9;
560  double energyClosest = 0;
561  closestMCParticle(genEvent, *aClus, dRClosest, energyClosest);
562 
563  if (dRClosest < 0.1)
564  {
565  hist_EE_CorSC_EoverTruth_->Fill(aClus->energy()/energyClosest);
566  hist_EE_CorSC_deltaR_->Fill(dRClosest);
567  }
568 
569  }
570 
571 }
572 
574 {
575  if (outputFile_.size() != 0 && dbe_) dbe_->save(outputFile_);
576 }
577 
578 //
579 // Closest MC Particle
580 //
582  double &dRClosest, double &energyClosest)
583 {
584 
585  // SuperCluster eta, phi
586  double scEta = sc.eta();
587  double scPhi = sc.phi();
588 
589  // initialize dRClosest to a large number
590  dRClosest = 999.9;
591 
592  // loop over the MC truth particles to find the
593  // closest to the superCluster in dR space
594  for(HepMC::GenEvent::particle_const_iterator currentParticle = genEvent->particles_begin();
595  currentParticle != genEvent->particles_end(); currentParticle++ )
596  {
597  if((*currentParticle)->status() == 1)
598  {
599  // need GenParticle in ECAL co-ordinates
600  HepMC::FourVector vtx = (*currentParticle)->production_vertex()->position();
601  double phiTrue = (*currentParticle)->momentum().phi();
602  double etaTrue = ecalEta((*currentParticle)->momentum().eta(), vtx.z()/10., vtx.perp()/10.);
603 
604  double dPhi = reco::deltaPhi(phiTrue, scPhi);
605  double dEta = scEta - etaTrue;
606  double deltaR = std::sqrt(dPhi*dPhi + dEta*dEta);
607 
608  if(deltaR < dRClosest)
609  {
610  dRClosest = deltaR;
611  energyClosest = (*currentParticle)->momentum().e();
612  }
613 
614  } // end if stable particle
615 
616  } // end loop on get particles
617 
618 }
619 
620 
621 //
622 // Compute Eta in the ECAL co-ordinate system
623 //
624 float EgammaSuperClusters::ecalEta(float EtaParticle , float Zvertex, float plane_Radius)
625 {
626  const float R_ECAL = 136.5;
627  const float Z_Endcap = 328.0;
628  const float etaBarrelEndcap = 1.479;
629 
630  if(EtaParticle != 0.)
631  {
632  float Theta = 0.0 ;
633  float ZEcal = (R_ECAL-plane_Radius)*sinh(EtaParticle)+Zvertex;
634 
635  if(ZEcal != 0.0) Theta = atan(R_ECAL/ZEcal);
636  if(Theta<0.0) Theta = Theta+Geom::pi() ;
637 
638  float ETA = - log(tan(0.5*Theta));
639 
640  if( fabs(ETA) > etaBarrelEndcap )
641  {
642  float Zend = Z_Endcap ;
643  if(EtaParticle<0.0 ) Zend = -Zend ;
644  float Zlen = Zend - Zvertex ;
645  float RR = Zlen/sinh(EtaParticle);
646  Theta = atan((RR+plane_Radius)/Zend);
647  if(Theta<0.0) Theta = Theta+Geom::pi() ;
648  ETA = - log(tan(0.5*Theta));
649  }
650 
651  return ETA;
652  }
653  else
654  {
655  edm::LogWarning("") << "[EgammaSuperClusters::ecalEta] Warning: Eta equals to zero, not correcting" ;
656  return EtaParticle;
657  }
658 }
659 
660 
T getParameter(std::string const &) const
MonitorElement * hist_EE_CorSC_ET_vs_Phi_
T getUntrackedParameter(std::string const &, T const &) const
MonitorElement * hist_EE_CorSC_ET_vs_R_
MonitorElement * hist_EE_PreSC_ET_
edm::InputTag barrelRawSuperClusterCollection_
MonitorElement * hist_EB_CorSC_ET_
edm::InputTag endcapCorSuperClusterCollection_
edm::InputTag MCTruthCollection_
MonitorElement * hist_EB_CorSC_NumBC_
MonitorElement * hist_EE_CorSC_etaWidth_
virtual void analyze(const edm::Event &, const edm::EventSetup &)
MonitorElement * hist_EE_PreSC_EoverTruth_
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:872
MonitorElement * hist_EE_CorSC_Eta_
MonitorElement * hist_EB_RawSC_Size_
MonitorElement * hist_EE_RawSC_NumBC_
MonitorElement * hist_EB_RawSC_S1toS9_
MonitorElement * hist_EB_RawSC_ET_
MonitorElement * hist_EE_CorSC_EoverTruth_
MonitorElement * hist_EE_CorSC_deltaR_
MonitorElement * hist_EB_RawSC_Phi_
MonitorElement * hist_EE_RawSC_ET_
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:161
void Fill(long long x)
MonitorElement * hist_EB_CorSC_phiWidth_
MonitorElement * hist_EE_PreSC_S1toS9_
MonitorElement * hist_EB_RawSC_Eta_
#define ETA
MonitorElement * hist_EB_CorSC_deltaR_
MonitorElement * hist_EE_RawSC_S1toS9_
edm::InputTag endcapRecHitCollection_
MonitorElement * hist_EB_CorSC_S25toE_
double dPhi(double phi1, double phi2)
Definition: JetUtil.h:30
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
MonitorElement * hist_EE_CorSC_S1toS9_
MonitorElement * hist_EE_RawSC_deltaR_
T sqrt(T t)
Definition: SSEVec.h:48
MonitorElement * hist_EE_PreSC_Phi_
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
MonitorElement * hist_EB_CorSC_ET_vs_Eta_
void save(const std::string &filename, const std::string &path="", const std::string &pattern="", const std::string &rewrite="", const uint32_t run=0, SaveReferenceTag ref=SaveWithReference, int minStatus=dqm::qstatus::STATUS_OK, const std::string &fileupdate="RECREATE")
Definition: DQMStore.cc:2296
MonitorElement * hist_EB_CorSC_EoverTruth_
tuple genEvent
Definition: MCTruth.py:33
MonitorElement * hist_EB_CorSC_etaWidth_
void setVerbose(unsigned level)
Definition: DQMStore.cc:548
edm::InputTag endcapRawSuperClusterCollection_
MonitorElement * hist_EE_CorSC_ET_
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:390
MonitorElement * hist_EE_PreSC_deltaR_
MonitorElement * hist_EE_CorSC_Size_
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:12
MonitorElement * hist_EE_PreSC_preshowerE_
MonitorElement * hist_EB_CorSC_Phi_
edm::InputTag barrelCorSuperClusterCollection_
float ecalEta(float EtaParticle, float Zvertex, float plane_Radius)
MonitorElement * hist_EB_CorSC_S1toS9_
MonitorElement * hist_EB_CorSC_Size_
MonitorElement * hist_EE_CorSC_Phi_
void closestMCParticle(const HepMC::GenEvent *genEvent, const reco::SuperCluster &sc, double &dRClosest, double &energyClosest)
MonitorElement * hist_EB_RawSC_EoverTruth_
static const float etaBarrelEndcap
T const * product() const
Definition: Handle.h:81
MonitorElement * hist_EE_RawSC_EoverTruth_
edm::InputTag barrelRecHitCollection_
EgammaSuperClusters(const edm::ParameterSet &)
static const float Z_Endcap
std::string const & label() const
Definition: InputTag.h:42
MonitorElement * hist_EB_RawSC_S25toE_
MonitorElement * hist_EE_PreSC_Eta_
MonitorElement * hist_EB_CorSC_Eta_
MonitorElement * hist_EB_RawSC_NumBC_
MonitorElement * hist_EE_CorSC_preshowerE_
double pi()
Definition: Pi.h:31
static const float R_ECAL
MonitorElement * hist_EE_RawSC_Size_
MonitorElement * hist_EE_CorSC_NumBC_
MonitorElement * hist_EE_CorSC_S25toE_
void showDirStructure(void) const
Definition: DQMStore.cc:2961
double phi() const
azimuthal angle of cluster centroid
Definition: CaloCluster.h:164
MonitorElement * hist_EE_RawSC_S25toE_
MonitorElement * book2D(const char *name, const char *title, int nchX, double lowX, double highX, int nchY, double lowY, double highY)
Book 2D histogram.
Definition: DQMStore.cc:1000
MonitorElement * hist_EE_CorSC_ET_vs_Eta_
MonitorElement * hist_EE_PreSC_NumBC_
edm::InputTag endcapPreSuperClusterCollection_
MonitorElement * hist_EE_PreSC_S25toE_
MonitorElement * hist_EE_RawSC_Eta_
MonitorElement * hist_EE_PreSC_Size_
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:584
MonitorElement * hist_EB_RawSC_deltaR_
MonitorElement * hist_EE_RawSC_Phi_
MonitorElement * hist_EB_CorSC_ET_vs_Phi_
MonitorElement * hist_EE_CorSC_phiWidth_