79 TSVDUnfold::TSVDUnfold(
const TH1D *bdat,
const TH1D *bini,
const TH1D *xini,
const TH2D *Adet )
98 if (bdat->GetNbinsX() != bini->GetNbinsX() ||
99 bdat->GetNbinsX() != xini->GetNbinsX() ||
100 bdat->GetNbinsX() != Adet->GetNbinsX() ||
101 bdat->GetNbinsX() != Adet->GetNbinsY()) {
102 TString msg =
"All histograms must have equal dimension.\n";
103 msg += Form(
" Found: dim(bdat)=%i\n", bdat->GetNbinsX() );
104 msg += Form(
" Found: dim(bini)=%i\n", bini->GetNbinsX() );
105 msg += Form(
" Found: dim(xini)=%i\n", xini->GetNbinsX() );
106 msg += Form(
" Found: dim(Adet)=%i,%i\n", Adet->GetNbinsX(), Adet->GetNbinsY() );
107 msg +=
"Please start again!";
109 Fatal(
"Init", msg,
"%s" );
112 fBcov = (TH2D*)fAdet->Clone(
"bcov");
114 for(
int i=1; i<=fBdat->GetNbinsX(); i++){
115 fBcov->SetBinContent(i, i, fBdat->GetBinError(i)*fBdat->GetBinError(i));
116 for(
int j=1; j<=fBdat->GetNbinsX(); j++){
118 fBcov->SetBinContent(i,j,0.);
122 fNdim = bdat->GetNbinsX();
132 TSVDUnfold::TSVDUnfold(
const TH1D *bdat, TH2D* Bcov,
const TH1D *bini,
const TH1D *xini,
const TH2D *Adet )
152 if (bdat->GetNbinsX() != bini->GetNbinsX() ||
153 bdat->GetNbinsX() != xini->GetNbinsX() ||
154 bdat->GetNbinsX() != Bcov->GetNbinsX() ||
155 bdat->GetNbinsX() != Bcov->GetNbinsY() ||
156 bdat->GetNbinsX() != Adet->GetNbinsX() ||
157 bdat->GetNbinsX() != Adet->GetNbinsY()) {
158 TString msg =
"All histograms must have equal dimension.\n";
159 msg += Form(
" Found: dim(bdat)=%i\n", bdat->GetNbinsX() );
160 msg += Form(
" Found: dim(Bcov)=%i,%i\n", Bcov->GetNbinsX(), Bcov->GetNbinsY() );
161 msg += Form(
" Found: dim(bini)=%i\n", bini->GetNbinsX() );
162 msg += Form(
" Found: dim(xini)=%i\n", xini->GetNbinsX() );
163 msg += Form(
" Found: dim(Adet)=%i,%i\n", Adet->GetNbinsX(), Adet->GetNbinsY() );
164 msg +=
"Please start again!";
166 Fatal(
"Init", msg,
"%s" );
170 fNdim = bdat->GetNbinsX();
177 TSVDUnfold::TSVDUnfold(
const TSVDUnfold& other )
181 fNormalize (other.fNormalize),
183 fDHist (other.fDHist),
184 fSVHist (other.fSVHist),
192 fToyhisto (other.fToyhisto),
193 fToymat (other.fToymat),
194 fToyMode (other.fToyMode),
195 fMatToyMode (other.fMatToyMode)
202 TSVDUnfold::~TSVDUnfold()
243 TH1D* TSVDUnfold::Unfold( Int_t kreg )
248 if (!fToyMode && !fMatToyMode) InitHistos( );
251 TVectorD vb(fNdim), vbini(fNdim), vxini(fNdim), vberr(fNdim);
252 TMatrixD mB(fNdim, fNdim), mA(fNdim, fNdim), mCurv(fNdim, fNdim), mC(fNdim, fNdim);
254 Double_t eps = 1e-12;
258 if (fToyMode) { H2V( fToyhisto, vb ); H2Verr( fToyhisto, vberr ); }
259 else { H2V( fBdat, vb ); H2Verr( fBdat, vberr ); }
264 if (fMatToyMode) H2M( fToymat, mA );
265 else H2M( fAdet, mA );
268 FillCurvatureMatrix( mCurv, mC );
272 TMatrixD CUort = CSVD.GetU();
273 TMatrixD CVort = CSVD.GetV();
274 TVectorD CSV = CSVD.GetSig();
276 TMatrixD CSVM(fNdim, fNdim);
277 for (Int_t i=0; i<fNdim; i++) CSVM(i,i) = 1/CSV(i);
279 CUort.Transpose( CUort );
280 TMatrixD mCinv = (CVort*CSVM)*CUort;
283 TDecompSVD BSVD( mB );
284 TMatrixD QT = BSVD.GetU();
286 TVectorD B2SV = BSVD.GetSig();
289 for(
int i=0; i<fNdim; i++){
290 BSV(i) = TMath::Sqrt(B2SV(i));
292 TMatrixD mAtmp(fNdim,fNdim);
293 TVectorD vbtmp(fNdim);
296 for(
int i=0; i<fNdim; i++){
297 for(
int j=0; j<fNdim; j++){
299 vbtmp(i) += QT(i,j)*vb(j)/BSV(i);
301 for(
int m=0; m<fNdim; m++){
303 mAtmp(i,j) += QT(i,m)*mA(m,j)/BSV(i);
312 TDecompSVD ASVD( mA*mCinv );
313 TMatrixD Uort = ASVD.GetU();
314 TMatrixD Vort = ASVD.GetV();
315 TVectorD ASV = ASVD.GetSig();
317 if (!fToyMode && !fMatToyMode) {
321 TMatrixD Vreg = mCinv*Vort;
323 Uort.Transpose(Uort);
324 TVectorD vd = Uort*vb;
326 if (!fToyMode && !fMatToyMode) {
331 Int_t k = GetKReg()-1;
337 TMatrixD Z(fNdim, fNdim);
338 for (Int_t i=0; i<fNdim; i++) {
339 if (ASV(i)<ASV(0)*eps) sreg = ASV(0)*eps;
341 vdz(i) = sreg/(sreg*sreg + ASV(k)*ASV(k));
342 Z(i,i) = vdz(i)*vdz(i);
344 TVectorD vz = CompProd( vd, vdz );
346 TMatrixD VortT(Vort);
347 VortT.Transpose(VortT);
348 TMatrixD W = mCinv*Vort*Z*VortT*mCinv;
350 TMatrixD Xtau(fNdim, fNdim);
351 TMatrixD Xinv(fNdim, fNdim);
354 for (Int_t i=0; i<fNdim; i++) {
355 for (Int_t j=0; j<fNdim; j++) {
356 Xtau(i,j) = vxini(i) * vxini(j) * W(i,j);
359 for (Int_t m=0; m<fNdim; m++) {
360 a += mA(m,i)*mA(m,j);
362 if(vxini(i) && vxini(j))
363 Xinv(i,j) = a/vxini(i)/vxini(j);
368 TVectorD vw = Vreg*vz;
371 vx = CompProd( vw, vxini );
374 Double_t scale = vx.Sum();
377 Xtau *= 1./scale/scale;
382 if (!fToyMode && !fMatToyMode) {
388 if (!fToyMode && !fMatToyMode) {
389 Info(
"Unfold",
"Unfolding param: %i",k+1 );
390 Info(
"Unfold",
"Curvature of weight distribution: %f", GetCurvature( vw, mCurv ) );
393 TH1D* h = (TH1D*)fBdat->Clone(
"unfoldingresult");
394 for(
int i=1; i<=fNdim; i++){
395 h->SetBinContent(i,0.);
396 h->SetBinError(i,0.);
411 TH2D* TSVDUnfold::GetUnfoldCovMatrix(
const TH2D* cov, Int_t ntoys, Int_t seed )
415 TH2D* unfcov = (TH2D*)fAdet->Clone(
"unfcovmat");
416 unfcov->SetTitle(
"Toy covariance matrix");
417 for(
int i=1; i<=fNdim; i++)
418 for(
int j=1; j<=fNdim; j++)
419 unfcov->SetBinContent(i,j,0.);
424 TMatrixD L(fNdim,fNdim); L *= 0;
426 for (Int_t iPar= 0; iPar < fNdim; iPar++) {
429 L(iPar,iPar) = cov->GetBinContent(iPar+1,iPar+1);
430 for (Int_t k=0; k<iPar; k++) L(iPar,iPar) -= TMath::Power( L(k,iPar), 2 );
431 if (L(iPar,iPar) > 0.0) L(iPar,iPar) = TMath::Sqrt(L(iPar,iPar));
432 else L(iPar,iPar) = 0.0;
435 for (Int_t jPar=iPar+1; jPar<fNdim; jPar++) {
436 L(iPar,jPar) = cov->GetBinContent(iPar+1,jPar+1);
437 for (Int_t k=0; k<iPar; k++) L(iPar,jPar) -= L(k,iPar)*L(k,jPar);
438 if (L(iPar,iPar)!=0.) L(iPar,jPar) /= L(iPar,iPar);
439 else L(iPar,jPar) = 0;
444 TMatrixD *Lt =
new TMatrixD(TMatrixD::kTransposed,L);
445 TRandom3 random(seed);
447 fToyhisto = (TH1D*)fBdat->Clone(
"toyhisto");
448 TH1D *toymean = (TH1D*)fBdat->Clone(
"toymean");
449 for (Int_t j=1; j<=fNdim; j++) toymean->SetBinContent(j,0.);
452 for (
int i=1; i<=ntoys; i++) {
456 for (Int_t k= 0; k < fNdim; k++) g(k) = random.Gaus(0.,1.);
462 for (
int j=1; j<=fNdim; j++) {
463 fToyhisto->SetBinContent(j,fBdat->GetBinContent(j)+g(j-1));
464 fToyhisto->SetBinError(j,fBdat->GetBinError(j));
467 unfres = Unfold(GetKReg());
469 for (Int_t j=1; j<=fNdim; j++) {
470 toymean->SetBinContent(j, toymean->GetBinContent(j) + unfres->GetBinContent(j)/ntoys);
477 random.SetSeed(seed);
480 for (
int i=1; i<=ntoys; i++) {
484 for (Int_t k= 0; k < fNdim; k++) g(k) = random.Gaus(0.,1.);
490 for (
int j=1; j<=fNdim; j++) {
491 fToyhisto->SetBinContent( j, fBdat->GetBinContent(j)+g(j-1) );
492 fToyhisto->SetBinError ( j, fBdat->GetBinError(j) );
494 unfres = Unfold(GetKReg());
496 for (Int_t j=1; j<=fNdim; j++) {
497 for (Int_t k=1; k<=fNdim; k++) {
498 unfcov->SetBinContent(j,k,unfcov->GetBinContent(j,k) + ( (unfres->GetBinContent(j) - toymean->GetBinContent(j))* (unfres->GetBinContent(k) - toymean->GetBinContent(k))/(ntoys-1)) );
517 TH2D* TSVDUnfold::GetAdetCovMatrix( Int_t ntoys, Int_t seed )
521 TH2D* unfcov = (TH2D*)fAdet->Clone(
"unfcovmat");
522 unfcov->SetTitle(
"Toy covariance matrix");
523 for(
int i=1; i<=fNdim; i++)
524 for(
int j=1; j<=fNdim; j++)
525 unfcov->SetBinContent(i,j,0.);
528 TRandom3 random(seed);
530 fToymat = (TH2D*)fAdet->Clone(
"toymat");
531 TH1D *toymean = (TH1D*)fXini->Clone(
"toymean");
532 for (Int_t j=1; j<=fNdim; j++) toymean->SetBinContent(j,0.);
534 for (
int i=1; i<=ntoys; i++) {
535 for (Int_t k=1; k<=fNdim; k++) {
536 for (Int_t m=1; m<=fNdim; m++) {
537 if (fAdet->GetBinContent(k,m)) {
538 fToymat->SetBinContent(k, m, random.Poisson(fAdet->GetBinContent(k,m)));
543 unfres = Unfold(GetKReg());
545 for (Int_t j=1; j<=fNdim; j++) {
546 toymean->SetBinContent(j, toymean->GetBinContent(j) + unfres->GetBinContent(j)/ntoys);
553 random.SetSeed(seed);
555 for (
int i=1; i<=ntoys; i++) {
556 for (Int_t k=1; k<=fNdim; k++) {
557 for (Int_t m=1; m<=fNdim; m++) {
558 if (fAdet->GetBinContent(k,m))
559 fToymat->SetBinContent(k, m, random.Poisson(fAdet->GetBinContent(k,m)));
563 unfres = Unfold(GetKReg());
565 for (Int_t j=1; j<=fNdim; j++) {
566 for (Int_t k=1; k<=fNdim; k++) {
567 unfcov->SetBinContent(j,k,unfcov->GetBinContent(j,k) + ( (unfres->GetBinContent(j) - toymean->GetBinContent(j))*(unfres->GetBinContent(k) - toymean->GetBinContent(k))/(ntoys-1)) );
574 fMatToyMode = kFALSE;
582 TH1D* TSVDUnfold::GetD()
const
584 for (
int i=1; i<=fDHist->GetNbinsX(); i++) {
585 if (fDHist->GetBinContent(i)<0.) fDHist->SetBinContent(i, TMath::Abs(fDHist->GetBinContent(i)));
593 TH1D* TSVDUnfold::GetSV()
const
602 TH2D* TSVDUnfold::GetXtau()
const
610 TH2D* TSVDUnfold::GetXinv()
const
618 TH2D* TSVDUnfold::GetBCov()
const
626 void TSVDUnfold::H2V(
const TH1D* histo, TVectorD& vec )
628 for (Int_t i=0; i<histo->GetNbinsX(); i++) vec(i) = histo->GetBinContent(i+1);
634 void TSVDUnfold::H2Verr(
const TH1D* histo, TVectorD& vec )
636 for (Int_t i=0; i<histo->GetNbinsX(); i++) vec(i) = histo->GetBinError(i+1);
642 void TSVDUnfold::V2H(
const TVectorD& vec, TH1D& histo )
644 for(Int_t i=0; i<vec.GetNrows(); i++) histo.SetBinContent(i+1, vec(i));
650 void TSVDUnfold::H2M(
const TH2D* histo, TMatrixD& mat )
652 for (Int_t j=0; j<histo->GetNbinsX(); j++) {
653 for (Int_t i=0; i<histo->GetNbinsY(); i++) {
654 mat(i,j) = histo->GetBinContent(i+1,j+1);
662 void TSVDUnfold::M2H(
const TMatrixD& mat, TH2D& histo )
664 for (Int_t j=0; j<mat.GetNcols(); j++) {
665 for (Int_t i=0; i<mat.GetNrows(); i++) {
666 histo.SetBinContent(i+1,j+1, mat(i,j));
674 TVectorD TSVDUnfold::VecDiv(
const TVectorD& vec1,
const TVectorD& vec2, Int_t zero )
676 TVectorD quot(vec1.GetNrows());
677 for (Int_t i=0; i<vec1.GetNrows(); i++) {
678 if (vec2(i) != 0) quot(i) = vec1(i) / vec2(i);
680 if (zero) quot(i) = 0;
681 else quot(i) = vec1(i);
690 TMatrixD TSVDUnfold::MatDivVec(
const TMatrixD& mat,
const TVectorD& vec, Int_t zero )
692 TMatrixD quotmat(mat.GetNrows(), mat.GetNcols());
693 for (Int_t i=0; i<mat.GetNrows(); i++) {
694 for (Int_t j=0; j<mat.GetNcols(); j++) {
695 if (vec(i) != 0) quotmat(i,j) = mat(i,j) / vec(i);
697 if (zero) quotmat(i,j) = 0;
698 else quotmat(i,j) = mat(i,j);
708 TVectorD TSVDUnfold::CompProd(
const TVectorD& vec1,
const TVectorD& vec2 )
710 TVectorD res(vec1.GetNrows());
711 for (Int_t i=0; i<vec1.GetNrows(); i++) res(i) = vec1(i) * vec2(i);
718 Double_t TSVDUnfold::GetCurvature(
const TVectorD& vec,
const TMatrixD& curv)
720 return vec*(curv*vec);
725 void TSVDUnfold::FillCurvatureMatrix( TMatrixD& tCurv, TMatrixD& tC )
const
727 Double_t eps = 0.00001;
729 Int_t ndim = tCurv.GetNrows();
735 if (fDdim == 0)
for (Int_t i=0; i<ndim; i++) tC(i,i) = 1;
736 else if (fDdim == 1) {
737 for (Int_t i=0; i<ndim; i++) {
738 if (i < ndim-1) tC(i,i+1) = 1.0;
742 else if (fDdim == 2) {
743 for (Int_t i=0; i<ndim; i++) {
744 if (i > 0) tC(i,i-1) = 1.0;
745 if (i < ndim-1) tC(i,i+1) = 1.0;
749 tC(ndim-1,ndim-1) = -1.0;
751 else if (fDdim == 3) {
752 for (Int_t i=1; i<ndim-2; i++) {
760 for (Int_t i=0; i<ndim; i++) {
761 if (i > 0) tC(i,i-1) = -4.0;
762 if (i < ndim-1) tC(i,i+1) = -4.0;
763 if (i > 1) tC(i,i-2) = 1.0;
764 if (i < ndim-2) tC(i,i+2) = 1.0;
768 tC(ndim-1,ndim-1) = 2.0;
770 tC(ndim-2,ndim-1) = -3.0;
772 tC(ndim-1,ndim-2) = -3.0;
774 tC(ndim-2,ndim-2) = 6.0;
776 else if (fDdim == 5) {
777 for (Int_t i=2; i < ndim-3; i++) {
786 else if (fDdim == 6) {
787 for (Int_t i = 3; i < ndim - 3; i++) {
799 for (Int_t i=0; i<ndim; i++) tC(i,i) = tC(i,i) + eps;
802 for (Int_t i=0; i<ndim; i++) {
803 for (Int_t j=0; j<ndim; j++) {
804 for (Int_t k=0; k<ndim; k++) {
805 tCurv(i,j) = tCurv(i,j) + tC(k,i)*tC(k,j);
813 void TSVDUnfold::InitHistos( )
815 fDHist =
new TH1D(
"dd",
"d vector after orthogonal transformation", fNdim, 0, fNdim );
818 fSVHist =
new TH1D(
"sv",
"Singular values of AC^-1", fNdim, 0, fNdim );
821 fXtau = (TH2D*)fAdet->Clone(
"Xtau");
822 fXtau->SetTitle(
"Regularized covariance matrix");
825 fXinv = (TH2D*)fAdet->Clone(
"Xinv");
826 fXinv->SetTitle(
"Inverse covariance matrix");
833 void TSVDUnfold::RegularisedSymMatInvert( TMatrixDSym& mat, Double_t eps )
836 const UInt_t n = mat.GetNrows();
839 UInt_t *ipos =
new UInt_t[n];
844 for (UInt_t i=0; i<n; i++) if (TMath::Abs(mat[i][i]) > ymax) ymax = TMath::Abs(mat[i][i]);
846 for (UInt_t i=0; i<n; i++) {
849 if (TMath::Abs(mat[i][i])/ymax > eps) ipos[nn++] = i;
853 TMatrixDSym matwork( nn );
854 for (UInt_t in=0; in<nn; in++)
for (UInt_t jn=0; jn<nn; jn++) matwork(in,jn) = 0;
857 for (UInt_t in=0; in<nn; in++) {
859 matwork[in][in] = mat[ipos[in]][ipos[in]];
860 for (UInt_t jn=in+1; jn<nn; jn++) {
861 matwork[in][jn] = mat[ipos[in]][ipos[jn]];
862 matwork[jn][in] = matwork[in][jn];
870 for (UInt_t i=0; i<n; i++)
for (UInt_t j=0; j<n; j++) mat[i][j] = 0;
873 for (UInt_t in=0; in<nn; in++) {
874 mat[ipos[in]][ipos[in]] = matwork[in][in];
875 for (UInt_t jn=in+1; jn<nn; jn++) {
876 mat[ipos[in]][ipos[jn]] = matwork[in][jn];
877 mat[ipos[jn]][ipos[in]] = mat[ipos[in]][ipos[jn]];
886 Double_t TSVDUnfold::ComputeChiSquared(
const TH1D& truspec,
const TH1D& unfspec)
888 UInt_t n = truspec.GetNbinsX();
892 for (UInt_t i=0; i<n; i++) {
893 for (UInt_t j=0; j<n; j++) {
894 chi2 += ( (truspec.GetBinContent( i+1 )-unfspec.GetBinContent( i+1 )) *
895 (truspec.GetBinContent( j+1 )-unfspec.GetBinContent( j+1 )) * fXinv->GetBinContent(i+1,j+1) );