CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 {
22  gROOT->SetStyle ("Plain") ;
23  gStyle->SetTextSize(0.5);
24  //gStyle->SetOptStat (1111111) ;
25  gStyle->SetOptStat (0) ;
26  //gStyle->SetOptFit (1111) ;
27  gStyle->SetOptFit (0) ;
28  gStyle->SetTitleBorderSize (0) ;
29  gStyle->SetTitleX (0.08) ;
30  gStyle->SetTitleY (0.97) ;
31  gStyle->SetPalette (1,0) ;
32  gStyle->SetStatColor (0) ;
33  gStyle->SetFrameFillStyle (0) ;
34  gStyle->SetFrameFillColor (0) ;
35  return ;
36 }
37 
38 
39 // -----------------------------------------------------------------
40 
41 
43 {
44  setStyle () ;
45  TCanvas * globalCanvas = static_cast<TCanvas*>
46  (gROOT->FindObject (name.c_str ())) ;
47  if (globalCanvas)
48  {
49  globalCanvas->Clear () ;
50  globalCanvas->UseCurrentStyle () ;
51  globalCanvas->SetWindowSize (700, 600) ;
52 
53  }
54  else
55  {
56  globalCanvas = new TCanvas (name.c_str (),name.c_str (), 700, 600) ;
57  }
58  return globalCanvas ;
59 }
60 
61 // -----------------------------------------------------------------
62 
63 
65 {
66 // std::cout << "writing " << name << std::endl ;
67 // setStyle () ;
68  TFile * globalTFile = (TFile*) gROOT->FindObject (name.c_str()) ;
69  if (!globalTFile)
70  {
71 // std::cout << "does not exist. creating it " << std::endl;
72  globalTFile = new TFile (name.c_str(),"RECREATE") ;
73  }
74 
75  return globalTFile ;
76 }
77 
78 
79 // -----------------------------------------------------------------
80 
81 
83 {
84  TFile * globalTFile = static_cast<TFile*>
85  (gROOT->FindObject (name.c_str ())) ;
86  if (!globalTFile) return 1 ;
87  globalTFile->Write () ;
88  globalTFile->Close () ;
89  delete globalTFile ;
90  return 0 ;
91 }
92 
93 
94 // -----------------------------------------------------------------
95 
96 
97 CLHEP::HepMatrix * getSavedMatrix (const std::string & name)
98 {
100  CLHEP::HepMatrix * savedMatrix ;
101  if (reader.touch (name))
102  {
103  savedMatrix = static_cast<CLHEP::HepMatrix *> (
104  reader.getMatrix (name)
105  );
106  }
107  else
108  {
109  savedMatrix = new CLHEP::HepMatrix (SCMaxEta,SCMaxPhi,0) ;
110  }
111 
112  return savedMatrix ;
113 }
114 
115 
116 
117 //========================================================================
118 //da usare
119 //HepGeom::Point3D<double> TBimpactPoint = TBposition (amplitude,m_beamEnergy) ;
120 
121 /*
122 dove trovare il codice di Chiara in CMSSW, per la ricostruzione
123 della posizione:
124 
125 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
126 
127 */
128 
130 HepGeom::Point3D<Float_t>
131 TBposition (const Float_t amplit[7][7],
132  const Float_t beamEne,
133  const Float_t w0,
134  const Float_t x0,
135  const Float_t a0,
136  const Float_t sideX, // crystal geometry, in mm
137  const Float_t sideY)
138 {
139  // variables
140  Float_t caloX = 0. ;
141  Float_t caloY = 0. ;
142  Float_t sumWeight = 0. ;
143  Float_t depth = x0 * (log (beamEne)+ a0) ; // shower depthh, in mm
144  Float_t sin3 = 0.052335956 ; // sin (3 degrees) , sin3 = sin(3.*3.141592654/180.)
145 
146  Float_t invE3x3 = 1. / get3x3 (amplit) ;
147 
148  // loop over 3x3 crystals
149  for (int eta = 2; eta <= 4; eta++)
150  {
151  for (int phi = 2; phi <= 4; phi++)
152  {
153  Float_t weight = log( amplit[eta][phi] * invE3x3) + w0 ;
154  if ( weight>0 )
155  {
156  caloX += (eta-3) * sideX * weight;
157  caloY -= (phi-3) * sideY * weight;
158  sumWeight += weight;
159  }
160  }
161  }
162 
163  caloX /=sumWeight;
164  caloY /=sumWeight;
165 
166  // correction for shower depthh
167  caloX -= depth*sin3;
168  caloY -= depth*sin3;
169 
170  // FIXME the z set to zero
171  HepGeom::Point3D<Float_t> TBposition (caloX, caloY, 0) ;
172 
173  return TBposition ;
174 
175 }
176 
177 // -----------------------------------------------------------------
178 
179 
182 double get5x5 (const Float_t energy[7][7])
183 {
184  double total = 0. ;
185 
186  for (int eta=1 ; eta<6 ; ++eta)
187  for (int phi=1 ; phi<6 ; ++phi)
188  total += energy[eta][phi] ;
189 
190  return total ;
191 }
192 
193 
194 // -----------------------------------------------------------------
195 
196 
199 double get3x3 (const Float_t energy[7][7])
200 {
201  double total = 0. ;
202 
203  for (int eta=2 ; eta<5 ; ++eta)
204  for (int phi=2 ; phi<5 ; ++phi)
205  total += energy[eta][phi] ;
206 
207  return total ;
208 }
209 
210 
211 // -----------------------------------------------------------------
212 
213 
215 /* int parseConfigFile (const TString& config)
216 {
217  if (gConfigParser) return 1 ;
218 
219  std::cout << "parsing "
220  << config << " file"
221  << std::endl ;
222 
223 
224  gConfigParser = new ConfigParser () ;
225  if ( !(gConfigParser->init(config)) )
226  {
227  std::cout << "Analysis::parseConfigFile: Could not open configuration file "
228  << config << std::endl;
229  perror ("Analysis::parseConfigFile: ") ;
230  exit (-1) ;
231  }
232  gConfigParser->print () ;
233  return 0 ;
234 }*/
235 
236 
237 // -----------------------------------------------------------------
238 
239 
240 // per un certo eta, il cristallo puo' essere qualsiasi intero
241 // Calibrationtra xmin e xmax
242 // lo posso fissare solo sfruttando anche phi
243 int xtalFromEtaPhi (const int & myEta, const int & myPhi)
244 {
245  int xMin = 20 * myEta + 1 ;
246  int xMax = 20 * (myEta + 1) + 1 ;
247 
248  int myCryst = 999999 ;
249 
250  for (int x = xMin ; x < xMax ; x++)
251  {
252  if (phiFromXtal (x) == myPhi)
253  myCryst = x ;
254  }
255  return myCryst ;
256 }
257 
258 
259 // -----------------------------------------------------------------
260 
261 
262 int xtalFromiEtaiPhi (const int & iEta, const int & iPhi)
263 {
264  assert (iEta >= 1) ;
265  assert (iEta <= 85) ;
266  assert (iPhi >= 1) ;
267  assert (iPhi <= 20) ;
268  return 20 * (iEta-1) + 21 - iPhi ;
269 }
270 
271 
272 // -----------------------------------------------------------------
273 
274 
275 int etaFromXtal (const int & xtal)
276 {
277 // return floor (static_cast<double> ((xtal-1) / 20)) ;
278  return int (floor ((xtal-1) / 20) );
279 }
280 
281 
282 // -----------------------------------------------------------------
283 
284 
285 int phiFromXtal (const int & xtal)
286 {
287  int phi = (xtal-1) - 20 * etaFromXtal (xtal) ;
288  return (20 - phi - 1) ;
289 }
290 
291 
292 // -----------------------------------------------------------------
293 
294 
295 int ietaFromXtal (const int & xtal)
296 {
297  return etaFromXtal (xtal) + 1 ;
298 }
299 
300 
301 // -----------------------------------------------------------------
302 
303 
304 int iphiFromXtal (const int & xtal)
305 {
306  return phiFromXtal (xtal) + 1 ;
307 }
308 
309 
310 // -----------------------------------------------------------------
311 
312 
313 int extract (std::vector<int> * output , const std::string & dati)
314  {
315  std::ifstream _dati (dati.c_str ()) ;
316  // loop over the file
317  while (!_dati.eof())
318  {
319  // get the line
320  std::string dataline ;
321  do { getline (_dati, dataline,'\n') ; }
322  while (*dataline.begin () == '#') ;
323  std::stringstream linea (dataline) ;
324  // loop over the line
325  while (!linea.eof ())
326  {
327  int buffer = -1 ;
328  linea >> buffer ;
329  if (buffer != -1) output->push_back (buffer) ;
330  } // loop over the line
331  } // loop over the file
332  return output->size () ;
333  }
334 
335 
336 // -----------------------------------------------------------------
337 
338 
339 // FIXME questi eta, phi sono quelli della matrice CLHEP,
340 // FIXME non quelli del super-modulo, giusto?
341 int writeCalibTxt (const CLHEP::HepMatrix & AmplitudeMatrix,
342  const CLHEP::HepMatrix & SigmaMatrix,
343  const CLHEP::HepMatrix & StatisticMatrix,
345 {
346  // look for the reference crystal
347  double reference = 0. ;
348  for (int eta = 0 ; eta<SCMaxEta ; ++eta)
349  for (int phi = 0 ; phi<SCMaxPhi ; ++phi)
350  {
351  if (AmplitudeMatrix[eta][phi] &&
352  SigmaMatrix[eta][phi] < 100 /*FIXME sigmaCut*/)
353  {
354  reference = AmplitudeMatrix[eta][phi] ;
355  std::cout << "[InvMatrixUtils][writeCalibTxt] reference crystal: "
356  << "(" << eta << "," << phi << ") -> "
357  << reference << "\n" ;
358  break ;
359  }
360  }
361  if (!reference)
362  {
363  std::cerr << "ERROR: no calibration coefficients found" << std::endl ;
364  return 1 ;
365  }
366 
367  // open the file for output
368  std::ofstream txt_outfile ;
369  txt_outfile.open (fileName.c_str ()) ;
370  txt_outfile << "# xtal\tcoeff\tsigma\tevt\tisGood\n" ;
371 
372  // loop over the crystals
373  for (int eta = 0 ; eta<SCMaxEta ; ++eta)
374  for (int phi = 0 ; phi<SCMaxPhi ; ++phi)
375  {
376  int isGood = 1 ;
377  if (AmplitudeMatrix[eta][phi] == 0) isGood = 0 ;
378  if (SigmaMatrix[eta][phi] > 100 /*FIXME sigmaCut*/) isGood = 0 ;
379  txt_outfile << xtalFromEtaPhi (eta,phi)
380  << "\t" << AmplitudeMatrix[eta][phi]/reference
381  << "\t" << SigmaMatrix[eta][phi]
382  << "\t" << StatisticMatrix[eta][phi]
383  << "\t" << isGood <<"\n" ;
384  }
385 
386  // save and close the file
387  txt_outfile.close () ;
388  return 0 ;
389 }
390 
391 
392 // -----------------------------------------------------------------
393 
394 
395 int writeCMSSWCoeff (const CLHEP::HepMatrix & amplMatrix,
396  double calibThres,
397  float ERef,
398  const CLHEP::HepMatrix & sigmaMatrix,
399  const CLHEP::HepMatrix & statisticMatrix,
401  std::string genTag,
404  std::string type)
405 {
406  // open the file for output
407  std::ofstream txt_outfile ;
408  txt_outfile.open (fileName.c_str ()) ;
409  txt_outfile << "1\n" ; // super-module number
410  txt_outfile << "-1\n" ; // number of events
411  txt_outfile << genTag << "\n" ;
412  txt_outfile << method << "\n" ;
413  txt_outfile << version << "\n" ;
414  txt_outfile << type << "\n" ;
415 
416  double reference = ERef ;
417 
418  // loop over crystals
419  for (int eta = 0 ; eta < SCMaxEta ; ++eta)
420  for (int phi = 0 ; phi < SCMaxPhi ; ++phi)
421  {
422  if (amplMatrix[eta][phi] <= calibThres)
423  txt_outfile << xtalFromiEtaiPhi (eta+1,phi+1)
424  << "\t" << 1
425  << "\t" << -1
426  << "\t" << -1
427  << "\t" << 0 <<"\n" ;
428  else
429  txt_outfile << xtalFromiEtaiPhi (eta+1,phi+1)
430  << "\t" << reference / amplMatrix[eta][phi]
431  << "\t" << sigmaMatrix[eta][phi]
432  << "\t" << statisticMatrix[eta][phi]
433  << "\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 
444 
445 int writeCMSSWCoeff (const CLHEP::HepMatrix & amplMatrix,
446  double calibThres,
447  int etaRef, int phiRef,
448  const CLHEP::HepMatrix & sigmaMatrix,
449  const CLHEP::HepMatrix & statisticMatrix,
451  std::string genTag,
454  std::string type)
455 {
456  // open the file for output
457  std::ofstream txt_outfile ;
458  txt_outfile.open (fileName.c_str ()) ;
459  txt_outfile << "1\n" ; // super-module number
460  txt_outfile << "-1\n" ; // number of events
461  txt_outfile << genTag << "\n" ;
462  txt_outfile << method << "\n" ;
463  txt_outfile << version << "\n" ;
464  txt_outfile << type << "\n" ;
465 
466  if (amplMatrix[etaRef-1][phiRef-1] == 0)
467  {
468  std::cerr << "The reference crystal: ("
469  << etaRef << "," << phiRef
470  << ") is out of range\n" ;
471  return 1 ;
472  }
473  double reference = amplMatrix[etaRef-1][phiRef-1] ;
474 
475  // loop over crystals
476  for (int eta = 0 ; eta < SCMaxEta ; ++eta)
477  for (int phi = 0 ; phi < SCMaxPhi ; ++phi)
478  {
479  if (amplMatrix[eta][phi] <= calibThres)
480  txt_outfile << xtalFromiEtaiPhi (eta+1,phi+1)
481  << "\t" << 1
482  << "\t" << -1
483  << "\t" << -1
484  << "\t" << 0 <<"\n" ;
485  else
486  txt_outfile << xtalFromiEtaiPhi (eta+1,phi+1)
487  << "\t" << reference / amplMatrix[eta][phi]
488  << "\t" << sigmaMatrix[eta][phi]
489  << "\t" << statisticMatrix[eta][phi]
490  << "\t" << 1 <<"\n" ;
491  } // loop over crystals
492 
493  // save and close the file
494  txt_outfile.close () ;
495  return 0 ;
496 }
497 
498 
499 // -----------------------------------------------------------------
500 
501 
502 int translateCoeff (const CLHEP::HepMatrix & calibcoeff,
503  const CLHEP::HepMatrix & sigmaMatrix,
504  const CLHEP::HepMatrix & statisticMatrix,
505  std::string SMnumber,
506  double calibThres,
508  std::string genTag,
511  std::string type)
512 {
513  // open the file for output
514  std::ofstream txt_outfile ;
515  txt_outfile.open (fileName.c_str ()) ;
516  txt_outfile << SMnumber << "\n" ; // super-module number
517  txt_outfile << "-1\n" ; // number of events
518  txt_outfile << genTag << "\n" ;
519  txt_outfile << method << "\n" ;
520  txt_outfile << version << "\n" ;
521  txt_outfile << type << "\n" ;
522 
523  // loop over crystals
524  for (int eta = 0 ; eta < SCMaxEta ; ++eta)
525  for (int phi = 0 ; phi < SCMaxPhi ; ++phi)
526  {
527  if (calibcoeff[eta][phi] < calibThres)
528  {
529  txt_outfile << xtalFromiEtaiPhi (eta+1,phi+1)
530  << "\t" << 1
531  << "\t" << -1
532  << "\t" << -1
533  << "\t" << 0 <<"\n" ;
534  std::cout << "[translateCoefff][" << SMnumber
535  << "]\t WARNING crystal " << xtalFromiEtaiPhi (eta+1,phi+1)
536  << " calib coeff below threshold: "
537  << "\t" << 1
538  << "\t" << -1
539  << "\t" << -1
540  << "\t" << 0 <<"\n" ;
541  }
542  else
543  txt_outfile << xtalFromiEtaiPhi (eta+1,phi+1)
544  << "\t" << calibcoeff[eta][phi]
545  << "\t" << sigmaMatrix[eta][phi]
546  << "\t" << statisticMatrix[eta][phi]
547  << "\t" << 1 <<"\n" ;
548  } // loop over crystals
549 
550  // save and close the file
551  txt_outfile.close () ;
552  return 0 ;
553 }
554 
555 
556 // -----------------------------------------------------------------
557 
558 
559 int readCMSSWcoeff (CLHEP::HepMatrix & calibcoeff,
560  const std::string & inputFileName,
561  double defaultVal)
562 {
563  std::ifstream CMSSWfile ;
564  CMSSWfile.open (inputFileName.c_str ()) ;
565  std::string buffer ;
566  CMSSWfile >> buffer ;
567  CMSSWfile >> buffer ;
568  CMSSWfile >> buffer ;
569  CMSSWfile >> buffer ;
570  CMSSWfile >> buffer ;
571  CMSSWfile >> buffer ;
572  while (!CMSSWfile.eof ())
573  {
574  int xtalnum ;
575  CMSSWfile >> xtalnum ;
576  double coeff ;
577  CMSSWfile >> coeff ;
578  double buffer ;
579  CMSSWfile >> buffer ;
580  int good ;
581  CMSSWfile >> good ;
582  CMSSWfile >> good ;
583  if (!good) coeff = defaultVal ; //FIXME 0 o 1?
584  calibcoeff[etaFromXtal (xtalnum)][phiFromXtal (xtalnum)] = coeff ;
585  }
586  return 0 ;
587 
588 }
589 
590 
591 // -----------------------------------------------------------------
592 
593 
594 int readCMSSWcoeffForComparison (CLHEP::HepMatrix & calibcoeff,
595  const std::string & inputFileName)
596 {
597  std::ifstream CMSSWfile ;
598  CMSSWfile.open (inputFileName.c_str ()) ;
599  std::string buffer ;
600  CMSSWfile >> buffer ;
601  CMSSWfile >> buffer ;
602  CMSSWfile >> buffer ;
603  CMSSWfile >> buffer ;
604  CMSSWfile >> buffer ;
605  CMSSWfile >> buffer ;
606  while (!CMSSWfile.eof ())
607  {
608  int xtalnum ;
609  CMSSWfile >> xtalnum ;
610  double coeff ;
611  CMSSWfile >> coeff ;
612  double buffer ;
613  CMSSWfile >> buffer ;
614  int good ;
615  CMSSWfile >> good ;
616  CMSSWfile >> good ;
617  if (!good) coeff = 0. ; //FIXME 0 o 1?
618  calibcoeff[etaFromXtal (xtalnum)][phiFromXtal (xtalnum)] = coeff ;
619  }
620  return 0 ;
621 
622 }
623 
624 
625 // -----------------------------------------------------------------
626 
627 
628 TH1D * smartProfile (TH2F * strip, double width)
629 {
630  TProfile * stripProfile = strip->ProfileX () ;
631 
632  // (from FitSlices of TH2.h)
633 
634  double xmin = stripProfile->GetXaxis ()->GetXmin () ;
635  double xmax = stripProfile->GetXaxis ()->GetXmax () ;
636  int profileBins = stripProfile->GetNbinsX () ;
637 
638  std::string name = strip->GetName () ;
639  name += "_smart" ;
640  TH1D * prof = new TH1D
641  (name.c_str (),strip->GetTitle (),profileBins,xmin,xmax) ;
642 
643  int cut = 0 ; // minimum number of entries per fitted bin
644  int nbins = strip->GetXaxis ()->GetNbins () ;
645  int binmin = 1 ;
646  int ngroup = 1 ; // bins per step
647  int binmax = nbins ;
648 
649  // loop over the strip bins
650  for (int bin=binmin ; bin<=binmax ; bin += ngroup)
651  {
652  TH1D *hpy = strip->ProjectionY ("_temp",bin,bin+ngroup-1,"e") ;
653  if (hpy == 0) continue ;
654  int nentries = Int_t (hpy->GetEntries ()) ;
655  if (nentries == 0 || nentries < cut) {delete hpy ; continue ;}
656 
657  Int_t biny = bin + ngroup/2 ;
658 
659  hpy->GetXaxis ()->SetRangeUser ( hpy->GetMean () - width * hpy->GetRMS (),
660  hpy->GetMean () + width * hpy->GetRMS ()) ;
661  prof->Fill (strip->GetXaxis ()->GetBinCenter (biny),
662  hpy->GetMean ()) ;
663  prof->SetBinError (biny,hpy->GetRMS()) ;
664 
665  delete hpy ;
666  } // loop over the bins
667 
668  delete stripProfile ;
669  return prof ;
670 }
671 
672 
673 // -----------------------------------------------------------------
674 
675 
676 TH1D * smartGausProfile (TH2F * strip, double width)
677 {
678  TProfile * stripProfile = strip->ProfileX () ;
679 
680  // (from FitSlices of TH2.h)
681 
682  double xmin = stripProfile->GetXaxis ()->GetXmin () ;
683  double xmax = stripProfile->GetXaxis ()->GetXmax () ;
684  int profileBins = stripProfile->GetNbinsX () ;
685 
686  std::string name = strip->GetName () ;
687  name += "_smartGaus" ;
688  TH1D * prof = new TH1D
689  (name.c_str (),strip->GetTitle (),profileBins,xmin,xmax) ;
690 
691  int cut = 0 ; // minimum number of entries per fitted bin
692  int nbins = strip->GetXaxis ()->GetNbins () ;
693  int binmin = 1 ;
694  int ngroup = 1 ; // bins per step
695  int binmax = nbins ;
696 
697  // loop over the strip bins
698  for (int bin=binmin ; bin<=binmax ; bin += ngroup)
699  {
700  TH1D *hpy = strip->ProjectionY ("_temp",bin,bin+ngroup-1,"e") ;
701  if (hpy == 0) continue ;
702  int nentries = Int_t (hpy->GetEntries ()) ;
703  if (nentries == 0 || nentries < cut) {delete hpy ; continue ;}
704 
705  Int_t biny = bin + ngroup/2 ;
706 
707  TF1 * gaussian = new TF1 ("gaussian","gaus", hpy->GetMean () - width * hpy->GetRMS (),
708  hpy->GetMean () + width * hpy->GetRMS ()) ;
709  gaussian->SetParameter (1,hpy->GetMean ()) ;
710  gaussian->SetParameter (2,hpy->GetRMS ()) ;
711  hpy->Fit ("gaussian","RQL") ;
712 
713  hpy->GetXaxis ()->SetRangeUser ( hpy->GetMean () - width * hpy->GetRMS (),
714  hpy->GetMean () + width * hpy->GetRMS ()) ;
715  prof->Fill (strip->GetXaxis ()->GetBinCenter (biny),
716  gaussian->GetParameter (1)) ;
717  prof->SetBinError (biny,gaussian->GetParameter (2)) ;
718 
719  delete gaussian ;
720  delete hpy ;
721  } // loop over the bins
722 
723  delete stripProfile ;
724  return prof ;
725 }
726 
727 
728 // -----------------------------------------------------------------
729 
730 
731 TH1D * smartError (TH1D * strip)
732 {
733 
734  double xmin = strip->GetXaxis ()->GetXmin () ;
735  double xmax = strip->GetXaxis ()->GetXmax () ;
736  int stripsBins = strip->GetNbinsX () ;
737 
738  std::string name = strip->GetName () ;
739  name += "_error" ;
740  TH1D * error = new TH1D
741  (name.c_str (),strip->GetTitle (),stripsBins,xmin,xmax) ;
742 
743  int binmin = 1 ;
744  int ngroup = 1 ; // bins per step
745  int binmax = stripsBins ;
746  for (int bin=binmin ; bin<=binmax ; bin += ngroup)
747  {
748  double dummyError = strip->GetBinError (bin) ;
749  error->SetBinContent (bin,dummyError) ;
750  }
751  return error;
752 }
753 
754 
755 // -----------------------------------------------------------------
756 
757 
758 double effectiveSigma (TH1F & histogram, int vSteps)
759 {
760  double totInt = histogram.Integral () ;
761  int maxBin = histogram.GetMaximumBin () ;
762  int maxBinVal = int(histogram.GetBinContent (maxBin)) ;
763  int totBins = histogram.GetNbinsX () ;
764  double area = totInt ;
765  double threshold = 0 ;
766  double vStep = maxBinVal / vSteps ;
767  int leftBin = 1 ;
768  int rightBin = totBins - 1 ;
769  //loop over the vertical range
770  while (area/totInt > 0.683)
771  {
772  threshold += vStep ;
773  // loop toward the left
774  for (int back = maxBin ; back > 0 ; --back)
775  {
776  if (histogram.GetBinContent (back) < threshold)
777  {
778  leftBin = back ;
779  break ;
780  }
781  } // loop toward the left
782 
783  // loop toward the right
784  for (int fwd = maxBin ; fwd < totBins ; ++fwd)
785  {
786  if (histogram.GetBinContent (fwd) < threshold)
787  {
788  rightBin = fwd ;
789  break ;
790  }
791  } // loop toward the right
792  area = histogram.Integral (leftBin,rightBin) ;
793  } //loop over the vertical range
794 
795  histogram.GetXaxis ()->SetRange (leftBin,rightBin) ;
796  // double sigmaEff = histogram.GetRMS () ;
797  double halfWidthRange = 0.5 * (histogram.GetBinCenter (rightBin) - histogram.GetBinCenter (leftBin)) ;
798  return halfWidthRange ;
799 }
800 
801 
802 // -----------------------------------------------------------------
803 
804 
805 std::pair<int,int> findSupport (TH1F & histogram, double thres)
806 {
807  int totBins = histogram.GetNbinsX () ;
808  if (thres >= histogram.GetMaximum ())
809  return std::pair<int,int> (0, totBins) ;
810 
811  int leftBin = totBins - 1 ;
812  // search from left for the minimum
813  for (int bin=1 ; bin<totBins ; ++bin)
814  {
815  if (histogram.GetBinContent (bin) > thres)
816  {
817  leftBin = bin ;
818  break ;
819  }
820  } // search from left for the minimum
821  int rightBin = 1 ;
822  // search from right for the maximum
823  for (int bin=totBins - 1 ; bin> 0 ; --bin)
824  {
825  if (histogram.GetBinContent (bin) > thres)
826  {
827  rightBin = bin ;
828  break ;
829  }
830  } // search from right for the maximum
831  return std::pair<int,int> (leftBin,rightBin) ;
832 }
833 
834 
835 // -----------------------------------------------------------------
836 
837 
838 void
840  CLHEP::HepMatrix * input,
841  double Default)
842 {
843  for (int eta = 0 ; eta < SCMaxEta ; ++eta)
844  for (int phi = 0 ; phi < SCMaxPhi ; ++phi)
845  {
846  if ((*input)[eta][phi])
847  output[eta][phi] = (*input)[eta][phi] ;
848  else output[eta][phi] = Default ;
849  }
850  return ;
851 }
852 
853 // -----------------------------------------------------------------
854 
855 double etaCorrE1E25 (int eta)
856 {
857  double p0 = 0.807883 ;
858  double p1 = 0.000182551 ;
859  double p2 = -5.76961e-06 ;
860  double p3 = 7.41903e-08 ;
861  double p4 = -2.25384e-10 ;
862 
863  double corr ;
864  if (eta < 6) corr = p0 ;
865  else corr = p0 + p1*eta + p2*eta*eta + p3*eta*eta*eta + p4*eta*eta*eta*eta;
866  return corr/p0 ;
867 }
868 // -----------------------------------------------------------------
869 
870 double etaCorrE1E49 (int eta)
871 {
872  double p0 = 0.799895 ;
873  double p1 = 0.000235487 ;
874  double p2 = -8.26496e-06 ;
875  double p3 = 1.21564e-07 ;
876  double p4 = -4.83286e-10 ;
877 
878  double corr ;
879  if (eta < 8) corr = p0 ;
880  else corr = p0 + p1*eta + p2*eta*eta + p3*eta*eta*eta + p4*eta*eta*eta*eta;
881  return corr/p0 ;
882 }
883 // -----------------------------------------------------------------
884 
885 double etaCorrE1E9 (int eta)
886 {
887  if (eta < 4) return 1.0 ;
888  // grazie Paolo
889  double p0 = 0.834629 ;
890  double p1 = 0.00015254 ;
891  double p2 = -4.91784e-06 ;
892  double p3 = 6.54652e-08 ;
893  double p4 = -2.4894e-10 ;
894 
895  double corr ;
896  if (eta < 6) corr = p0 ;
897  else corr = p0 + p1*eta + p2*eta*eta + p3*eta*eta*eta + p4*eta*eta*eta*eta;
898  return corr/p0 ;
899 }
900 
type
Definition: HCALResponse.h:21
TH1D * smartError(TH1D *strip)
#define SCMaxEta
#define Default
Definition: vmac.h:108
int xtalFromiEtaiPhi(const int &iEta, const int &iPhi)
int readCMSSWcoeff(CLHEP::HepMatrix &calibcoeff, const std::string &inputFileName, double defaultVal=1.)
#define SCMaxPhi
double etaCorrE1E49(int eta)
assert(m_qm.get())
void mtrTransfer(double output[85][20], CLHEP::HepMatrix *input, double Default)
int readCMSSWcoeffForComparison(CLHEP::HepMatrix &calibcoeff, const std::string &inputFileName)
double effectiveSigma(TH1F &histogram, int vSteps=100)
T eta() const
void setStyle()
static std::string const input
Definition: EdmProvDump.cc:43
int iphiFromXtal(const int &xtal)
return((rh^lh)&mask)
int writeCMSSWCoeff(const CLHEP::HepMatrix &amplMatrix, double calibThres, float ERef, const CLHEP::HepMatrix &sigmaMatrix, const CLHEP::HepMatrix &statisticMatrix, std::string fileName="calibOutput.txt", std::string genTag="CAL_GENTAG", std::string method="CAL_METHOD", std::string version="CAL_VERSION", std::string type="CAL_TYPE")
TH1D * smartProfile(TH2F *strip, double width)
int phiFromXtal(const int &xtal)
TH1D * smartGausProfile(TH2F *strip, double width)
int writeCalibTxt(const CLHEP::HepMatrix &AmplitudeMatrix, const CLHEP::HepMatrix &SigmaMatrix, const CLHEP::HepMatrix &StatisticMatrix, std::string fileName="calibOutput.txt")
save (read) CLHEP::HepMatrix to (from) text files
Definition: matrixSaver.h:22
double p4[4]
Definition: TauolaWrapper.h:92
TFile * getGlobalTFile(std::string name="Inv MatrixTFile.root")
int touch(std::string inputFileName)
Definition: matrixSaver.cc:138
TCanvas * getGlobalCanvas(std::string name="Inv MatrixCanvas")
double p2[4]
Definition: TauolaWrapper.h:90
int translateCoeff(const CLHEP::HepMatrix &calibcoeff, const CLHEP::HepMatrix &sigmaMatrix, const CLHEP::HepMatrix &statisticMatrix, std::string SMnumber="1", double calibThres=0.01, std::string fileName="calibOutput.txt", std::string genTag="CAL_GENTAG", std::string method="CAL_METHOD", std::string version="CAL_VERSION", std::string type="CAL_TYPE")
JetCorrectorParameters corr
Definition: classes.h:5
int saveGlobalTFile(std::string name="Inv MatrixFile.root")
CLHEP::HepGenMatrix * getMatrix(std::string inputFileName)
Definition: matrixSaver.cc:147
int extract(std::vector< int > *output, const std::string &dati)
double get3x3(const Float_t energy[7][7])
int ietaFromXtal(const int &xtal)
double get5x5(const Float_t energy[7][7])
double p1[4]
Definition: TauolaWrapper.h:89
std::pair< int, int > findSupport(TH1F &histogram, double thres=0.)
double etaCorrE1E9(int eta)
tuple cout
Definition: gather_cfg.py:121
CLHEP::HepMatrix * getSavedMatrix(const std::string &name)
Definition: DDAxes.h:10
int weight
Definition: histoStyle.py:50
int xtalFromEtaPhi(const int &myEta, const int &myPhi)
double etaCorrE1E25(int eta)
HepGeom::Point3D< Float_t > TBposition(const Float_t amplit[7][7], const Float_t beamEne, const Float_t w0=4.0, const Float_t x0=8.9, const Float_t a0=6.2, const Float_t sideX=24.06, const Float_t sideY=22.02)
int etaFromXtal(const int &xtal)
double p3[4]
Definition: TauolaWrapper.h:91
tuple log
Definition: cmsBatch.py:347
Definition: DDAxes.h:10