CMS 3D CMS Logo

ElectronCalibration.cc
Go to the documentation of this file.
1 
3 
22 #include "TFile.h"
23 #include "TH1.h"
24 #include "TH2.h"
25 #include "TF1.h"
26 #include "TRandom.h"
27 
28 #include <iostream>
29 #include <string>
30 #include <stdexcept>
31 #include <vector>
32 
34  rootfile_ = iConfig.getParameter<std::string>("rootfile");
35  recHitToken_ = consumes<EBRecHitCollection>(iConfig.getParameter<edm::InputTag>("ebRecHitsLabel"));
36  electronToken_ = consumes<reco::GsfElectronCollection>(iConfig.getParameter<edm::InputTag>("electronLabel"));
37  trackLabel_ = iConfig.getParameter<edm::InputTag>("trackLabel");
38  calibAlgo_ = iConfig.getParameter<std::string>("CALIBRATION_ALGO");
39  edm::LogVerbatim("EcalCalib") << " The used Algorithm is " << calibAlgo_;
40  keventweight_ = iConfig.getParameter<int>("keventweight");
41  ClusterSize_ = iConfig.getParameter<int>("Clustersize");
42  ElePt_ = iConfig.getParameter<double>("ElePt");
43  maxeta_ = iConfig.getParameter<int>("maxeta");
44  mineta_ = iConfig.getParameter<int>("mineta");
45  maxphi_ = iConfig.getParameter<int>("maxphi");
46  minphi_ = iConfig.getParameter<int>("minphi");
47  cut1_ = iConfig.getParameter<double>("cut1");
48  cut2_ = iConfig.getParameter<double>("cut2");
49  cut3_ = iConfig.getParameter<double>("cut3");
50  elecclass_ = iConfig.getParameter<int>("elecclass");
51  edm::LogVerbatim("EcalCalib") << " The electronclass is " << elecclass_;
52  numevent_ = iConfig.getParameter<int>("numevent");
53  miscalibfile_ = iConfig.getParameter<std::string>("miscalibfile");
54 
55  cutEPCalo1_ = iConfig.getParameter<double>("cutEPCaloMin");
56  cutEPCalo2_ = iConfig.getParameter<double>("cutEPCaloMax");
57  cutEPin1_ = iConfig.getParameter<double>("cutEPinMin");
58  cutEPin2_ = iConfig.getParameter<double>("cutEPinMax");
59  cutCalo1_ = iConfig.getParameter<double>("cutCaloMin");
60  cutCalo2_ = iConfig.getParameter<double>("cutCaloMax");
61 
62  cutESeed_ = iConfig.getParameter<double>("cutESeed");
63 }
64 
66 
67 //========================================================================
69  //========================================================================
70  f = new TFile(rootfile_.c_str(), "RECREATE");
71 
72  // Book histograms
73  e9 = new TH1F("e9", "E9 energy", 300, 0., 150.);
74  e25 = new TH1F("e25", "E25 energy", 300, 0., 150.);
75  scE = new TH1F("scE", "SC energy", 300, 0., 150.);
76  trP = new TH1F("trP", "Trk momentum", 300, 0., 150.);
77  EoP = new TH1F("EoP", "EoP", 600, 0., 3.);
78  EoP_all = new TH1F("EoP_all", "EoP_all", 600, 0., 3.);
79 
80  if (elecclass_ == 0 || elecclass_ == -1) {
81  calibs = new TH1F("calib", "Calibration constants", 4000, 0.5, 2.);
82  } else {
83  calibs = new TH1F("calib", "Calibration constants", 800, 0.5, 2.);
84  }
85 
86  e25OverScE = new TH1F("e25OverscE", "E25 / SC energy", 400, 0., 2.);
87  E25oP = new TH1F("E25oP", "E25 / P", 1200, 0., 1.5);
88 
89  Map = new TH2F("Map", "Nb Events in Crystal", 85, 1, 85, 70, 5, 75);
90  e9Overe25 = new TH1F("e9Overe25", "E9 / E25", 400, 0., 2.);
91  Map3Dcalib = new TH2F("3Dcalib", "3Dcalib", 85, 1, 85, 70, 5, 75);
92 
93  MapCor1 = new TH2F("MapCor1", "Correlation E25/Pcalo versus E25/Pin", 100, 0., 5., 100, 0., 5.);
94  MapCor2 = new TH2F("MapCor2", "Correlation E25/Pcalo versus E/P", 100, 0., 5., 100, 0., 5.);
95  MapCor3 = new TH2F("MapCor3", "Correlation E25/Pcalo versus Pout/Pin", 100, 0., 5., 100, 0., 5.);
96  MapCor4 = new TH2F("MapCor4", "Correlation E25/Pcalo versus E25/highestP", 100, 0., 5., 100, 0., 5.);
97  MapCor5 = new TH2F("MapCor5", "Correlation E25/Pcalo versus Pcalo/Pout", 100, 0., 5., 100, 0., 5.);
98  MapCor6 = new TH2F("MapCor6", "Correlation Pout/Pin versus E25/Pin", 100, 0., 5., 100, 0., 5.);
99  MapCor7 = new TH2F("MapCor7", "Correlation Pout/Pin versus Pcalo/Pout", 100, 0., 5., 100, 0., 5.);
100  MapCor8 = new TH2F("MapCor8", "Correlation E25/Pin versus Pcalo/Pout", 100, 0., 5., 100, 0., 5.);
101  MapCor9 = new TH2F("MapCor9", "Correlation E25/Pcalo versus Eseed/Pout", 100, 0., 5., 100, 0., 5.);
102  MapCor10 = new TH2F("MapCor10", "Correlation Eseed/Pout versus Pout/Pin", 100, 0., 5., 100, 0., 5.);
103  MapCor11 = new TH2F("MapCor11", "Correlation Eseed/Pout versus E25/Pin", 100, 0., 5., 100, 0., 5.);
104  MapCorCalib =
105  new TH2F("MapCorCalib", "Correlation Miscalibration versus Calibration constants", 100, 0.5, 1.5, 100, 0.5, 1.5);
106 
107  PinMinPout = new TH1F("PinMinPout", "(Pin - Pout)/Pin", 600, -2.0, 2.0);
108 
109  if (elecclass_ == 0 || elecclass_ == -1) {
110  calibinter = new TH1F("calibinter", "internal calibration constants", 2000, 0.5, 2.);
111  PinOverPout = new TH1F("PinOverPout", "pinOverpout", 600, 0., 3.);
112  eSeedOverPout = new TH1F("eSeedOverPout", "eSeedOverpout ", 600, 0., 3.);
113  MisCalibs = new TH1F("MisCalibs", "Miscalibration constants", 4000, 0.5, 2.);
114  RatioCalibs = new TH1F("RatioCalibs", "Ratio in Calibration Constants", 4000, 0.5, 2.0);
115  DiffCalibs = new TH1F("DiffCalibs", "Difference in Calibration constants", 4000, -1.0, 1.0);
116  } else {
117  calibinter = new TH1F("calibinter", "internal calibration constants", 400, 0.5, 2.);
118  PinOverPout = new TH1F("PinOverPout", "pinOverpout", 600, 0., 3.);
119  eSeedOverPout = new TH1F("eSeedOverPout", "eSeedOverpout ", 600, 0., 3.);
120  MisCalibs = new TH1F("MisCalibs", "Miscalibration constants", 800, 0.5, 2.);
121  RatioCalibs = new TH1F("RatioCalibs", "Ratio in Calibration Constants", 800, 0.5, 2.0);
122  DiffCalibs = new TH1F("DiffCalibs", "Difference in Calibration constants", 800, -1.0, 1.0);
123  }
124  Error1 = new TH1F("Error1", "DeltaP/Pin", 800, -1.0, 1.0);
125  Error2 = new TH1F("Error2", "DeltaP/Pout", 800, -1.0, 1.0);
126  Error3 = new TH1F("Error3", "DeltaP/Pcalo", 800, -1.0, 1.0);
127  eSeedOverPout2 = new TH1F("eSeedOverPout2", "eSeedOverpout (No Supercluster)", 600, 0., 4.);
128  hadOverEm = new TH1F("hadOverEm", "Had/EM distribution", 600, -2., 2.);
129 
130  // Book histograms
131  Map3DcalibNoCuts = new TH2F("3DcalibNoCuts", "3Dcalib (Before Cuts)", 85, 1, 85, 70, 5, 75);
132  e9NoCuts = new TH1F("e9NoCuts", "E9 energy (Before Cuts)", 300, 0., 150.);
133  e25NoCuts = new TH1F("e25NoCuts", "E25 energy (Before Cuts)", 300, 0., 150.);
134  scENoCuts = new TH1F("scENoCuts", "SC energy (Before Cuts)", 300, 0., 150.);
135  trPNoCuts = new TH1F("trPNoCuts", "Trk momentum (Before Cuts)", 300, 0., 150.);
136  EoPNoCuts = new TH1F("EoPNoCuts", "EoP (Before Cuts)", 600, 0., 3.);
137  if (elecclass_ == 0 || elecclass_ == -1) {
138  calibsNoCuts = new TH1F("calibNoCuts", "Calibration constants (Before Cuts)", 4000, 0., 2.);
139  } else {
140  calibsNoCuts = new TH1F("calibNoCuts", "Calibration constants (Before Cuts)", 800, 0., 2.);
141  }
142  e25OverScENoCuts = new TH1F("e25OverscENoCuts", "E25 / SC energy (Before Cuts)", 400, 0., 2.);
143  E25oPNoCuts = new TH1F("E25oPNoCuts", "E25 / P (Before Cuts)", 1200, 0., 1.5);
144  MapNoCuts = new TH2F("MapNoCuts", "Nb Events in Crystal (Before Cuts)", 85, 1, 85, 70, 5, 75);
145  e9Overe25NoCuts = new TH1F("e9Overe25NoCuts", "E9 / E25 (Before Cuts)", 400, 0., 2.);
146  PinOverPoutNoCuts = new TH1F("PinOverPoutNoCuts", "pinOverpout (Before Cuts)", 600, 0., 3.);
147  eSeedOverPoutNoCuts = new TH1F(" eSeedOverPoutNoCuts", "eSeedOverpout (Before Cuts) ", 600, 0., 4.);
148  PinMinPoutNoCuts = new TH1F("PinMinPoutNoCuts", "(Pin - Pout)/Pin (Before Cuts)", 600, -2.0, 2.0);
149 
150  RatioCalibsNoCuts = new TH1F("RatioCalibsNoCuts", "Ratio in Calibration Constants (Before Cuts)", 4000, 0.5, 2.0);
151  DiffCalibsNoCuts = new TH1F("DiffCalibsNoCuts", "Difference in Calibration constants (Before Cuts)", 4000, -1.0, 1.0);
152  calibinterNoCuts = new TH1F("calibinterNoCuts", "internal calibration constants", 2000, 0.5, 2.);
153 
154  MapCor1NoCuts =
155  new TH2F("MapCor1NoCuts", "Correlation E25/PatCalo versus E25/Pin (Before Cuts)", 100, 0., 5., 100, 0., 5.);
156  MapCor2NoCuts =
157  new TH2F("MapCor2NoCuts", "Correlation E25/PatCalo versus E/P (Before Cuts)", 100, 0., 5., 100, 0., 5.);
158  MapCor3NoCuts =
159  new TH2F("MapCor3NoCuts", "Correlation E25/PatCalo versus Pout/Pin (Before Cuts)", 100, 0., 5., 100, 0., 5.);
160  MapCor4NoCuts =
161  new TH2F("MapCor4NoCuts", "Correlation E25/PatCalo versus E25/highestP (Before Cuts)", 100, 0., 5., 100, 0., 5.);
162  MapCor5NoCuts =
163  new TH2F("MapCor5NoCuts", "Correlation E25/Pcalo versus Pcalo/Pout (Before Cuts)", 100, 0., 5., 100, 0., 5.);
164  MapCor6NoCuts =
165  new TH2F("MapCor6NoCuts", "Correlation Pout/Pin versus E25/Pin (Before Cuts)", 100, 0., 5., 100, 0., 5.);
166  MapCor7NoCuts =
167  new TH2F("MapCor7NoCuts", "Correlation Pout/Pin versus Pcalo/Pout (Before Cuts)", 100, 0., 5., 100, 0., 5.);
168  MapCor8NoCuts =
169  new TH2F("MapCor8NoCuts", "Correlation E25/Pin versus Pcalo/Pout (Before Cuts)", 100, 0., 5., 100, 0., 5.);
170  MapCor9NoCuts =
171  new TH2F("MapCor9NoCuts", "Correlation E25/Pcalo versus Eseed/Pout (Before Cuts)", 100, 0., 5., 100, 0., 5.);
173  new TH2F("MapCor10NoCuts", "Correlation Eseed/Pout versus Pout/Pin (Before Cuts)", 100, 0., 5., 100, 0., 5.);
175  new TH2F("MapCor11NoCuts", "Correlation Eseed/Pout versus E25/Pin (Before Cuts)", 100, 0., 5., 100, 0., 5.);
176  MapCorCalibNoCuts = new TH2F("MapCorCalibNoCuts",
177  "Correlation Miscalibration versus Calibration constants (Before Cuts)",
178  100,
179  0.,
180  3.,
181  100,
182  0.,
183  3.);
184 
185  Error1NoCuts = new TH1F("Eror1NoCuts", "DeltaP/Pin (Before Cuts)", 800, -1.0, 1.0);
186  Error2NoCuts = new TH1F("Error2NoCuts", "DeltaP/Pout (Before Cuts)", 800, -1.0, 1.0);
187  Error3NoCuts = new TH1F("Error3NoCuts", "DeltaP/Pcalo (Before Cuts)", 800, -1.0, 1.0);
188  eSeedOverPout2NoCuts = new TH1F("eSeedOverPout2NoCuts", "eSeedOverpout (No Supercluster, Before Cuts)", 600, 0., 4.);
189  hadOverEmNoCuts = new TH1F("hadOverEmNoCuts", "Had/EM distribution (Before Cuts)", 600, -2., 2.);
190 
191  //Book histograms after ESeed cut
192  MapCor1ESeed =
193  new TH2F("MapCor1ESeed", "Correlation E25/Pcalo versus E25/Pin (after Eseed/Pout cut)", 100, 0., 5., 100, 0., 5.);
194  MapCor2ESeed =
195  new TH2F("MapCor2ESeed", "Correlation E25/Pcalo versus E/P (after Eseed/Pout cut)", 100, 0., 5., 100, 0., 5.);
196  MapCor3ESeed = new TH2F(
197  "MapCor3ESeed", "Correlation E25/Pcalo versus Pout/Pin (after Eseed/Pout cut)", 100, 0., 5., 100, 0., 5.);
198  MapCor4ESeed = new TH2F(
199  "MapCor4ESeed", "Correlation E25/Pcalo versus E25/highestP (after Eseed/Pout cut)", 100, 0., 5., 100, 0., 5.);
200  MapCor5ESeed = new TH2F(
201  "MapCor5ESeed", "Correlation E25/Pcalo versus Pcalo/Pout (after Eseed/Pout cut)", 100, 0., 5., 100, 0., 5.);
202  MapCor6ESeed =
203  new TH2F("MapCor6ESeed", "Correlation Pout/Pin versus E25/Pin (after Eseed/Pout cut)", 100, 0., 5., 100, 0., 5.);
204  MapCor7ESeed = new TH2F(
205  "MapCor7ESeed", "Correlation Pout/Pin versus Pcalo/Pout (after Eseed/Pout cut)", 100, 0., 5., 100, 0., 5.);
206  MapCor8ESeed = new TH2F(
207  "MapCor8ESeed", "Correlation E25/Pin versus Pcalo/Pout (after Eseed/Pout cut)", 100, 0., 5., 100, 0., 5.);
208  MapCor9ESeed = new TH2F(
209  "MapCor9ESeed", "Correlation E25/Pcalo versus Eseed/Pout (after Eseed/Pout cut)", 100, 0., 5., 100, 0., 5.);
210  MapCor10ESeed = new TH2F(
211  "MapCor10ESeed", "Correlation Eseed/Pout versus Pout/Pin (after Eseed/Pout cut)", 100, 0., 5., 100, 0., 5.);
212  MapCor11ESeed = new TH2F(
213  "MapCor11ESeed", "Correlation Eseed/Pout versus E25/Pin (after Eseed/Pout cut)", 100, 0., 5., 100, 0., 5.);
214 
216  new TH1F("eSeedOverPout2ESeed", "eSeedOverpout (No Supercluster, after Eseed/Pout cut)", 600, 0., 4.);
217 
218  hadOverEmESeed = new TH1F("hadOverEmESeed", "Had/EM distribution (after Eseed/Pout cut)", 600, -2., 2.);
219 
220  //Book histograms without any cut
221  GeneralMap = new TH2F("GeneralMap", "Map without any cuts", 85, 1, 85, 70, 5, 75);
222 
224  etaMin = mineta_;
225  etaMax = maxeta_;
226  phiMin = minphi_;
227  phiMax = maxphi_;
228  if (calibAlgo_ == "L3") {
230  } else {
231  if (calibAlgo_ == "HH" || calibAlgo_ == "HHReg") {
233  } else {
234  edm::LogError("EcalCalib") << " Name of Algorithm is not recognize " << calibAlgo_
235  << " Should be either L3, HH or HHReg. Abort! ";
236  }
237  }
238  read_events = 0;
239 
240  // get Region to be calibrated
242 
243  oldCalibs.resize(ReducedMap.size(), 0.);
244 
245  // table is set to zero
246  for (int phi = 0; phi < 360; phi++) {
247  for (int eta = 0; eta < 171; eta++) {
248  eventcrystal[eta][phi] = 0;
249  }
250  }
251 
252  edm::LogVerbatim("EcalCalib") << " Begin JOB ";
253 }
254 
255 //========================================================================
256 
258  //========================================================================
259 
260  int nIterations = 10;
261  if (calibAlgo_ == "L3") {
263  } else {
264  if (calibAlgo_ == "HH") {
266  } else {
267  if (calibAlgo_ == "HHReg") {
269  } else {
270  edm::LogError("EcalCalib") << " Calibration not run due to problem in Algo Choice...";
271  return;
272  }
273  }
274  }
275  for (int ii = 0; ii < (int)solution.size(); ii++) {
276  edm::LogVerbatim("EcalCalib") << "solution[" << ii << "] = " << solution[ii];
277  calibs->Fill(solution[ii]);
278  }
279 
280  newCalibs.resize(ReducedMap.size(), 0.);
281 
282  calibXMLwriter write_calibrations;
283 
284  FILE* MisCalib;
285  MisCalib = fopen(miscalibfile_.c_str(), "r");
286  int fileStatus = 1;
287  int eta = -1;
288  int phi = -1;
289  float coeff = -1;
290 
291  std::map<EBDetId, float> OldCoeff;
292 
293  while (fileStatus != EOF) {
294  fileStatus = fscanf(MisCalib, "%d %d %f\n", &eta, &phi, &coeff);
295  if (eta != -1 && phi != -1 && coeff != -1) {
296  OldCoeff.insert(std::make_pair(EBDetId(eta, phi, EBDetId::ETAPHIMODE), coeff));
297  }
298  }
299 
300  fclose(MisCalib);
301 
302  int icry = 0;
303  CalibrationCluster::CalibMap::iterator itmap;
304  for (itmap = ReducedMap.begin(); itmap != ReducedMap.end(); itmap++) {
305  newCalibs[icry] = solution[icry];
306 
307  write_calibrations.writeLine(itmap->first, newCalibs[icry]);
308  float Compare = 1.;
309  std::map<EBDetId, float>::iterator iter = OldCoeff.find(itmap->first);
310  if (iter != OldCoeff.end())
311  Compare = iter->second;
312 
313  if ((itmap->first).ieta() > mineta_ && (itmap->first).ieta() < maxeta_ && (itmap->first).iphi() > minphi_ &&
314  (itmap->first).iphi() < maxphi_) {
315  Map3Dcalib->Fill((itmap->first).ieta(), (itmap->first).iphi(), newCalibs[icry] * Compare);
316  MisCalibs->Fill(Compare);
317  }
318  if ((itmap->first).ieta() < mineta_ + 2) {
319  icry++;
320  continue;
321  }
322  if ((itmap->first).ieta() > maxeta_ - 2) {
323  icry++;
324  continue;
325  }
326  if ((itmap->first).iphi() < minphi_ + 2) {
327  icry++;
328  continue;
329  }
330  if ((itmap->first).iphi() > maxphi_ - 2) {
331  icry++;
332  continue;
333  }
334 
335  calibinter->Fill(newCalibs[icry]);
336  DiffCalibs->Fill(newCalibs[icry] - 1. / Compare);
337  RatioCalibs->Fill(newCalibs[icry] * Compare);
338  MapCorCalib->Fill(1. / Compare, newCalibs[icry]);
339  icry++;
340  }
341 
342  if (calibAlgo_ == "L3") {
345  } else {
346  if (calibAlgo_ == "HH") {
348  } else {
349  if (calibAlgo_ == "HHReg") {
351  } else {
352  edm::LogError("EcalCalib") << " Calibration not run due to problem in AlgoChoice...";
353  return;
354  }
355  }
356  }
357  for (int ii = 0; ii < (int)solutionNoCuts.size(); ii++) {
359  }
360  int icryp = 0;
361  CalibrationCluster::CalibMap::iterator itmapp;
362  for (itmapp = ReducedMap.begin(); itmapp != ReducedMap.end(); itmapp++) {
363  newCalibs[icryp] = solutionNoCuts[icryp];
364  float Compare2 = 1.;
365  std::map<EBDetId, float>::iterator iter2 = OldCoeff.find(itmapp->first);
366  if (iter2 != OldCoeff.end())
367  Compare2 = iter2->second;
368 
369  if ((itmapp->first).ieta() > mineta_ && (itmapp->first).ieta() < maxeta_ && (itmapp->first).iphi() > minphi_ &&
370  (itmapp->first).iphi() < maxphi_)
371  Map3DcalibNoCuts->Fill((itmapp->first).ieta(), (itmapp->first).iphi(), newCalibs[icryp] * Compare2);
372  if ((itmapp->first).ieta() < mineta_ + 2) {
373  icryp++;
374  continue;
375  }
376  if ((itmapp->first).ieta() > maxeta_ - 2) {
377  icryp++;
378  continue;
379  }
380  if ((itmapp->first).iphi() < minphi_ + 2) {
381  icryp++;
382  continue;
383  }
384  if ((itmapp->first).iphi() > maxphi_ - 2) {
385  icryp++;
386  continue;
387  }
388  calibinterNoCuts->Fill(newCalibs[icryp]);
389  DiffCalibsNoCuts->Fill(newCalibs[icryp] - 1. / (Compare2));
390  RatioCalibsNoCuts->Fill(newCalibs[icryp] * Compare2);
391  MapCorCalibNoCuts->Fill(1. / Compare2, newCalibs[icryp]);
392  icryp++;
393  }
394 
396 
397  edm::LogVerbatim("EcalCalib") << "\n************* STATISTICS **************";
398  edm::LogVerbatim("EcalCalib") << " Events Studied " << read_events;
399 
401 
402  f->Write();
403 
404  f->Close();
405 }
406 
407 //=================================================================================
409  //=================================================================================
410 
411  EcalRecHitCollection ecrh = *phits;
413  EBDetId save;
414  float en_save = 0;
415  for (it = ecrh.begin(); it != ecrh.end(); it++) {
416  EBDetId p = EBDetId(it->id().rawId());
417  if (it->energy() > en_save) {
418  en_save = it->energy();
419  save = p;
420  }
421  }
422  return save;
423 }
424 
425 //=================================================================================
426 EBDetId ElectronCalibration::findMaxHit2(const std::vector<DetId>& v1, const EBRecHitCollection* hits) {
427  //=================================================================================
428 
429  double currEnergy = 0.;
430  EBDetId maxHit;
431 
432  for (std::vector<DetId>::const_iterator idsIt = v1.begin(); idsIt != v1.end(); ++idsIt) {
433  if (idsIt->subdetId() != 1)
434  continue;
436  itrechit = hits->find(*idsIt);
437 
438  if (itrechit == hits->end()) {
439  edm::LogVerbatim("EcalCalib") << "ElectronCalibration::findMaxHit2: rechit not found! ";
440  continue;
441  }
442  if (itrechit->energy() > currEnergy) {
443  currEnergy = itrechit->energy();
444  maxHit = *idsIt;
445  }
446  }
447 
448  return maxHit;
449 }
450 
451 //=================================================================================
453  //=================================================================================
454  using namespace edm;
455 
456  // Get EBRecHits
458  if (!phits.isValid()) {
459  edm::LogError("EcalCalib") << "Error! can't get the product EBRecHitCollection: ";
460  }
461 
462  const EBRecHitCollection* hits = phits.product(); // get a ptr to the product
463 
464  // Get pixelElectrons
465  const Handle<reco::GsfElectronCollection>& pElectrons = iEvent.getHandle(electronToken_);
466  if (!pElectrons.isValid()) {
467  edm::LogError("EcalCalib") << "Error! can't get the product ElectronCollection: ";
468  }
469 
471  read_events++;
472  if (read_events % 1000 == 0)
473  edm::LogVerbatim("EcalCalib") << "read_events = " << read_events << std::endl;
474 
475  if (!hits)
476  return;
477  if (hits->empty())
478  return;
479  if (!electronCollection)
480  return;
481  if (electronCollection->empty())
482  return;
483 
485  // START HERE....
487  reco::GsfElectronCollection::const_iterator eleIt = electronCollection->begin();
488 
489  reco::GsfElectron highPtElectron;
490 
491  float highestElePt = 0.;
492  bool found = false;
493  for (eleIt = electronCollection->begin(); eleIt != electronCollection->end(); eleIt++) {
494  //Comments
495  if (fabs(eleIt->eta()) > (maxeta_ + 3) * 0.0175)
496  continue;
497  if (eleIt->eta() < (mineta_ - 3) * 0.0175)
498  continue;
499 
500  if (eleIt->pt() > highestElePt) {
501  highestElePt = eleIt->pt();
502  highPtElectron = *eleIt;
503  found = true;
504  }
505  }
506  if (highestElePt < ElePt_)
507  return;
508  if (!found)
509  return;
510  const reco::SuperCluster& sc = *(highPtElectron.superCluster());
511  if (fabs(sc.eta()) > (maxeta_ + 3) * 0.0175) {
512  edm::LogVerbatim("EcalCalib") << "++++ Problem with electron, electron eta is " << highPtElectron.eta()
513  << " while SC is " << sc.eta() << std::endl;
514  return;
515  }
516 
517  std::vector<DetId> v1;
518  //Loop to fill the vector of DetIds
519  for (std::vector<std::pair<DetId, float> >::const_iterator idsIt = sc.hitsAndFractions().begin();
520  idsIt != sc.hitsAndFractions().end();
521  ++idsIt) {
522  v1.push_back(idsIt->first);
523  }
524 
525  //Change function name
526  EBDetId maxHitId;
527 
528  maxHitId = findMaxHit2(v1, hits);
529 
530  if (maxHitId.null()) {
531  edm::LogVerbatim("EcalCalib") << " Null " << std::endl;
532  return;
533  }
534 
535  int maxCC_Eta = maxHitId.ieta();
536  int maxCC_Phi = maxHitId.iphi();
537 
538  if (maxCC_Eta > maxeta_)
539  return;
540  if (maxCC_Eta < mineta_)
541  return;
542  if (maxCC_Phi > maxphi_)
543  return;
544  if (maxCC_Phi < minphi_)
545  return;
546 
547  // number of events per crystal is set
548  if (numevent_ > 0) {
549  eventcrystal[maxCC_Eta + 85][maxCC_Phi - 1] += 1;
550  if (eventcrystal[maxCC_Eta + 85][maxCC_Phi - 1] > numevent_)
551  return;
552  }
553 
554  std::vector<EBDetId> Xtals5x5 = calibCluster.get5x5Id(maxHitId);
555 
556  if ((int)Xtals5x5.size() != ClusterSize_ * ClusterSize_)
557  return;
558 
559  // fill cluster energy
560  std::vector<float> energy;
561  float energy3x3 = 0.;
562  float energy5x5 = 0.;
563 
564  for (int icry = 0; icry < ClusterSize_ * ClusterSize_; icry++) {
566  if (Xtals5x5[icry].subdetId() != 1)
567  continue;
568  itrechit = hits->find(Xtals5x5[icry]);
569  if (itrechit == hits->end()) {
570  edm::LogVerbatim("EcalCalib") << "DetId not is e25";
571  continue;
572  }
573 
574  if (edm::isNotFinite(itrechit->energy()))
575  return;
576  energy.push_back(itrechit->energy());
577  energy5x5 += energy[icry];
578 
579  if (icry == 6 || icry == 7 || icry == 8 || icry == 11 || icry == 12 || icry == 13 || icry == 16 || icry == 17 ||
580  icry == 18) {
581  energy3x3 += energy[icry];
582  }
583  }
584  if ((int)energy.size() != ClusterSize_ * ClusterSize_)
585  return;
586  // Once we have the matrix 5x5, we have to correct for gaps/cracks/umbrella and maincontainement
587 
588  GeneralMap->Fill(maxCC_Eta, maxCC_Phi);
589 
590  EoP_all->Fill(highPtElectron.eSuperClusterOverP());
591 
592  if (highPtElectron.classification() == elecclass_ || elecclass_ == -1) {
593  float Ptrack_in =
594  sqrt(pow(highPtElectron.trackMomentumAtVtx().X(), 2) + pow(highPtElectron.trackMomentumAtVtx().Y(), 2) +
595  pow(highPtElectron.trackMomentumAtVtx().Z(), 2));
596 
597  float UncorrectedPatCalo =
598  sqrt(pow(highPtElectron.trackMomentumAtCalo().X(), 2) + pow(highPtElectron.trackMomentumAtCalo().Y(), 2) +
599  pow(highPtElectron.trackMomentumAtCalo().Z(), 2));
600 
601  float Ptrack_out =
602  sqrt(pow(highPtElectron.trackMomentumOut().X(), 2) + pow(highPtElectron.trackMomentumOut().Y(), 2) +
603  pow(highPtElectron.trackMomentumOut().Z(), 2));
604 
605  EventMatrixNoCuts.push_back(energy);
606  EnergyVectorNoCuts.push_back(UncorrectedPatCalo);
607 
608  MaxCCetaNoCuts.push_back(maxCC_Eta);
609  MaxCCphiNoCuts.push_back(maxCC_Phi);
610 
611  WeightVectorNoCuts.push_back(energy5x5 / UncorrectedPatCalo);
612 
613  //---------------------------------------------------No Cuts-------------------------------------------------------
614  e9NoCuts->Fill(energy3x3);
615  e25NoCuts->Fill(energy5x5);
616  e9Overe25NoCuts->Fill(energy3x3 / energy5x5);
617  scENoCuts->Fill(sc.energy());
618 
619  trPNoCuts->Fill(UncorrectedPatCalo);
620 
621  EoPNoCuts->Fill(highPtElectron.eSuperClusterOverP());
622  e25OverScENoCuts->Fill(energy5x5 / sc.energy());
623 
624  E25oPNoCuts->Fill(energy5x5 / UncorrectedPatCalo);
625 
626  MapNoCuts->Fill(maxCC_Eta, maxCC_Phi);
627  PinOverPoutNoCuts->Fill(
628  sqrt(pow(highPtElectron.trackMomentumAtVtx().X(), 2) + pow(highPtElectron.trackMomentumAtVtx().Y(), 2) +
629  pow(highPtElectron.trackMomentumAtVtx().Z(), 2)) /
630  sqrt(pow(highPtElectron.trackMomentumOut().X(), 2) + pow(highPtElectron.trackMomentumOut().Y(), 2) +
631  pow(highPtElectron.trackMomentumOut().Z(), 2)));
632  eSeedOverPoutNoCuts->Fill(highPtElectron.eSuperClusterOverP());
633 
634  MapCor1NoCuts->Fill(energy5x5 / UncorrectedPatCalo, energy5x5 / Ptrack_in);
635  MapCor2NoCuts->Fill(energy5x5 / UncorrectedPatCalo, highPtElectron.eSuperClusterOverP());
636  MapCor3NoCuts->Fill(energy5x5 / UncorrectedPatCalo, Ptrack_out / Ptrack_in);
637  MapCor4NoCuts->Fill(energy5x5 / UncorrectedPatCalo, energy5x5 / highPtElectron.p());
638  MapCor5NoCuts->Fill(energy5x5 / UncorrectedPatCalo, UncorrectedPatCalo / Ptrack_out);
639  MapCor6NoCuts->Fill(Ptrack_out / Ptrack_in, energy5x5 / Ptrack_in);
640  MapCor7NoCuts->Fill(Ptrack_out / Ptrack_in, UncorrectedPatCalo / Ptrack_out);
641  MapCor8NoCuts->Fill(energy5x5 / Ptrack_in, UncorrectedPatCalo / Ptrack_out);
642  MapCor9NoCuts->Fill(energy5x5 / UncorrectedPatCalo, highPtElectron.eSeedClusterOverPout());
643  MapCor10NoCuts->Fill(highPtElectron.eSeedClusterOverPout(), Ptrack_out / Ptrack_in);
644  MapCor11NoCuts->Fill(highPtElectron.eSeedClusterOverPout(), energy5x5 / Ptrack_in);
645 
646  PinMinPoutNoCuts->Fill((Ptrack_in - Ptrack_out) / Ptrack_in);
647 
648  Error1NoCuts->Fill(highPtElectron.trackMomentumError() / Ptrack_in);
649  Error2NoCuts->Fill(highPtElectron.trackMomentumError() / Ptrack_out);
650  Error3NoCuts->Fill(highPtElectron.trackMomentumError() / UncorrectedPatCalo);
651  eSeedOverPout2NoCuts->Fill(highPtElectron.eSeedClusterOverPout());
652 
653  hadOverEmNoCuts->Fill(highPtElectron.hadronicOverEm());
654 
655  //------------------------------------------------Cuts-----------------------------------------------------
656  //Cuts!
657  if ((energy3x3 / energy5x5) < cut1_)
658  return;
659 
660  if ((Ptrack_out / Ptrack_in) < cut2_ || (Ptrack_out / Ptrack_in) > cut3_)
661  return;
662  if ((energy5x5 / Ptrack_in) < cutEPin1_ || (energy5x5 / Ptrack_in) > cutEPin2_)
663  return;
664 
665  e9->Fill(energy3x3);
666  e25->Fill(energy5x5);
667  e9Overe25->Fill(energy3x3 / energy5x5);
668  scE->Fill(sc.energy());
669  trP->Fill(UncorrectedPatCalo);
670 
671  EoP->Fill(highPtElectron.eSuperClusterOverP());
672  e25OverScE->Fill(energy5x5 / sc.energy());
673 
674  E25oP->Fill(energy5x5 / UncorrectedPatCalo);
675 
676  Map->Fill(maxCC_Eta, maxCC_Phi);
677  PinOverPout->Fill(
678  sqrt(pow(highPtElectron.trackMomentumAtVtx().X(), 2) + pow(highPtElectron.trackMomentumAtVtx().Y(), 2) +
679  pow(highPtElectron.trackMomentumAtVtx().Z(), 2)) /
680  sqrt(pow(highPtElectron.trackMomentumOut().X(), 2) + pow(highPtElectron.trackMomentumOut().Y(), 2) +
681  pow(highPtElectron.trackMomentumOut().Z(), 2)));
682  eSeedOverPout->Fill(highPtElectron.eSuperClusterOverP());
683 
684  MapCor1->Fill(energy5x5 / UncorrectedPatCalo, energy5x5 / Ptrack_in);
685  MapCor2->Fill(energy5x5 / UncorrectedPatCalo, highPtElectron.eSuperClusterOverP());
686  MapCor3->Fill(energy5x5 / UncorrectedPatCalo, Ptrack_out / Ptrack_in);
687  MapCor4->Fill(energy5x5 / UncorrectedPatCalo, energy5x5 / highPtElectron.p());
688  MapCor5->Fill(energy5x5 / UncorrectedPatCalo, UncorrectedPatCalo / Ptrack_out);
689  MapCor6->Fill(Ptrack_out / Ptrack_in, energy5x5 / Ptrack_in);
690  MapCor7->Fill(Ptrack_out / Ptrack_in, UncorrectedPatCalo / Ptrack_out);
691  MapCor8->Fill(energy5x5 / Ptrack_in, UncorrectedPatCalo / Ptrack_out);
692  MapCor9->Fill(energy5x5 / UncorrectedPatCalo, highPtElectron.eSeedClusterOverPout());
693  MapCor10->Fill(highPtElectron.eSeedClusterOverPout(), Ptrack_out / Ptrack_in);
694  MapCor11->Fill(highPtElectron.eSeedClusterOverPout(), energy5x5 / Ptrack_in);
695 
696  PinMinPout->Fill((Ptrack_in - Ptrack_out) / Ptrack_in);
697 
698  Error1->Fill(highPtElectron.trackMomentumError() / Ptrack_in);
699  Error2->Fill(highPtElectron.trackMomentumError() / Ptrack_out);
700  Error3->Fill(highPtElectron.trackMomentumError() / UncorrectedPatCalo);
701 
702  eSeedOverPout2->Fill(highPtElectron.eSeedClusterOverPout());
703  hadOverEm->Fill(highPtElectron.hadronicOverEm());
704 
705  EventMatrix.push_back(energy);
706  EnergyVector.push_back(UncorrectedPatCalo);
707  MaxCCeta.push_back(maxCC_Eta);
708  MaxCCphi.push_back(maxCC_Phi);
709 
710  WeightVector.push_back(energy5x5 / UncorrectedPatCalo);
711 
712  //-------------------------------------------------------Extra Cut-----------------------------------------------------
713  if (highPtElectron.eSeedClusterOverPout() < cutESeed_)
714  return;
715 
716  MapCor1ESeed->Fill(energy5x5 / UncorrectedPatCalo, energy5x5 / Ptrack_in);
717  MapCor2ESeed->Fill(energy5x5 / UncorrectedPatCalo, highPtElectron.eSuperClusterOverP());
718  MapCor3ESeed->Fill(energy5x5 / UncorrectedPatCalo, Ptrack_out / Ptrack_in);
719  MapCor4ESeed->Fill(energy5x5 / UncorrectedPatCalo, energy5x5 / highPtElectron.p());
720  MapCor5ESeed->Fill(energy5x5 / UncorrectedPatCalo, UncorrectedPatCalo / Ptrack_out);
721  MapCor6ESeed->Fill(Ptrack_out / Ptrack_in, energy5x5 / Ptrack_in);
722  MapCor7ESeed->Fill(Ptrack_out / Ptrack_in, UncorrectedPatCalo / Ptrack_out);
723  MapCor8ESeed->Fill(energy5x5 / Ptrack_in, UncorrectedPatCalo / Ptrack_out);
724  MapCor9ESeed->Fill(energy5x5 / UncorrectedPatCalo, highPtElectron.eSeedClusterOverPout());
725  MapCor10ESeed->Fill(highPtElectron.eSeedClusterOverPout(), Ptrack_out / Ptrack_in);
726  MapCor11ESeed->Fill(highPtElectron.eSeedClusterOverPout(), energy5x5 / Ptrack_in);
727 
728  eSeedOverPout2ESeed->Fill(highPtElectron.eSeedClusterOverPout());
729 
730  hadOverEmESeed->Fill(highPtElectron.hadronicOverEm());
731 
732  } else {
733  return;
734  }
735 }
Log< level::Info, true > LogVerbatim
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
CalibrationCluster::CalibMap ReducedMap
edm::EDGetTokenT< EBRecHitCollection > recHitToken_
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:209
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
std::vector< float > WeightVector
int iphi() const
get the crystal iphi
Definition: EBDetId.h:51
float trackMomentumError() const
Definition: GsfElectron.h:884
void analyze(const edm::Event &, const edm::EventSetup &) override
HouseholderDecomposition * MyHH
T const * product() const
Definition: Handle.h:70
std::vector< EcalRecHit >::const_iterator const_iterator
edm::EDGetTokenT< reco::GsfElectronCollection > electronToken_
CalibMap getMap(int, int, int, int)
std::vector< float > solutionNoCuts
MinL3Algorithm * MyL3Algo1
void beginJob() override
std::vector< int > MaxCCphiNoCuts
std::vector< int > MaxCCetaNoCuts
Log< level::Error, false > LogError
std::vector< GsfElectron > GsfElectronCollection
collection of GsfElectron objects
float eSuperClusterOverP() const
Definition: GsfElectron.h:221
Classification classification() const
Definition: GsfElectron.h:805
std::vector< float > solution
int ieta() const
get the crystal ieta
Definition: EBDetId.h:49
math::XYZVectorF trackMomentumAtCalo() const
Definition: GsfElectron.h:269
float eSeedClusterOverPout() const
Definition: GsfElectron.h:223
constexpr bool null() const
is this a null id ?
Definition: DetId.h:59
std::vector< float > newCalibs
math::XYZVectorF trackMomentumOut() const
Definition: GsfElectron.h:270
std::vector< float > EnergyVector
int iEvent
Definition: GenABIO.cc:224
double p() const final
magnitude of momentum vector
ElectronCalibration(const edm::ParameterSet &)
TupleMultiplicity< TrackerTraits > const *__restrict__ TrackingRecHitSoAConstView< TrackerTraits > TrackerTraits::tindex_type *__restrict__ double *__restrict__ phits
EBDetId findMaxHit(edm::Handle< EBRecHitCollection > &)
std::vector< float > iterate(const std::vector< std::vector< float > > &eventMatrix, const std::vector< int > &VmaxCeta, const std::vector< int > &VmaxCphi, const std::vector< float > &energyVector, const int &nIter, const bool &normalizeFlag=false)
std::vector< int > MaxCCphi
T sqrt(T t)
Definition: SSEVec.h:19
std::vector< float > WeightVectorNoCuts
std::vector< std::vector< float > > EventMatrixNoCuts
CalibrationCluster calibCluster
math::XYZVectorF trackMomentumAtVtx() const
Definition: GsfElectron.h:268
const_iterator begin() const
static const int ETAPHIMODE
Definition: EBDetId.h:158
ii
Definition: cuy.py:589
std::vector< EcalRecHit >::iterator iterator
const_iterator end() const
double energy() const
cluster energy
Definition: CaloCluster.h:148
std::vector< std::vector< float > > EventMatrix
std::vector< float > runRegional(const std::vector< std::vector< float > > &eventMatrix, const std::vector< int > &VmaxCeta, const std::vector< int > &VmaxCphi, const std::vector< float > &energyVector, const int &nIter, const int &regLength=5)
bool isValid() const
Definition: HandleBase.h:70
EBDetId findMaxHit2(const std::vector< DetId > &, const EBRecHitCollection *)
void writeLine(EBDetId const &, float)
HLT enums.
std::vector< float > EnergyVectorNoCuts
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:180
std::vector< EBDetId > get5x5Id(EBDetId const &)
save
Definition: cuy.py:1164
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
float hadronicOverEm() const
Definition: GsfElectron.h:500
std::vector< float > oldCalibs
SuperClusterRef superCluster() const override
reference to a SuperCluster
Definition: GsfElectron.h:155
std::vector< float > iterate(const std::vector< std::vector< float > > &eventMatrix, const std::vector< int > &VmaxCeta, const std::vector< int > &VmaxCphi, const std::vector< float > &energyVector, const int &nIter, const bool &normalizeFlag=false)
double eta() const final
momentum pseudorapidity
std::vector< int > MaxCCeta