00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070 #include <iostream>
00071
00072 #include "TSVDUnfold.h"
00073 #include "TH1D.h"
00074 #include "TH2D.h"
00075 #include "TDecompSVD.h"
00076 #include "TRandom3.h"
00077 #include "TMath.h"
00078
00079 ClassImp(TSVDUnfold)
00080
00081 using namespace std;
00082
00083
00084 TSVDUnfold::TSVDUnfold( const TH1D *bdat, const TH1D *bini, const TH1D *xini, const TH2D *Adet )
00085 : TObject (),
00086 fNdim (0),
00087 fDdim (2),
00088 fNormalize (kFALSE),
00089 fKReg (-1),
00090 fDHist (NULL),
00091 fSVHist (NULL),
00092 fBdat (bdat),
00093 fBini (bini),
00094 fXini (xini),
00095 fAdet (Adet),
00096 fToyhisto (NULL),
00097 fToymat (NULL),
00098 fToyMode (kFALSE),
00099 fMatToyMode (kFALSE)
00100 {
00101
00102
00103
00104
00105 if (bdat->GetNbinsX() != bini->GetNbinsX() ||
00106 bdat->GetNbinsX() != xini->GetNbinsX() ||
00107 bdat->GetNbinsX() != Adet->GetNbinsX() ||
00108 bdat->GetNbinsX() != Adet->GetNbinsY()) {
00109 TString msg = "All histograms must have equal dimension.\n";
00110 msg += Form( " Found: dim(bdat)=%i\n", bdat->GetNbinsX() );
00111 msg += Form( " Found: dim(bini)=%i\n", bini->GetNbinsX() );
00112 msg += Form( " Found: dim(xini)=%i\n", xini->GetNbinsX() );
00113 msg += Form( " Found: dim(Adet)=%i,%i\n", Adet->GetNbinsX(), Adet->GetNbinsY() );
00114 msg += "Please start again!";
00115
00116 Fatal( "Init", msg, "%s" );
00117 }
00118
00119
00120 fNdim = bdat->GetNbinsX();
00121 fDdim = 2;
00122 }
00123
00124
00125 TSVDUnfold::TSVDUnfold( const TSVDUnfold& other )
00126 : TObject ( other ),
00127 fNdim (other.fNdim),
00128 fDdim (other.fDdim),
00129 fNormalize (other.fNormalize),
00130 fKReg (other.fKReg),
00131 fDHist (other.fDHist),
00132 fSVHist (other.fSVHist),
00133 fBdat (other.fBdat),
00134 fBini (other.fBini),
00135 fXini (other.fXini),
00136 fAdet (other.fAdet),
00137 fToyhisto (other.fToyhisto),
00138 fToymat (other.fToymat),
00139 fToyMode (other.fToyMode),
00140 fMatToyMode (other.fMatToyMode)
00141 {
00142
00143 }
00144
00145
00146 TSVDUnfold::~TSVDUnfold()
00147 {
00148
00149 }
00150
00151
00152 TH1D* TSVDUnfold::Unfold( Int_t kreg )
00153 {
00154
00155 fKReg = kreg;
00156
00157
00158 if (!fToyMode && !fMatToyMode) InitHistos( );
00159
00160
00161 TVectorD vb(fNdim), vbini(fNdim), vxini(fNdim), vberr(fNdim);
00162 TMatrixD mA(fNdim, fNdim), mCurv(fNdim, fNdim), mC(fNdim, fNdim);
00163
00164 Double_t eps = 1e-12;
00165 Double_t sreg;
00166
00167
00168 if (fToyMode) { H2V( fToyhisto, vb ); H2Verr( fToyhisto, vberr ); }
00169 else { H2V( fBdat, vb ); H2Verr( fBdat, vberr ); }
00170
00171 H2V( fBini, vbini );
00172 H2V( fXini, vxini );
00173 if (fMatToyMode) H2M( fToymat, mA );
00174 else H2M( fAdet, mA );
00175
00176
00177 Double_t scale = vb.Sum()/vbini.Sum();
00178 vbini *= scale;
00179 vxini *= scale;
00180 mA *= scale;
00181
00182
00183 FillCurvatureMatrix( mCurv, mC );
00184
00185
00186 TDecompSVD CSVD(mC);
00187 TMatrixD CUort = CSVD.GetU();
00188 TMatrixD CVort = CSVD.GetV();
00189 TVectorD CSV = CSVD.GetSig();
00190
00191 TMatrixD CSVM(fNdim, fNdim);
00192 for (Int_t i=0; i<fNdim; i++) CSVM(i,i) = 1/CSV(i);
00193
00194 CUort.Transpose( CUort );
00195 TMatrixD mCinv = (CVort*CSVM)*CUort;
00196
00197
00198 vbini = VecDiv ( vbini, vberr );
00199 vb = VecDiv ( vb, vberr, 1 );
00200 mA = MatDivVec( mA, vberr, 1 );
00201 vberr = VecDiv ( vberr, vberr, 1 );
00202
00203
00204 TDecompSVD ASVD( mA*mCinv );
00205 TMatrixD Uort = ASVD.GetU();
00206 TMatrixD Vort = ASVD.GetV();
00207 TVectorD ASV = ASVD.GetSig();
00208
00209 if (!fToyMode && !fMatToyMode) {
00210 V2H(ASV, *fSVHist);
00211 }
00212
00213 TMatrixD Vreg = mCinv*Vort;
00214
00215 Uort.Transpose(Uort);
00216 TVectorD vdini = Uort*vbini;
00217 TVectorD vd = Uort*vb;
00218
00219 if (!fToyMode && !fMatToyMode) {
00220 V2H(vd, *fDHist);
00221 }
00222
00223
00224 Int_t k = GetKReg()-1;
00225
00226 TVectorD vx(fNdim);
00227
00228
00229 TVectorD vdz(fNdim);
00230 for (Int_t i=0; i<fNdim; i++) {
00231 if (ASV(i)<ASV(0)*eps) sreg = ASV(0)*eps;
00232 else sreg = ASV(i);
00233 vdz(i) = sreg/(sreg*sreg + ASV(k)*ASV(k));
00234 }
00235 TVectorD vz = CompProd( vd, vdz );
00236
00237
00238 TVectorD vw = Vreg*vz;
00239
00240
00241 vx = CompProd( vw, vxini );
00242
00243 if(fNormalize){
00244 scale = vx.Sum();
00245 if (scale > 0) vx *= 1.0/scale;
00246 }
00247
00248
00249 if (!fToyMode && !fMatToyMode) {
00250 Info( "Unfold", "Unfolding param: %i",k+1 );
00251 Info( "Unfold", "Curvature of weight distribution: %f", GetCurvature( vw, mCurv ) );
00252 }
00253
00254 TH1D* h = (TH1D*)fBdat->Clone("unfoldingresult");
00255 for(int i=1; i<=fNdim; i++){
00256 h->SetBinContent(i,0.);
00257 h->SetBinError(i,0.);
00258 }
00259 V2H( vx, *h );
00260
00261 return h;
00262 }
00263
00264
00265 TH2D* TSVDUnfold::GetUnfoldCovMatrix( const TH2D* cov, Int_t ntoys, Int_t seed )
00266 {
00267
00268 fToyMode = true;
00269 TH1D* unfres = 0;
00270 TH2D* unfcov = (TH2D*)fAdet->Clone("unfcovmat");
00271 for(int i=1; i<=fNdim; i++)
00272 for(int j=1; j<=fNdim; j++)
00273 unfcov->SetBinContent(i,j,0.);
00274
00275
00276
00277
00278 TMatrixD L(fNdim,fNdim); L *= 0;
00279
00280 for (Int_t iPar= 0; iPar < fNdim; iPar++) {
00281
00282
00283 L(iPar,iPar) = cov->GetBinContent(iPar+1,iPar+1);
00284 for (Int_t k=0; k<iPar; k++) L(iPar,iPar) -= TMath::Power( L(k,iPar), 2 );
00285 if (L(iPar,iPar) > 0.0) L(iPar,iPar) = TMath::Sqrt(L(iPar,iPar));
00286 else L(iPar,iPar) = 0.0;
00287
00288
00289 for (Int_t jPar=iPar+1; jPar<fNdim; jPar++) {
00290 L(iPar,jPar) = cov->GetBinContent(iPar+1,jPar+1);
00291 for (Int_t k=0; k<iPar; k++) L(iPar,jPar) -= L(k,iPar)*L(k,jPar);
00292 if (L(iPar,iPar)!=0.) L(iPar,jPar) /= L(iPar,iPar);
00293 else L(iPar,jPar) = 0;
00294 }
00295 }
00296
00297
00298 TMatrixD *Lt = new TMatrixD(TMatrixD::kTransposed,L);
00299 TRandom3 random(seed);
00300
00301 fToyhisto = (TH1D*)fBdat->Clone("toyhisto");
00302 TH1D *toymean = (TH1D*)fBdat->Clone("toymean");
00303 for (Int_t j=1; j<=fNdim; j++) toymean->SetBinContent(j,0.);
00304
00305
00306 for (int i=1; i<=ntoys; i++) {
00307
00308
00309 TVectorD g(fNdim);
00310 for (Int_t k= 0; k < fNdim; k++) g(k) = random.Gaus(0.,1.);
00311
00312
00313 g *= (*Lt);
00314
00315
00316 for (int j=1; j<=fNdim; j++) {
00317 fToyhisto->SetBinContent(j,fBdat->GetBinContent(j)+g(j-1));
00318 fToyhisto->SetBinError(j,fBdat->GetBinError(j));
00319 }
00320
00321 unfres = Unfold(GetKReg());
00322
00323 for (Int_t j=1; j<=fNdim; j++) {
00324 toymean->SetBinContent(j, toymean->GetBinContent(j) + unfres->GetBinContent(j)/ntoys);
00325 }
00326 }
00327
00328
00329 random.SetSeed(seed);
00330
00331
00332 for (int i=1; i<=ntoys; i++) {
00333
00334
00335 TVectorD g(fNdim);
00336 for (Int_t k= 0; k < fNdim; k++) g(k) = random.Gaus(0.,1.);
00337
00338
00339 g *= (*Lt);
00340
00341
00342 for (int j=1; j<=fNdim; j++) {
00343 fToyhisto->SetBinContent( j, fBdat->GetBinContent(j)+g(j-1) );
00344 fToyhisto->SetBinError ( j, fBdat->GetBinError(j) );
00345 }
00346 unfres = Unfold(GetKReg());
00347
00348 for (Int_t j=1; j<=fNdim; j++) {
00349 for (Int_t k=1; k<=fNdim; k++) {
00350 unfcov->SetBinContent(j,k,unfcov->GetBinContent(j,k) + ( (unfres->GetBinContent(j) - toymean->GetBinContent(j))* (unfres->GetBinContent(k) - toymean->GetBinContent(k))/(ntoys-1)) );
00351 }
00352 }
00353 }
00354 delete Lt;
00355 delete unfres;
00356 fToyMode = kFALSE;
00357
00358 return unfcov;
00359 }
00360
00361
00362 TH2D* TSVDUnfold::GetAdetCovMatrix( Int_t ntoys, Int_t seed )
00363 {
00364
00365 fMatToyMode = true;
00366 TH1D* unfres = 0;
00367 TH2D* unfcov = (TH2D*)fAdet->Clone("unfcovmat");
00368 for(int i=1; i<=fNdim; i++)
00369 for(int j=1; j<=fNdim; j++)
00370 unfcov->SetBinContent(i,j,0.);
00371
00372
00373 TRandom3 random(seed);
00374
00375 fToymat = (TH2D*)fAdet->Clone("toymat");
00376 TH1D *toymean = (TH1D*)fXini->Clone("toymean");
00377 for (Int_t j=1; j<=fNdim; j++) toymean->SetBinContent(j,0.);
00378
00379 for (int i=1; i<=ntoys; i++) {
00380 for (Int_t k=1; k<=fNdim; k++) {
00381 for (Int_t m=1; m<=fNdim; m++) {
00382 if (fAdet->GetBinContent(k,m)) {
00383 fToymat->SetBinContent(k, m, random.Poisson(fAdet->GetBinContent(k,m)));
00384 }
00385 }
00386 }
00387
00388 unfres = Unfold(GetKReg());
00389
00390 for (Int_t j=1; j<=fNdim; j++) {
00391 toymean->SetBinContent(j, toymean->GetBinContent(j) + unfres->GetBinContent(j)/ntoys);
00392 }
00393 }
00394
00395
00396 random.SetSeed(seed);
00397
00398 for (int i=1; i<=ntoys; i++) {
00399 for (Int_t k=1; k<=fNdim; k++) {
00400 for (Int_t m=1; m<=fNdim; m++) {
00401 if (fAdet->GetBinContent(k,m))
00402 fToymat->SetBinContent(k, m, random.Poisson(fAdet->GetBinContent(k,m)));
00403 }
00404 }
00405
00406 unfres = Unfold(GetKReg());
00407
00408 for (Int_t j=1; j<=fNdim; j++) {
00409 for (Int_t k=1; k<=fNdim; k++) {
00410 unfcov->SetBinContent(j,k,unfcov->GetBinContent(j,k) + ( (unfres->GetBinContent(j) - toymean->GetBinContent(j))*(unfres->GetBinContent(k) - toymean->GetBinContent(k))/(ntoys-1)) );
00411 }
00412 }
00413 }
00414 delete unfres;
00415 fMatToyMode = kFALSE;
00416
00417 return unfcov;
00418 }
00419
00420
00421 TH1D* TSVDUnfold::GetD() const
00422 {
00423
00424 for (int i=1; i<=fDHist->GetNbinsX(); i++) {
00425 if (fDHist->GetBinContent(i)<0.) fDHist->SetBinContent(i, TMath::Abs(fDHist->GetBinContent(i)));
00426 }
00427 return fDHist;
00428 }
00429
00430
00431 TH1D* TSVDUnfold::GetSV() const
00432 {
00433
00434 return fSVHist;
00435 }
00436
00437
00438 void TSVDUnfold::H2V( const TH1D* histo, TVectorD& vec )
00439 {
00440
00441 for (Int_t i=0; i<histo->GetNbinsX(); i++) vec(i) = histo->GetBinContent(i+1);
00442 }
00443
00444
00445 void TSVDUnfold::H2Verr( const TH1D* histo, TVectorD& vec )
00446 {
00447
00448 for (Int_t i=0; i<histo->GetNbinsX(); i++) vec(i) = histo->GetBinError(i+1);
00449 }
00450
00451
00452 void TSVDUnfold::V2H( const TVectorD& vec, TH1D& histo )
00453 {
00454
00455 for(Int_t i=0; i<vec.GetNrows(); i++) histo.SetBinContent(i+1, vec(i));
00456 }
00457
00458
00459 void TSVDUnfold::H2M( const TH2D* histo, TMatrixD& mat )
00460 {
00461
00462 for (Int_t j=0; j<histo->GetNbinsX(); j++) {
00463 for (Int_t i=0; i<histo->GetNbinsY(); i++) {
00464 mat(i,j) = histo->GetBinContent(i+1,j+1);
00465 }
00466 }
00467 }
00468
00469
00470 TVectorD TSVDUnfold::VecDiv( const TVectorD& vec1, const TVectorD& vec2, Int_t zero )
00471 {
00472
00473 TVectorD quot(vec1.GetNrows());
00474 for (Int_t i=0; i<vec1.GetNrows(); i++) {
00475 if (vec2(i) != 0) quot(i) = vec1(i) / vec2(i);
00476 else {
00477 if (zero) quot(i) = 0;
00478 else quot(i) = vec1(i);
00479 }
00480 }
00481 return quot;
00482 }
00483
00484
00485 TMatrixD TSVDUnfold::MatDivVec( const TMatrixD& mat, const TVectorD& vec, Int_t zero )
00486 {
00487
00488 TMatrixD quotmat(mat.GetNrows(), mat.GetNcols());
00489 for (Int_t i=0; i<mat.GetNrows(); i++) {
00490 for (Int_t j=0; j<mat.GetNcols(); j++) {
00491 if (vec(i) != 0) quotmat(i,j) = mat(i,j) / vec(i);
00492 else {
00493 if (zero) quotmat(i,j) = 0;
00494 else quotmat(i,j) = mat(i,j);
00495 }
00496 }
00497 }
00498 return quotmat;
00499 }
00500
00501
00502 TVectorD TSVDUnfold::CompProd( const TVectorD& vec1, const TVectorD& vec2 )
00503 {
00504
00505 TVectorD res(vec1.GetNrows());
00506 for (Int_t i=0; i<vec1.GetNrows(); i++) res(i) = vec1(i) * vec2(i);
00507 return res;
00508 }
00509
00510
00511 Double_t TSVDUnfold::GetCurvature(const TVectorD& vec, const TMatrixD& curv)
00512 {
00513
00514 return vec*(curv*vec);
00515 }
00516
00517
00518 void TSVDUnfold::FillCurvatureMatrix( TMatrixD& tCurv, TMatrixD& tC ) const
00519 {
00520 Double_t eps = 0.00001;
00521
00522 Int_t ndim = tCurv.GetNrows();
00523
00524
00525 tCurv *= 0;
00526 tC *= 0;
00527
00528 if (fDdim == 0) for (Int_t i=0; i<ndim; i++) tC(i,i) = 1;
00529 else if (ndim == 1) {
00530 for (Int_t i=0; i<ndim; i++) {
00531 if (i < ndim-1) tC(i,i+1) = 1.0;
00532 tC(i,i) = 1.0;
00533 }
00534 }
00535 else if (fDdim == 2) {
00536 for (Int_t i=0; i<ndim; i++) {
00537 if (i > 0) tC(i,i-1) = 1.0;
00538 if (i < ndim-1) tC(i,i+1) = 1.0;
00539 tC(i,i) = -2.0;
00540 }
00541 tC(0,0) = -1.0;
00542 tC(ndim-1,ndim-1) = -1.0;
00543 }
00544 else if (fDdim == 3) {
00545 for (Int_t i=1; i<ndim-2; i++) {
00546 tC(i,i-1) = 1.0;
00547 tC(i,i) = -3.0;
00548 tC(i,i+1) = 3.0;
00549 tC(i,i+2) = -1.0;
00550 }
00551 }
00552 else if (fDdim==4) {
00553 for (Int_t i=0; i<ndim; i++) {
00554 if (i > 0) tC(i,i-1) = -4.0;
00555 if (i < ndim-1) tC(i,i+1) = -4.0;
00556 if (i > 1) tC(i,i-2) = 1.0;
00557 if (i < ndim-2) tC(i,i+2) = 1.0;
00558 tC(i,i) = 6.0;
00559 }
00560 tC(0,0) = 2.0;
00561 tC(ndim-1,ndim-1) = 2.0;
00562 tC(0,1) = -3.0;
00563 tC(ndim-2,ndim-1) = -3.0;
00564 tC(1,0) = -3.0;
00565 tC(ndim-1,ndim-2) = -3.0;
00566 tC(1,1) = 6.0;
00567 tC(ndim-2,ndim-2) = 6.0;
00568 }
00569 else if (fDdim == 5) {
00570 for (Int_t i=2; i < ndim-3; i++) {
00571 tC(i,i-2) = 1.0;
00572 tC(i,i-1) = -5.0;
00573 tC(i,i) = 10.0;
00574 tC(i,i+1) = -10.0;
00575 tC(i,i+2) = 5.0;
00576 tC(i,i+3) = -1.0;
00577 }
00578 }
00579 else if (fDdim == 6) {
00580 for (Int_t i = 3; i < ndim - 3; i++) {
00581 tC(i,i-3) = 1.0;
00582 tC(i,i-2) = -6.0;
00583 tC(i,i-1) = 15.0;
00584 tC(i,i) = -20.0;
00585 tC(i,i+1) = 15.0;
00586 tC(i,i+2) = -6.0;
00587 tC(i,i+3) = 1.0;
00588 }
00589 }
00590
00591
00592 for (Int_t i=0; i<ndim; i++) tC(i,i) = tC(i,i) + eps;
00593
00594
00595 for (Int_t i=0; i<ndim; i++) {
00596 for (Int_t j=0; j<ndim; j++) {
00597 for (Int_t k=0; k<ndim; k++) {
00598 tCurv(i,j) = tCurv(i,j) + tC(k,i)*tC(k,j);
00599 }
00600 }
00601 }
00602 }
00603
00604
00605 void TSVDUnfold::InitHistos( )
00606 {
00607
00608 fDHist = new TH1D( "dd", "d vector after orthogonal transformation", fNdim, 0, fNdim );
00609 fDHist->Sumw2();
00610
00611 fSVHist = new TH1D( "sv", "Singular values of AC^-1", fNdim, 0, fNdim );
00612 fSVHist->Sumw2();
00613 }
00614
00615
00616 void TSVDUnfold::RegularisedSymMatInvert( TMatrixDSym& mat, Double_t eps )
00617 {
00618
00619
00620
00621 const UInt_t n = mat.GetNrows();
00622 UInt_t nn = 0;
00623
00624 UInt_t *ipos = new UInt_t[n];
00625
00626
00627
00628 Double_t ymax = 0;
00629 for (UInt_t i=0; i<n; i++) if (TMath::Abs(mat[i][i]) > ymax) ymax = TMath::Abs(mat[i][i]);
00630
00631 for (UInt_t i=0; i<n; i++) {
00632
00633
00634 if (TMath::Abs(mat[i][i])/ymax > eps) ipos[nn++] = i;
00635 }
00636
00637
00638 TMatrixDSym matwork( nn );
00639 for (UInt_t in=0; in<nn; in++) for (UInt_t jn=0; jn<nn; jn++) matwork(in,jn) = 0;
00640
00641
00642 for (UInt_t in=0; in<nn; in++) {
00643
00644 matwork[in][in] = mat[ipos[in]][ipos[in]];
00645 for (UInt_t jn=in+1; jn<nn; jn++) {
00646 matwork[in][jn] = mat[ipos[in]][ipos[jn]];
00647 matwork[jn][in] = matwork[in][jn];
00648 }
00649 }
00650
00651
00652 matwork.Invert();
00653
00654
00655 for (UInt_t i=0; i<n; i++) for (UInt_t j=0; j<n; j++) mat[i][j] = 0;
00656
00657
00658 for (UInt_t in=0; in<nn; in++) {
00659 mat[ipos[in]][ipos[in]] = matwork[in][in];
00660 for (UInt_t jn=in+1; jn<nn; jn++) {
00661 mat[ipos[in]][ipos[jn]] = matwork[in][jn];
00662 mat[ipos[jn]][ipos[in]] = mat[ipos[in]][ipos[jn]];
00663 }
00664 }
00665 delete [] ipos;
00666 }
00667
00668
00669 Double_t TSVDUnfold::ComputeChiSquared( const TH1D& truspec, const TH1D& unfspec, const TH2D& covmat, Double_t regpar )
00670 {
00671
00672 UInt_t n = truspec.GetNbinsX();
00673 TMatrixDSym mat( n );
00674 for (UInt_t i=0; i<n; i++) {
00675 for (UInt_t j=i; j<n; j++) {
00676 mat[i][j] = covmat.GetBinContent( i+1, j+1 );
00677 mat[j][i] = mat[i][j];
00678 }
00679 }
00680
00681 RegularisedSymMatInvert( mat, regpar );
00682
00683
00684 Double_t chi2 = 0;
00685 for (UInt_t i=0; i<n; i++) {
00686 for (UInt_t j=0; j<n; j++) {
00687 chi2 += ( (truspec.GetBinContent( i+1 )-unfspec.GetBinContent( i+1 )) *
00688 (truspec.GetBinContent( j+1 )-unfspec.GetBinContent( j+1 )) * mat[i][j] );
00689 }
00690 }
00691
00692 return chi2;
00693 }
00694