CMS 3D CMS Logo

InvMatrixUtils.cc
Go to the documentation of this file.
1 
6 #include "TStyle.h"
7 #include "TROOT.h"
8 #include "CLHEP/Geometry/Point3D.h"
9 //#include "ConfigParser.h"
11 
12 #include <iostream>
13 #include <string>
14 #include <fstream>
15 #include <sstream>
16 #include <vector>
17 #include <cassert>
18 
20 void setStyle() {
21  gROOT->SetStyle("Plain");
22  gStyle->SetTextSize(0.5);
23  //gStyle->SetOptStat (1111111) ;
24  gStyle->SetOptStat(0);
25  //gStyle->SetOptFit (1111) ;
26  gStyle->SetOptFit(0);
27  gStyle->SetTitleBorderSize(0);
28  gStyle->SetTitleX(0.08);
29  gStyle->SetTitleY(0.97);
30  gStyle->SetPalette(1, nullptr);
31  gStyle->SetStatColor(0);
32  gStyle->SetFrameFillStyle(0);
33  gStyle->SetFrameFillColor(0);
34  return;
35 }
36 
37 // -----------------------------------------------------------------
38 
40  setStyle();
41  TCanvas *globalCanvas = static_cast<TCanvas *>(gROOT->FindObject(name.c_str()));
42  if (globalCanvas) {
43  globalCanvas->Clear();
44  globalCanvas->UseCurrentStyle();
45  globalCanvas->SetWindowSize(700, 600);
46 
47  } else {
48  globalCanvas = new TCanvas(name.c_str(), name.c_str(), 700, 600);
49  }
50  return globalCanvas;
51 }
52 
53 // -----------------------------------------------------------------
54 
56  // std::cout << "writing " << name << std::endl ;
57  // setStyle () ;
58  TFile *globalTFile = (TFile *)gROOT->FindObject(name.c_str());
59  if (!globalTFile) {
60  // std::cout << "does not exist. creating it " << std::endl;
61  globalTFile = new TFile(name.c_str(), "RECREATE");
62  }
63 
64  return globalTFile;
65 }
66 
67 // -----------------------------------------------------------------
68 
70  TFile *globalTFile = static_cast<TFile *>(gROOT->FindObject(name.c_str()));
71  if (!globalTFile)
72  return 1;
73  globalTFile->Write();
74  globalTFile->Close();
75  delete globalTFile;
76  return 0;
77 }
78 
79 // -----------------------------------------------------------------
80 
81 CLHEP::HepMatrix *getSavedMatrix(const std::string &name) {
83  CLHEP::HepMatrix *savedMatrix;
84  if (reader.touch(name)) {
85  savedMatrix = static_cast<CLHEP::HepMatrix *>(reader.getMatrix(name));
86  } else {
87  savedMatrix = new CLHEP::HepMatrix(SCMaxEta, SCMaxPhi, 0);
88  }
89 
90  return savedMatrix;
91 }
92 
93 //========================================================================
94 //da usare
95 //HepGeom::Point3D<double> TBimpactPoint = TBposition (amplitude,m_beamEnergy) ;
96 
97 /*
98 dove trovare il codice di Chiara in CMSSW, per la ricostruzione
99 della posizione:
100 
101 http://cmssw.cvs.cern.ch/cgi-bin/cmssw.cgi/CMSSW/RecoTBCalo/EcalTBAnalysisCoreTools/src/TBPositionCalc.cc?rev=1.1&cvsroot=CMSSW&content-type=text/vnd.viewcvs-markup
102 
103 */
104 
106 HepGeom::Point3D<Float_t> TBposition(const Float_t amplit[7][7],
107  const Float_t beamEne,
108  const Float_t w0,
109  const Float_t x0,
110  const Float_t a0,
111  const Float_t sideX, // crystal geometry, in mm
112  const Float_t sideY) {
113  // variables
114  Float_t caloX = 0.;
115  Float_t caloY = 0.;
116  Float_t sumWeight = 0.;
117  Float_t depth = x0 * (log(beamEne) + a0); // shower depthh, in mm
118  Float_t sin3 = 0.052335956; // sin (3 degrees) , sin3 = sin(3.*3.141592654/180.)
119 
120  Float_t invE3x3 = 1. / get3x3(amplit);
121 
122  // loop over 3x3 crystals
123  for (int eta = 2; eta <= 4; eta++) {
124  for (int phi = 2; phi <= 4; phi++) {
125  Float_t weight = log(amplit[eta][phi] * invE3x3) + w0;
126  if (weight > 0) {
127  caloX += (eta - 3) * sideX * weight;
128  caloY -= (phi - 3) * sideY * weight;
129  sumWeight += weight;
130  }
131  }
132  }
133 
134  caloX /= sumWeight;
135  caloY /= sumWeight;
136 
137  // correction for shower depthh
138  caloX -= depth * sin3;
139  caloY -= depth * sin3;
140 
141  // FIXME the z set to zero
142  HepGeom::Point3D<Float_t> TBposition(caloX, caloY, 0);
143 
144  return TBposition;
145 }
146 
147 // -----------------------------------------------------------------
148 
151 double get5x5(const Float_t energy[7][7]) {
152  double total = 0.;
153 
154  for (int eta = 1; eta < 6; ++eta)
155  for (int phi = 1; phi < 6; ++phi)
156  total += energy[eta][phi];
157 
158  return total;
159 }
160 
161 // -----------------------------------------------------------------
162 
165 double get3x3(const Float_t energy[7][7]) {
166  double total = 0.;
167 
168  for (int eta = 2; eta < 5; ++eta)
169  for (int phi = 2; phi < 5; ++phi)
170  total += energy[eta][phi];
171 
172  return total;
173 }
174 
175 // -----------------------------------------------------------------
176 
178 /* int parseConfigFile (const TString& config)
179 {
180  if (gConfigParser) return 1 ;
181 
182  std::cout << "parsing "
183  << config << " file"
184  << std::endl ;
185 
186 
187  gConfigParser = new ConfigParser () ;
188  if ( !(gConfigParser->init(config)) )
189  {
190  std::cout << "Analysis::parseConfigFile: Could not open configuration file "
191  << config << std::endl;
192  perror ("Analysis::parseConfigFile: ") ;
193  exit (-1) ;
194  }
195  gConfigParser->print () ;
196  return 0 ;
197 }*/
198 
199 // -----------------------------------------------------------------
200 
201 // per un certo eta, il cristallo puo' essere qualsiasi intero
202 // Calibrationtra xmin e xmax
203 // lo posso fissare solo sfruttando anche phi
204 int xtalFromEtaPhi(const int &myEta, const int &myPhi) {
205  int xMin = 20 * myEta + 1;
206  int xMax = 20 * (myEta + 1) + 1;
207 
208  int myCryst = 999999;
209 
210  for (int x = xMin; x < xMax; x++) {
211  if (phiFromXtal(x) == myPhi)
212  myCryst = x;
213  }
214  return myCryst;
215 }
216 
217 // -----------------------------------------------------------------
218 
219 int xtalFromiEtaiPhi(const int &iEta, const int &iPhi) {
220  assert(iEta >= 1);
221  assert(iEta <= 85);
222  assert(iPhi >= 1);
223  assert(iPhi <= 20);
224  return 20 * (iEta - 1) + 21 - iPhi;
225 }
226 
227 // -----------------------------------------------------------------
228 
229 int etaFromXtal(const int &xtal) {
230  // return floor (static_cast<double> ((xtal-1) / 20)) ;
231  return int(floor((xtal - 1) / 20));
232 }
233 
234 // -----------------------------------------------------------------
235 
236 int phiFromXtal(const int &xtal) {
237  int phi = (xtal - 1) - 20 * etaFromXtal(xtal);
238  return (20 - phi - 1);
239 }
240 
241 // -----------------------------------------------------------------
242 
243 int ietaFromXtal(const int &xtal) { return etaFromXtal(xtal) + 1; }
244 
245 // -----------------------------------------------------------------
246 
247 int iphiFromXtal(const int &xtal) { return phiFromXtal(xtal) + 1; }
248 
249 // -----------------------------------------------------------------
250 
251 int extract(std::vector<int> *output, const std::string &dati) {
252  std::ifstream _dati(dati.c_str());
253  // loop over the file
254  while (!_dati.eof()) {
255  // get the line
256  std::string dataline;
257  do {
258  getline(_dati, dataline, '\n');
259  } while (*dataline.begin() == '#');
260  std::stringstream linea(dataline);
261  // loop over the line
262  while (!linea.eof()) {
263  int buffer = -1;
264  linea >> buffer;
265  if (buffer != -1)
266  output->push_back(buffer);
267  } // loop over the line
268  } // loop over the file
269  return output->size();
270 }
271 
272 // -----------------------------------------------------------------
273 
274 // FIXME questi eta, phi sono quelli della matrice CLHEP,
275 // FIXME non quelli del super-modulo, giusto?
276 int writeCalibTxt(const CLHEP::HepMatrix &AmplitudeMatrix,
277  const CLHEP::HepMatrix &SigmaMatrix,
278  const CLHEP::HepMatrix &StatisticMatrix,
280  // look for the reference crystal
281  double reference = 0.;
282  for (int eta = 0; eta < SCMaxEta; ++eta)
283  for (int phi = 0; phi < SCMaxPhi; ++phi) {
284  if (AmplitudeMatrix[eta][phi] && SigmaMatrix[eta][phi] < 100 /*FIXME sigmaCut*/) {
285  reference = AmplitudeMatrix[eta][phi];
286  std::cout << "[InvMatrixUtils][writeCalibTxt] reference crystal: "
287  << "(" << eta << "," << phi << ") -> " << reference << "\n";
288  break;
289  }
290  }
291  if (!reference) {
292  std::cerr << "ERROR: no calibration coefficients found" << std::endl;
293  return 1;
294  }
295 
296  // open the file for output
297  std::ofstream txt_outfile;
298  txt_outfile.open(fileName.c_str());
299  txt_outfile << "# xtal\tcoeff\tsigma\tevt\tisGood\n";
300 
301  // loop over the crystals
302  for (int eta = 0; eta < SCMaxEta; ++eta)
303  for (int phi = 0; phi < SCMaxPhi; ++phi) {
304  int isGood = 1;
305  if (AmplitudeMatrix[eta][phi] == 0)
306  isGood = 0;
307  if (SigmaMatrix[eta][phi] > 100 /*FIXME sigmaCut*/)
308  isGood = 0;
309  txt_outfile << xtalFromEtaPhi(eta, phi) << "\t" << AmplitudeMatrix[eta][phi] / reference << "\t"
310  << SigmaMatrix[eta][phi] << "\t" << StatisticMatrix[eta][phi] << "\t" << isGood << "\n";
311  }
312 
313  // save and close the file
314  txt_outfile.close();
315  return 0;
316 }
317 
318 // -----------------------------------------------------------------
319 
320 int writeCMSSWCoeff(const CLHEP::HepMatrix &amplMatrix,
321  double calibThres,
322  float ERef,
323  const CLHEP::HepMatrix &sigmaMatrix,
324  const CLHEP::HepMatrix &statisticMatrix,
326  std::string genTag,
329  std::string type) {
330  // open the file for output
331  std::ofstream txt_outfile;
332  txt_outfile.open(fileName.c_str());
333  txt_outfile << "1\n"; // super-module number
334  txt_outfile << "-1\n"; // number of events
335  txt_outfile << genTag << "\n";
336  txt_outfile << method << "\n";
337  txt_outfile << version << "\n";
338  txt_outfile << type << "\n";
339 
340  double reference = ERef;
341 
342  // loop over crystals
343  for (int eta = 0; eta < SCMaxEta; ++eta)
344  for (int phi = 0; phi < SCMaxPhi; ++phi) {
345  if (amplMatrix[eta][phi] <= calibThres)
346  txt_outfile << xtalFromiEtaiPhi(eta + 1, phi + 1) << "\t" << 1 << "\t" << -1 << "\t" << -1 << "\t" << 0 << "\n";
347  else
348  txt_outfile << xtalFromiEtaiPhi(eta + 1, phi + 1) << "\t" << reference / amplMatrix[eta][phi] << "\t"
349  << sigmaMatrix[eta][phi] << "\t" << statisticMatrix[eta][phi] << "\t" << 1 << "\n";
350  } // loop over crystals
351 
352  // save and close the file
353  txt_outfile.close();
354  return 0;
355 }
356 
357 // -----------------------------------------------------------------
358 
359 int writeCMSSWCoeff(const CLHEP::HepMatrix &amplMatrix,
360  double calibThres,
361  int etaRef,
362  int phiRef,
363  const CLHEP::HepMatrix &sigmaMatrix,
364  const CLHEP::HepMatrix &statisticMatrix,
366  std::string genTag,
369  std::string type) {
370  // open the file for output
371  std::ofstream txt_outfile;
372  txt_outfile.open(fileName.c_str());
373  txt_outfile << "1\n"; // super-module number
374  txt_outfile << "-1\n"; // number of events
375  txt_outfile << genTag << "\n";
376  txt_outfile << method << "\n";
377  txt_outfile << version << "\n";
378  txt_outfile << type << "\n";
379 
380  if (amplMatrix[etaRef - 1][phiRef - 1] == 0) {
381  std::cerr << "The reference crystal: (" << etaRef << "," << phiRef << ") is out of range\n";
382  return 1;
383  }
384  double reference = amplMatrix[etaRef - 1][phiRef - 1];
385 
386  // loop over crystals
387  for (int eta = 0; eta < SCMaxEta; ++eta)
388  for (int phi = 0; phi < SCMaxPhi; ++phi) {
389  if (amplMatrix[eta][phi] <= calibThres)
390  txt_outfile << xtalFromiEtaiPhi(eta + 1, phi + 1) << "\t" << 1 << "\t" << -1 << "\t" << -1 << "\t" << 0 << "\n";
391  else
392  txt_outfile << xtalFromiEtaiPhi(eta + 1, phi + 1) << "\t" << reference / amplMatrix[eta][phi] << "\t"
393  << sigmaMatrix[eta][phi] << "\t" << statisticMatrix[eta][phi] << "\t" << 1 << "\n";
394  } // loop over crystals
395 
396  // save and close the file
397  txt_outfile.close();
398  return 0;
399 }
400 
401 // -----------------------------------------------------------------
402 
403 int translateCoeff(const CLHEP::HepMatrix &calibcoeff,
404  const CLHEP::HepMatrix &sigmaMatrix,
405  const CLHEP::HepMatrix &statisticMatrix,
406  std::string SMnumber,
407  double calibThres,
409  std::string genTag,
412  std::string type) {
413  // open the file for output
414  std::ofstream txt_outfile;
415  txt_outfile.open(fileName.c_str());
416  txt_outfile << SMnumber << "\n"; // super-module number
417  txt_outfile << "-1\n"; // number of events
418  txt_outfile << genTag << "\n";
419  txt_outfile << method << "\n";
420  txt_outfile << version << "\n";
421  txt_outfile << type << "\n";
422 
423  // loop over crystals
424  for (int eta = 0; eta < SCMaxEta; ++eta)
425  for (int phi = 0; phi < SCMaxPhi; ++phi) {
426  if (calibcoeff[eta][phi] < calibThres) {
427  txt_outfile << xtalFromiEtaiPhi(eta + 1, phi + 1) << "\t" << 1 << "\t" << -1 << "\t" << -1 << "\t" << 0 << "\n";
428  std::cout << "[translateCoefff][" << SMnumber << "]\t WARNING crystal " << xtalFromiEtaiPhi(eta + 1, phi + 1)
429  << " calib coeff below threshold: "
430  << "\t" << 1 << "\t" << -1 << "\t" << -1 << "\t" << 0 << "\n";
431  } else
432  txt_outfile << xtalFromiEtaiPhi(eta + 1, phi + 1) << "\t" << calibcoeff[eta][phi] << "\t"
433  << sigmaMatrix[eta][phi] << "\t" << statisticMatrix[eta][phi] << "\t" << 1 << "\n";
434  } // loop over crystals
435 
436  // save and close the file
437  txt_outfile.close();
438  return 0;
439 }
440 
441 // -----------------------------------------------------------------
442 
443 int readCMSSWcoeff(CLHEP::HepMatrix &calibcoeff, const std::string &inputFileName, double defaultVal) {
444  std::ifstream CMSSWfile;
445  CMSSWfile.open(inputFileName.c_str());
447  CMSSWfile >> buffer;
448  CMSSWfile >> buffer;
449  CMSSWfile >> buffer;
450  CMSSWfile >> buffer;
451  CMSSWfile >> buffer;
452  CMSSWfile >> buffer;
453  while (!CMSSWfile.eof()) {
454  int xtalnum;
455  CMSSWfile >> xtalnum;
456  double coeff;
457  CMSSWfile >> coeff;
458  double buffer;
459  CMSSWfile >> buffer;
460  int good;
461  CMSSWfile >> good;
462  CMSSWfile >> good;
463  if (!good)
464  coeff = defaultVal; //FIXME 0 o 1?
465  calibcoeff[etaFromXtal(xtalnum)][phiFromXtal(xtalnum)] = coeff;
466  }
467  return 0;
468 }
469 
470 // -----------------------------------------------------------------
471 
472 int readCMSSWcoeffForComparison(CLHEP::HepMatrix &calibcoeff, const std::string &inputFileName) {
473  std::ifstream CMSSWfile;
474  CMSSWfile.open(inputFileName.c_str());
476  CMSSWfile >> buffer;
477  CMSSWfile >> buffer;
478  CMSSWfile >> buffer;
479  CMSSWfile >> buffer;
480  CMSSWfile >> buffer;
481  CMSSWfile >> buffer;
482  while (!CMSSWfile.eof()) {
483  int xtalnum;
484  CMSSWfile >> xtalnum;
485  double coeff;
486  CMSSWfile >> coeff;
487  double buffer;
488  CMSSWfile >> buffer;
489  int good;
490  CMSSWfile >> good;
491  CMSSWfile >> good;
492  if (!good)
493  coeff = 0.; //FIXME 0 o 1?
494  calibcoeff[etaFromXtal(xtalnum)][phiFromXtal(xtalnum)] = coeff;
495  }
496  return 0;
497 }
498 
499 // -----------------------------------------------------------------
500 
501 TH1D *smartProfile(TH2F *strip, double width) {
502  TProfile *stripProfile = strip->ProfileX();
503 
504  // (from FitSlices of TH2.h)
505 
506  double xmin = stripProfile->GetXaxis()->GetXmin();
507  double xmax = stripProfile->GetXaxis()->GetXmax();
508  int profileBins = stripProfile->GetNbinsX();
509 
510  std::string name = strip->GetName();
511  name += "_smart";
512  TH1D *prof = new TH1D(name.c_str(), strip->GetTitle(), profileBins, xmin, xmax);
513 
514  int cut = 0; // minimum number of entries per fitted bin
515  int nbins = strip->GetXaxis()->GetNbins();
516  int binmin = 1;
517  int ngroup = 1; // bins per step
518  int binmax = nbins;
519 
520  // loop over the strip bins
521  for (int bin = binmin; bin <= binmax; bin += ngroup) {
522  TH1D *hpy = strip->ProjectionY("_temp", bin, bin + ngroup - 1, "e");
523  if (hpy == nullptr)
524  continue;
525  int nentries = Int_t(hpy->GetEntries());
526  if (nentries == 0 || nentries < cut) {
527  delete hpy;
528  continue;
529  }
530 
531  Int_t biny = bin + ngroup / 2;
532 
533  hpy->GetXaxis()->SetRangeUser(hpy->GetMean() - width * hpy->GetRMS(), hpy->GetMean() + width * hpy->GetRMS());
534  prof->Fill(strip->GetXaxis()->GetBinCenter(biny), hpy->GetMean());
535  prof->SetBinError(biny, hpy->GetRMS());
536 
537  delete hpy;
538  } // loop over the bins
539 
540  delete stripProfile;
541  return prof;
542 }
543 
544 // -----------------------------------------------------------------
545 
546 TH1D *smartGausProfile(TH2F *strip, double width) {
547  TProfile *stripProfile = strip->ProfileX();
548 
549  // (from FitSlices of TH2.h)
550 
551  double xmin = stripProfile->GetXaxis()->GetXmin();
552  double xmax = stripProfile->GetXaxis()->GetXmax();
553  int profileBins = stripProfile->GetNbinsX();
554 
555  std::string name = strip->GetName();
556  name += "_smartGaus";
557  TH1D *prof = new TH1D(name.c_str(), strip->GetTitle(), profileBins, xmin, xmax);
558 
559  int cut = 0; // minimum number of entries per fitted bin
560  int nbins = strip->GetXaxis()->GetNbins();
561  int binmin = 1;
562  int ngroup = 1; // bins per step
563  int binmax = nbins;
564 
565  // loop over the strip bins
566  for (int bin = binmin; bin <= binmax; bin += ngroup) {
567  TH1D *hpy = strip->ProjectionY("_temp", bin, bin + ngroup - 1, "e");
568  if (hpy == nullptr)
569  continue;
570  int nentries = Int_t(hpy->GetEntries());
571  if (nentries == 0 || nentries < cut) {
572  delete hpy;
573  continue;
574  }
575 
576  Int_t biny = bin + ngroup / 2;
577 
578  TF1 *gaussian =
579  new TF1("gaussian", "gaus", hpy->GetMean() - width * hpy->GetRMS(), hpy->GetMean() + width * hpy->GetRMS());
580  gaussian->SetParameter(1, hpy->GetMean());
581  gaussian->SetParameter(2, hpy->GetRMS());
582  hpy->Fit("gaussian", "RQL");
583 
584  hpy->GetXaxis()->SetRangeUser(hpy->GetMean() - width * hpy->GetRMS(), hpy->GetMean() + width * hpy->GetRMS());
585  prof->Fill(strip->GetXaxis()->GetBinCenter(biny), gaussian->GetParameter(1));
586  prof->SetBinError(biny, gaussian->GetParameter(2));
587 
588  delete gaussian;
589  delete hpy;
590  } // loop over the bins
591 
592  delete stripProfile;
593  return prof;
594 }
595 
596 // -----------------------------------------------------------------
597 
598 TH1D *smartError(TH1D *strip) {
599  double xmin = strip->GetXaxis()->GetXmin();
600  double xmax = strip->GetXaxis()->GetXmax();
601  int stripsBins = strip->GetNbinsX();
602 
603  std::string name = strip->GetName();
604  name += "_error";
605  TH1D *error = new TH1D(name.c_str(), strip->GetTitle(), stripsBins, xmin, xmax);
606 
607  int binmin = 1;
608  int ngroup = 1; // bins per step
609  int binmax = stripsBins;
610  for (int bin = binmin; bin <= binmax; bin += ngroup) {
611  double dummyError = strip->GetBinError(bin);
612  error->SetBinContent(bin, dummyError);
613  }
614  return error;
615 }
616 
617 // -----------------------------------------------------------------
618 
619 double effectiveSigma(TH1F &histogram, int vSteps) {
620  double totInt = histogram.Integral();
621  int maxBin = histogram.GetMaximumBin();
622  int maxBinVal = int(histogram.GetBinContent(maxBin));
623  int totBins = histogram.GetNbinsX();
624  double area = totInt;
625  double threshold = 0;
626  double vStep = maxBinVal / vSteps;
627  int leftBin = 1;
628  int rightBin = totBins - 1;
629  //loop over the vertical range
630  while (area / totInt > 0.683) {
631  threshold += vStep;
632  // loop toward the left
633  for (int back = maxBin; back > 0; --back) {
634  if (histogram.GetBinContent(back) < threshold) {
635  leftBin = back;
636  break;
637  }
638  } // loop toward the left
639 
640  // loop toward the right
641  for (int fwd = maxBin; fwd < totBins; ++fwd) {
642  if (histogram.GetBinContent(fwd) < threshold) {
643  rightBin = fwd;
644  break;
645  }
646  } // loop toward the right
647  area = histogram.Integral(leftBin, rightBin);
648  } //loop over the vertical range
649 
650  histogram.GetXaxis()->SetRange(leftBin, rightBin);
651  // double sigmaEff = histogram.GetRMS () ;
652  double halfWidthRange = 0.5 * (histogram.GetBinCenter(rightBin) - histogram.GetBinCenter(leftBin));
653  return halfWidthRange;
654 }
655 
656 // -----------------------------------------------------------------
657 
658 std::pair<int, int> findSupport(TH1F &histogram, double thres) {
659  int totBins = histogram.GetNbinsX();
660  if (thres >= histogram.GetMaximum())
661  return std::pair<int, int>(0, totBins);
662 
663  int leftBin = totBins - 1;
664  // search from left for the minimum
665  for (int bin = 1; bin < totBins; ++bin) {
666  if (histogram.GetBinContent(bin) > thres) {
667  leftBin = bin;
668  break;
669  }
670  } // search from left for the minimum
671  int rightBin = 1;
672  // search from right for the maximum
673  for (int bin = totBins - 1; bin > 0; --bin) {
674  if (histogram.GetBinContent(bin) > thres) {
675  rightBin = bin;
676  break;
677  }
678  } // search from right for the maximum
679  return std::pair<int, int>(leftBin, rightBin);
680 }
681 
682 // -----------------------------------------------------------------
683 
684 void mtrTransfer(double output[SCMaxEta][SCMaxPhi], CLHEP::HepMatrix *input, double Default) {
685  for (int eta = 0; eta < SCMaxEta; ++eta)
686  for (int phi = 0; phi < SCMaxPhi; ++phi) {
687  if ((*input)[eta][phi])
688  output[eta][phi] = (*input)[eta][phi];
689  else
690  output[eta][phi] = Default;
691  }
692  return;
693 }
694 
695 // -----------------------------------------------------------------
696 
697 double etaCorrE1E25(int eta) {
698  double p0 = 0.807883;
699  double p1 = 0.000182551;
700  double p2 = -5.76961e-06;
701  double p3 = 7.41903e-08;
702  double p4 = -2.25384e-10;
703 
704  double corr;
705  if (eta < 6)
706  corr = p0;
707  else
708  corr = p0 + p1 * eta + p2 * eta * eta + p3 * eta * eta * eta + p4 * eta * eta * eta * eta;
709  return corr / p0;
710 }
711 // -----------------------------------------------------------------
712 
713 double etaCorrE1E49(int eta) {
714  double p0 = 0.799895;
715  double p1 = 0.000235487;
716  double p2 = -8.26496e-06;
717  double p3 = 1.21564e-07;
718  double p4 = -4.83286e-10;
719 
720  double corr;
721  if (eta < 8)
722  corr = p0;
723  else
724  corr = p0 + p1 * eta + p2 * eta * eta + p3 * eta * eta * eta + p4 * eta * eta * eta * eta;
725  return corr / p0;
726 }
727 // -----------------------------------------------------------------
728 
729 double etaCorrE1E9(int eta) {
730  if (eta < 4)
731  return 1.0;
732  // grazie Paolo
733  double p0 = 0.834629;
734  double p1 = 0.00015254;
735  double p2 = -4.91784e-06;
736  double p3 = 6.54652e-08;
737  double p4 = -2.4894e-10;
738 
739  double corr;
740  if (eta < 6)
741  corr = p0;
742  else
743  corr = p0 + p1 * eta + p2 * eta * eta + p3 * eta * eta * eta + p4 * eta * eta * eta * eta;
744  return corr / p0;
745 }
double etaCorrE1E9(int eta)
int translateCoeff(const CLHEP::HepMatrix &calibcoeff, const CLHEP::HepMatrix &sigmaMatrix, const CLHEP::HepMatrix &statisticMatrix, std::string SMnumber, double calibThres, std::string fileName, std::string genTag, std::string method, std::string version, std::string type)
TCanvas * getGlobalCanvas(std::string name)
#define SCMaxEta
int iphiFromXtal(const int &xtal)
TH1D * smartError(TH1D *strip)
TH1D * smartProfile(TH2F *strip, double width)
int ietaFromXtal(const int &xtal)
int readCMSSWcoeffForComparison(CLHEP::HepMatrix &calibcoeff, const std::string &inputFileName)
#define SCMaxPhi
int extract(std::vector< int > *output, const std::string &dati)
HepGeom::Point3D< Float_t > TBposition(const Float_t amplit[7][7], const Float_t beamEne, const Float_t w0, const Float_t x0, const Float_t a0, const Float_t sideX, const Float_t sideY)
double effectiveSigma(TH1F &histogram, int vSteps)
Definition: weight.py:1
reader
Definition: DQM.py:105
constexpr unsigned int maxBin
int phiFromXtal(const int &xtal)
assert(be >=bs)
static std::string const input
Definition: EdmProvDump.cc:50
void mtrTransfer(double output[85][20], CLHEP::HepMatrix *input, double Default)
TFile * getGlobalTFile(std::string name)
TH1D * smartGausProfile(TH2F *strip, double width)
std::pair< int, int > findSupport(TH1F &histogram, double thres)
dictionary corr
save (read) CLHEP::HepMatrix to (from) text files
Definition: matrixSaver.h:22
double etaCorrE1E25(int eta)
int xtalFromEtaPhi(const int &myEta, const int &myPhi)
CLHEP::HepMatrix * getSavedMatrix(const std::string &name)
void setStyle()
double etaCorrE1E49(int eta)
int etaFromXtal(const int &xtal)
int writeCalibTxt(const CLHEP::HepMatrix &AmplitudeMatrix, const CLHEP::HepMatrix &SigmaMatrix, const CLHEP::HepMatrix &StatisticMatrix, std::string fileName)
static constexpr float a0
int xtalFromiEtaiPhi(const int &iEta, const int &iPhi)
constexpr bool Default
Definition: SoACommon.h:75
int readCMSSWcoeff(CLHEP::HepMatrix &calibcoeff, const std::string &inputFileName, double defaultVal)
float x
Definition: output.py:1
double get5x5(const Float_t energy[7][7])
int writeCMSSWCoeff(const CLHEP::HepMatrix &amplMatrix, double calibThres, float ERef, const CLHEP::HepMatrix &sigmaMatrix, const CLHEP::HepMatrix &statisticMatrix, std::string fileName, std::string genTag, std::string method, std::string version, std::string type)
double get3x3(const Float_t energy[7][7])
int saveGlobalTFile(std::string name)