8 #include "CLHEP/Geometry/Point3D.h"
22 gROOT->SetStyle (
"Plain") ;
23 gStyle->SetTextSize(0.5);
25 gStyle->SetOptStat (0) ;
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) ;
45 TCanvas * globalCanvas =
static_cast<TCanvas*
>
46 (gROOT->FindObject (name.c_str ())) ;
49 globalCanvas->Clear () ;
50 globalCanvas->UseCurrentStyle () ;
51 globalCanvas->SetWindowSize (700, 600) ;
56 globalCanvas =
new TCanvas (name.c_str (),name.c_str (), 700, 600) ;
68 TFile * globalTFile = (TFile*) gROOT->FindObject (name.c_str()) ;
72 globalTFile =
new TFile (name.c_str(),
"RECREATE") ;
84 TFile * globalTFile =
static_cast<TFile*
>
85 (gROOT->FindObject (name.c_str ())) ;
86 if (!globalTFile)
return 1 ;
87 globalTFile->Write () ;
88 globalTFile->Close () ;
100 CLHEP::HepMatrix * savedMatrix ;
101 if (reader.
touch (name))
103 savedMatrix =
static_cast<CLHEP::HepMatrix *
> (
130 HepGeom::Point3D<Float_t>
132 const Float_t beamEne,
142 Float_t sumWeight = 0. ;
143 Float_t
depth = x0 * (
log (beamEne)+ a0) ;
144 Float_t sin3 = 0.052335956 ;
146 Float_t invE3x3 = 1. /
get3x3 (amplit) ;
156 caloX += (
eta-3) * sideX * weight;
157 caloY -= (
phi-3) * sideY * weight;
171 HepGeom::Point3D<Float_t>
TBposition (caloX, caloY, 0) ;
188 total += energy[
eta][
phi] ;
205 total += energy[
eta][
phi] ;
245 int xMin = 20 * myEta + 1 ;
246 int xMax = 20 * (myEta + 1) + 1 ;
248 int myCryst = 999999 ;
250 for (
int x = xMin ;
x < xMax ;
x++)
268 return 20 * (iEta-1) + 21 - iPhi ;
278 return int (floor ((xtal-1) / 20) );
288 return (20 - phi - 1) ;
315 std::ifstream _dati (dati.c_str ()) ;
321 do { getline (_dati, dataline,
'\n') ; }
322 while (*dataline.begin () ==
'#') ;
323 std::stringstream linea (dataline) ;
325 while (!linea.eof ())
329 if (buffer != -1) output->push_back (buffer) ;
332 return output->size () ;
342 const CLHEP::HepMatrix & SigmaMatrix,
343 const CLHEP::HepMatrix & StatisticMatrix,
351 if (AmplitudeMatrix[
eta][
phi] &&
352 SigmaMatrix[
eta][
phi] < 100 )
354 reference = AmplitudeMatrix[
eta][
phi] ;
355 std::cout <<
"[InvMatrixUtils][writeCalibTxt] reference crystal: "
356 <<
"(" <<
eta <<
"," <<
phi <<
") -> "
357 << reference <<
"\n" ;
363 std::cerr <<
"ERROR: no calibration coefficients found" << std::endl ;
368 std::ofstream txt_outfile ;
369 txt_outfile.open (fileName.c_str ()) ;
370 txt_outfile <<
"# xtal\tcoeff\tsigma\tevt\tisGood\n" ;
377 if (AmplitudeMatrix[
eta][
phi] == 0) isGood = 0 ;
378 if (SigmaMatrix[
eta][
phi] > 100 ) isGood = 0 ;
380 <<
"\t" << AmplitudeMatrix[
eta][
phi]/reference
381 <<
"\t" << SigmaMatrix[
eta][
phi]
382 <<
"\t" << StatisticMatrix[
eta][
phi]
383 <<
"\t" << isGood <<
"\n" ;
387 txt_outfile.close () ;
398 const CLHEP::HepMatrix & sigmaMatrix,
399 const CLHEP::HepMatrix & statisticMatrix,
407 std::ofstream txt_outfile ;
408 txt_outfile.open (fileName.c_str ()) ;
409 txt_outfile <<
"1\n" ;
410 txt_outfile <<
"-1\n" ;
411 txt_outfile << genTag <<
"\n" ;
412 txt_outfile << method <<
"\n" ;
413 txt_outfile << version <<
"\n" ;
414 txt_outfile << type <<
"\n" ;
422 if (amplMatrix[
eta][
phi] <= calibThres)
427 <<
"\t" << 0 <<
"\n" ;
430 <<
"\t" << reference / amplMatrix[
eta][
phi]
431 <<
"\t" << sigmaMatrix[
eta][
phi]
432 <<
"\t" << statisticMatrix[
eta][
phi]
433 <<
"\t" << 1 <<
"\n" ;
437 txt_outfile.close () ;
447 int etaRef,
int phiRef,
448 const CLHEP::HepMatrix & sigmaMatrix,
449 const CLHEP::HepMatrix & statisticMatrix,
457 std::ofstream txt_outfile ;
458 txt_outfile.open (fileName.c_str ()) ;
459 txt_outfile <<
"1\n" ;
460 txt_outfile <<
"-1\n" ;
461 txt_outfile << genTag <<
"\n" ;
462 txt_outfile << method <<
"\n" ;
463 txt_outfile << version <<
"\n" ;
464 txt_outfile << type <<
"\n" ;
466 if (amplMatrix[etaRef-1][phiRef-1] == 0)
469 << etaRef <<
"," << phiRef
470 <<
") is out of range\n" ;
473 double reference = amplMatrix[etaRef-1][phiRef-1] ;
479 if (amplMatrix[
eta][
phi] <= calibThres)
484 <<
"\t" << 0 <<
"\n" ;
487 <<
"\t" << reference / amplMatrix[
eta][
phi]
488 <<
"\t" << sigmaMatrix[
eta][
phi]
489 <<
"\t" << statisticMatrix[
eta][
phi]
490 <<
"\t" << 1 <<
"\n" ;
494 txt_outfile.close () ;
503 const CLHEP::HepMatrix & sigmaMatrix,
504 const CLHEP::HepMatrix & statisticMatrix,
514 std::ofstream txt_outfile ;
515 txt_outfile.open (fileName.c_str ()) ;
516 txt_outfile << SMnumber <<
"\n" ;
517 txt_outfile <<
"-1\n" ;
518 txt_outfile << genTag <<
"\n" ;
519 txt_outfile << method <<
"\n" ;
520 txt_outfile << version <<
"\n" ;
521 txt_outfile << type <<
"\n" ;
527 if (calibcoeff[
eta][
phi] < calibThres)
533 <<
"\t" << 0 <<
"\n" ;
534 std::cout <<
"[translateCoefff][" << SMnumber
536 <<
" calib coeff below threshold: "
540 <<
"\t" << 0 <<
"\n" ;
544 <<
"\t" << calibcoeff[
eta][
phi]
545 <<
"\t" << sigmaMatrix[
eta][
phi]
546 <<
"\t" << statisticMatrix[
eta][
phi]
547 <<
"\t" << 1 <<
"\n" ;
551 txt_outfile.close () ;
563 std::ifstream CMSSWfile ;
564 CMSSWfile.open (inputFileName.c_str ()) ;
566 CMSSWfile >> buffer ;
567 CMSSWfile >> buffer ;
568 CMSSWfile >> buffer ;
569 CMSSWfile >> buffer ;
570 CMSSWfile >> buffer ;
571 CMSSWfile >> buffer ;
572 while (!CMSSWfile.eof ())
575 CMSSWfile >> xtalnum ;
579 CMSSWfile >> buffer ;
583 if (!good) coeff = defaultVal ;
597 std::ifstream CMSSWfile ;
598 CMSSWfile.open (inputFileName.c_str ()) ;
600 CMSSWfile >> buffer ;
601 CMSSWfile >> buffer ;
602 CMSSWfile >> buffer ;
603 CMSSWfile >> buffer ;
604 CMSSWfile >> buffer ;
605 CMSSWfile >> buffer ;
606 while (!CMSSWfile.eof ())
609 CMSSWfile >> xtalnum ;
613 CMSSWfile >> buffer ;
617 if (!good) coeff = 0. ;
630 TProfile * stripProfile = strip->ProfileX () ;
634 double xmin = stripProfile->GetXaxis ()->GetXmin () ;
635 double xmax = stripProfile->GetXaxis ()->GetXmax () ;
636 int profileBins = stripProfile->GetNbinsX () ;
640 TH1D * prof =
new TH1D
641 (name.c_str (),strip->GetTitle (),profileBins,
xmin,
xmax) ;
644 int nbins = strip->GetXaxis ()->GetNbins () ;
650 for (
int bin=binmin ;
bin<=binmax ;
bin += ngroup)
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 ;}
657 Int_t biny =
bin + ngroup/2 ;
659 hpy->GetXaxis ()->SetRangeUser ( hpy->GetMean () - width * hpy->GetRMS (),
660 hpy->GetMean () + width * hpy->GetRMS ()) ;
661 prof->Fill (strip->GetXaxis ()->GetBinCenter (biny),
663 prof->SetBinError (biny,hpy->GetRMS()) ;
668 delete stripProfile ;
678 TProfile * stripProfile = strip->ProfileX () ;
682 double xmin = stripProfile->GetXaxis ()->GetXmin () ;
683 double xmax = stripProfile->GetXaxis ()->GetXmax () ;
684 int profileBins = stripProfile->GetNbinsX () ;
687 name +=
"_smartGaus" ;
688 TH1D * prof =
new TH1D
689 (name.c_str (),strip->GetTitle (),profileBins,
xmin,
xmax) ;
692 int nbins = strip->GetXaxis ()->GetNbins () ;
698 for (
int bin=binmin ;
bin<=binmax ;
bin += ngroup)
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 ;}
705 Int_t biny =
bin + ngroup/2 ;
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") ;
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)) ;
723 delete stripProfile ;
734 double xmin = strip->GetXaxis ()->GetXmin () ;
735 double xmax = strip->GetXaxis ()->GetXmax () ;
736 int stripsBins = strip->GetNbinsX () ;
740 TH1D *
error =
new TH1D
741 (name.c_str (),strip->GetTitle (),stripsBins,
xmin,
xmax) ;
745 int binmax = stripsBins ;
746 for (
int bin=binmin ;
bin<=binmax ;
bin += ngroup)
748 double dummyError = strip->GetBinError (
bin) ;
749 error->SetBinContent (
bin,dummyError) ;
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 ;
766 double vStep = maxBinVal / vSteps ;
768 int rightBin = totBins - 1 ;
770 while (area/totInt > 0.683)
774 for (
int back = maxBin ; back > 0 ; --back)
776 if (histogram.GetBinContent (back) <
threshold)
784 for (
int fwd = maxBin ; fwd < totBins ; ++fwd)
786 if (histogram.GetBinContent (fwd) <
threshold)
792 area = histogram.Integral (leftBin,rightBin) ;
795 histogram.GetXaxis ()->SetRange (leftBin,rightBin) ;
797 double halfWidthRange = 0.5 * (histogram.GetBinCenter (rightBin) - histogram.GetBinCenter (leftBin)) ;
798 return halfWidthRange ;
807 int totBins = histogram.GetNbinsX () ;
808 if (thres >= histogram.GetMaximum ())
809 return std::pair<int,int> (0, totBins) ;
811 int leftBin = totBins - 1 ;
815 if (histogram.GetBinContent (
bin) > thres)
823 for (
int bin=totBins - 1 ;
bin> 0 ; --
bin)
825 if (histogram.GetBinContent (
bin) > thres)
831 return std::pair<int,int> (leftBin,rightBin) ;
840 CLHEP::HepMatrix *
input,
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 ;
864 if (eta < 6) corr = p0 ;
865 else corr = p0 + p1*eta + p2*eta*eta + p3*eta*eta*eta + p4*eta*eta*eta*
eta;
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 ;
879 if (eta < 8) corr = p0 ;
880 else corr = p0 + p1*eta + p2*eta*eta + p3*eta*eta*eta + p4*eta*eta*eta*
eta;
887 if (eta < 4)
return 1.0 ;
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 ;
896 if (eta < 6) corr = p0 ;
897 else corr = p0 + p1*eta + p2*eta*eta + p3*eta*eta*eta + p4*eta*eta*eta*
eta;
TH1D * smartError(TH1D *strip)
int xtalFromiEtaiPhi(const int &iEta, const int &iPhi)
int readCMSSWcoeff(CLHEP::HepMatrix &calibcoeff, const std::string &inputFileName, double defaultVal=1.)
double etaCorrE1E49(int eta)
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)
static std::string const input
int iphiFromXtal(const int &xtal)
int writeCMSSWCoeff(const CLHEP::HepMatrix &lMatrix, 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)
T x() const
Cartesian x coordinate.
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
TFile * getGlobalTFile(std::string name="Inv MatrixTFile.root")
bool touch(std::string inputFileName)
TCanvas * getGlobalCanvas(std::string name="Inv MatrixCanvas")
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")
int saveGlobalTFile(std::string name="Inv MatrixFile.root")
CLHEP::HepGenMatrix * getMatrix(std::string inputFileName)
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])
Geom::Phi< T > phi() const
return(e1-e2)*(e1-e2)+dp *dp
std::pair< int, int > findSupport(TH1F &histogram, double thres=0.)
double etaCorrE1E9(int eta)
CLHEP::HepMatrix * getSavedMatrix(const std::string &name)
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)