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
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118
00119
00120
00121
00122
00123
00124
00125
00126
00127
00128
00129
00130
00131
00132
00133
00134
00135
00136
00137
00138
00139
00140
00141
00142
00143
00144
00145
00146
00147
00148
00149
00150
00151
00152
00153
00154
00155
00156
00157
00158
00159
00160
00161
00162
00163
00164
00165
00166
00167
00168
00169
00170
00171
00172
00173
00174
00175
00176
00177
00178
00179
00180
00181
00182
00183
00184
00185
00186
00187
00188
00189
00190
00191
00192
00193
00194
00195
00196
00197
00198
00199
00200
00201
00202
00203
00204
00205
00206
00207
00208
00209
00210
00211
00212
00213
00214
00215
00216
00217
00218
00219
00220
00221
00222
00223
00224
00225
00226
00227
00228
00229
00230
00231
00232
00233
00234
00235
00236
00237
00238
00239
00240
00241
00242
00243
00244
00245
00246
00247
00248
00249
00250
00251
00252
00253
00254
00255
00256
00257
00258
00259
00260
00261
00262
00263
00264
00265
00266
00267
00268
00269
00270
00271
00272
00273
00274
00275
00276
00277
00278
00279
00280
00281
00282
00283
00284
00285
00286
00287
00288
00289
00290
00291
00292
00293
00294
00295
00296
00297
00298
00299
00300
00301
00302
00303
00304
00305
00306
00307
00308
00309
00310
00311
00312
00313
00314
00315
00316
00317
00318
00319
00320
00321
00322
00323 #include <iostream>
00324 #include <TMatrixD.h>
00325 #include <TMatrixDSparse.h>
00326 #include <TMatrixDSym.h>
00327 #include <TMath.h>
00328 #include <TDecompBK.h>
00329
00330 #include <TUnfold.h>
00331
00332 #include <map>
00333 #include <vector>
00334
00335
00336
00337
00338 #define SCHUR_COMPLEMENT_MATRIX_INVERSION
00339
00340
00341
00342
00343
00344
00345 #define INVERT_WITH_CONDITIONING
00346
00347
00348
00349
00350
00351
00352 #ifdef DEBUG_LCURVE
00353 #include <TCanvas.h>
00354 #endif
00355
00356 ClassImp(TUnfold)
00357
00358
00359 const char *TUnfold::GetTUnfoldVersion(void) {
00360 return TUnfold_VERSION;
00361 }
00362
00363
00364 void TUnfold::InitTUnfold(void)
00365 {
00366
00367 fXToHist.Set(0);
00368 fHistToX.Set(0);
00369 fSumOverY.Set(0);
00370 fA = 0;
00371 fLsquared = 0;
00372 fVyy = 0;
00373 fY = 0;
00374 fX0 = 0;
00375 fTauSquared = 0.0;
00376 fBiasScale = 0.0;
00377 fNdf = 0;
00378 fConstraint = kEConstraintNone;
00379 fRegMode = kRegModeNone;
00380
00381 fVxx = 0;
00382 fX = 0;
00383 fAx = 0;
00384 fChi2A = 0.0;
00385 fLXsquared = 0.0;
00386 fRhoMax = 999.0;
00387 fRhoAvg = -1.0;
00388 fDXDAM[0] = 0;
00389 fDXDAZ[0] = 0;
00390 fDXDAM[1] = 0;
00391 fDXDAZ[1] = 0;
00392 fDXDtauSquared = 0;
00393 fDXDY = 0;
00394 fEinv = 0;
00395 fE = 0;
00396 fVxxInv = 0;
00397 }
00398
00399 void TUnfold::DeleteMatrix(TMatrixD **m) {
00400 if(*m) delete *m;
00401 *m=0;
00402 }
00403
00404 void TUnfold::DeleteMatrix(TMatrixDSparse **m) {
00405 if(*m) delete *m;
00406 *m=0;
00407 }
00408
00409 void TUnfold::ClearResults(void) {
00410
00411
00412
00413 DeleteMatrix(&fVxx);
00414 DeleteMatrix(&fX);
00415 DeleteMatrix(&fAx);
00416 for(Int_t i=0;i<2;i++) {
00417 DeleteMatrix(fDXDAM+i);
00418 DeleteMatrix(fDXDAZ+i);
00419 }
00420 DeleteMatrix(&fDXDtauSquared);
00421 DeleteMatrix(&fDXDY);
00422 DeleteMatrix(&fEinv);
00423 DeleteMatrix(&fE);
00424 DeleteMatrix(&fVxxInv);
00425 fChi2A = 0.0;
00426 fLXsquared = 0.0;
00427 fRhoMax = 999.0;
00428 fRhoAvg = -1.0;
00429 }
00430
00431 TUnfold::TUnfold(void)
00432 {
00433
00434 InitTUnfold();
00435 }
00436
00437 Double_t TUnfold::DoUnfold(void)
00438 {
00439
00440
00441
00442
00443
00444
00445
00446
00447
00448
00449
00450
00451
00452
00453
00454
00455
00456
00457
00458
00459
00460
00461
00462
00463
00464
00465
00466
00467
00468
00469 ClearResults();
00470
00471
00472
00473 Int_t ny=fVyy->GetNrows();
00474 const Int_t *vyy_rows=fVyy->GetRowIndexArray();
00475 const Int_t *vyy_cols=fVyy->GetColIndexArray();
00476 const Double_t *vyy_data=fVyy->GetMatrixArray();
00477 Int_t *usedBin=new Int_t[ny];
00478 for(Int_t i=0;i<ny;i++) {
00479 usedBin[i]=0;
00480 }
00481 for(Int_t i=0;i<ny;i++) {
00482 for(Int_t k=vyy_rows[i];k<vyy_rows[i+1];k++) {
00483 if(vyy_data[k]>0.0) {
00484 usedBin[i]++;
00485 usedBin[vyy_cols[k]]++;
00486 }
00487 }
00488 }
00489 Int_t n=0;
00490 Int_t *yToI=new Int_t[ny];
00491 Int_t *iToY=new Int_t[ny];
00492 for(Int_t i=0;i<ny;i++) {
00493 if(usedBin[i]) {
00494 yToI[i]=n;
00495 iToY[n]=i;
00496 n++;
00497 } else {
00498 yToI[i]=-1;
00499 }
00500 }
00501 delete[] usedBin;
00502
00503
00504
00505
00506
00507
00508 Int_t nn=vyy_rows[ny];
00509 Int_t *vyyred_rows=new Int_t[nn];
00510 Int_t *vyyred_cols=new Int_t[nn];
00511 Double_t *vyyred_data=new Double_t[nn];
00512 nn=0;
00513 for(Int_t i=0;i<ny;i++) {
00514 for(Int_t k=vyy_rows[i];k<vyy_rows[i+1];k++) {
00515 if(vyy_data[k]>0.0) {
00516 vyyred_rows[nn]=yToI[i];
00517 vyyred_cols[nn]=yToI[vyy_cols[k]];
00518 vyyred_data[nn]=vyy_data[k];
00519 nn++;
00520 }
00521 }
00522 }
00523 TMatrixDSparse *vyyred=CreateSparseMatrix
00524 (n,n,nn,vyyred_rows,vyyred_cols,vyyred_data);
00525 delete[] vyyred_rows;
00526 delete[] vyyred_cols;
00527 delete[] vyyred_data;
00528
00529 TMatrixD *vyyredinv=InvertMSparse(vyyred);
00530 DeleteMatrix(&vyyred);
00531
00532 nn=ny*ny;
00533 Int_t *vyyinv_rows=new Int_t[nn];
00534 Int_t *vyyinv_cols=new Int_t[nn];
00535 Double_t *vyyinv_data=new Double_t[nn];
00536 nn=0;
00537 for(Int_t ix=0;ix<n;ix++) {
00538 for(Int_t iy=0;iy<n;iy++) {
00539 Double_t d=(*vyyredinv)(ix,iy);
00540 if(d != 0.0) {
00541 vyyinv_rows[nn]=iToY[ix];
00542 vyyinv_cols[nn]=iToY[iy];
00543 vyyinv_data[nn]=d;
00544 nn++;
00545 }
00546 }
00547 }
00548 DeleteMatrix(&vyyredinv);
00549 TMatrixDSparse *Vyyinv=CreateSparseMatrix
00550 (ny,ny,nn,vyyinv_rows,vyyinv_cols,vyyinv_data);
00551
00552 delete[] vyyinv_rows;
00553 delete[] vyyinv_cols;
00554 delete[] vyyinv_data;
00555
00556 delete[] yToI;
00557 delete[] iToY;
00558
00559
00560
00561
00562
00563
00564 TMatrixDSparse *AtVyyinv=MultiplyMSparseTranspMSparse(fA,Vyyinv);
00565
00566
00567
00568
00569
00570 TMatrixDSparse *rhs=MultiplyMSparseM(AtVyyinv,fY);
00571 if (fBiasScale != 0.0) {
00572 TMatrixDSparse *rhs2=MultiplyMSparseM(fLsquared,fX0);
00573 AddMSparse(rhs, fTauSquared * fBiasScale ,rhs2);
00574 DeleteMatrix(&rhs2);
00575 }
00576
00577
00578
00579
00580
00581 fEinv=MultiplyMSparseMSparse(AtVyyinv,fA);
00582 AddMSparse(fEinv,fTauSquared,fLsquared);
00583
00584
00585
00586
00587
00588
00589 TMatrixD *EE = InvertMSparse(fEinv);
00590 fE=new TMatrixDSparse(*EE);
00591
00592
00593
00594
00595
00596 fX = new TMatrixD(*EE, TMatrixD::kMult, *rhs);
00597 DeleteMatrix(&rhs);
00598
00599
00600 Double_t lambda_half=0.0;
00601 Double_t one_over_epsEeps=0.0;
00602 TMatrixDSparse *epsilon=0;
00603 TMatrixDSparse *Eepsilon=0;
00604 if(fConstraint != kEConstraintNone) {
00605
00606 const Int_t *A_rows=fA->GetRowIndexArray();
00607 const Int_t *A_cols=fA->GetColIndexArray();
00608 const Double_t *A_data=fA->GetMatrixArray();
00609 TMatrixD epsilonNosparse(fA->GetNcols(),1);
00610 for(Int_t i=0;i<A_rows[fA->GetNrows()];i++) {
00611 epsilonNosparse(A_cols[i],0) += A_data[i];
00612 }
00613 epsilon=new TMatrixDSparse(epsilonNosparse);
00614
00615 Eepsilon=MultiplyMSparseMSparse(fE,epsilon);
00616
00617 TMatrixDSparse *epsilonEepsilon=MultiplyMSparseTranspMSparse(epsilon,
00618 Eepsilon);
00619
00620 if(epsilonEepsilon->GetRowIndexArray()[1]==1) {
00621 one_over_epsEeps=1./epsilonEepsilon->GetMatrixArray()[0];
00622 } else {
00623 Fatal("TUnfold::Unfold","epsilon#Eepsilon has dimension %d != 1",
00624 epsilonEepsilon->GetRowIndexArray()[1]);
00625 }
00626 DeleteMatrix(&epsilonEepsilon);
00627
00628 Double_t y_minus_epsx=0.0;
00629 for(Int_t iy=0;iy<fY->GetNrows();iy++) {
00630 y_minus_epsx += (*fY)(iy,0);
00631 }
00632
00633 for(Int_t ix=0;ix<epsilonNosparse.GetNrows();ix++) {
00634 y_minus_epsx -= epsilonNosparse(ix,0) * (*fX)(ix,0);
00635 }
00636
00637 lambda_half=y_minus_epsx*one_over_epsEeps;
00638
00639 const Int_t *EEpsilon_rows=Eepsilon->GetRowIndexArray();
00640 const Double_t *EEpsilon_data=Eepsilon->GetMatrixArray();
00641 for(Int_t ix=0;ix<Eepsilon->GetNrows();ix++) {
00642 if(EEpsilon_rows[ix]<EEpsilon_rows[ix+1]) {
00643 (*fX)(ix,0) += lambda_half * EEpsilon_data[EEpsilon_rows[ix]];
00644 }
00645 }
00646 }
00647
00648
00649
00650
00651 fDXDY = MultiplyMSparseMSparse(fE,AtVyyinv);
00652
00653
00654 if(fConstraint != kEConstraintNone) {
00655
00656 Int_t *rows=new Int_t[GetNy()];
00657 Int_t *cols=new Int_t[GetNy()];
00658 Double_t *data=new Double_t[GetNy()];
00659 for(Int_t i=0;i<GetNy();i++) {
00660 rows[i]=0;
00661 cols[i]=i;
00662 data[i]=one_over_epsEeps;
00663 }
00664 TMatrixDSparse *temp=CreateSparseMatrix
00665 (1,GetNy(),GetNy(),rows,cols,data);
00666 delete[] data;
00667 delete[] rows;
00668 delete[] cols;
00669
00670 TMatrixDSparse *epsilonB=MultiplyMSparseTranspMSparse(epsilon,fDXDY);
00671
00672 AddMSparse(temp, -one_over_epsEeps, epsilonB);
00673 DeleteMatrix(&epsilonB);
00674
00675 TMatrixDSparse *corr=MultiplyMSparseMSparse(Eepsilon,temp);
00676 DeleteMatrix(&temp);
00677
00678 AddMSparse(fDXDY,1.0,corr);
00679 DeleteMatrix(&corr);
00680 }
00681
00682 DeleteMatrix(&AtVyyinv);
00683
00684
00685
00686
00687 TMatrixDSparse *AE = MultiplyMSparseMSparse(fA,fE);
00688 fVxx = MultiplyMSparseMSparse(fDXDY,AE);
00689
00690 DeleteMatrix(&AE);
00691
00692
00693
00694
00695
00696 fAx = MultiplyMSparseM(fA,fX);
00697
00698
00699
00700
00701
00702 TMatrixD dy(*fY, TMatrixD::kMinus, *fAx);
00703 TMatrixDSparse *VyyinvDy=MultiplyMSparseM(Vyyinv,&dy);
00704 DeleteMatrix(&Vyyinv);
00705
00706 const Int_t *VyyinvDy_rows=VyyinvDy->GetRowIndexArray();
00707 const Double_t *VyyinvDy_data=VyyinvDy->GetMatrixArray();
00708 fChi2A=0.0;
00709 for(Int_t iy=0;iy<VyyinvDy->GetNrows();iy++) {
00710 if(VyyinvDy_rows[iy]<VyyinvDy_rows[iy+1]) {
00711 fChi2A += VyyinvDy_data[VyyinvDy_rows[iy]]*dy(iy,0);
00712 }
00713 }
00714 TMatrixD dx( fBiasScale * (*fX0), TMatrixD::kMinus,(*fX));
00715 TMatrixDSparse *LsquaredDx=MultiplyMSparseM(fLsquared,&dx);
00716 const Int_t *LsquaredDx_rows=LsquaredDx->GetRowIndexArray();
00717 const Double_t *LsquaredDx_data=LsquaredDx->GetMatrixArray();
00718 fLXsquared = 0.0;
00719 for(Int_t ix=0;ix<LsquaredDx->GetNrows();ix++) {
00720 if(LsquaredDx_rows[ix]<LsquaredDx_rows[ix+1]) {
00721 fLXsquared += LsquaredDx_data[LsquaredDx_rows[ix]]*dx(ix,0);
00722 }
00723 }
00724
00725
00726
00727 fDXDtauSquared=MultiplyMSparseMSparse(fE,LsquaredDx);
00728
00729 if(fConstraint != kEConstraintNone) {
00730 TMatrixDSparse *temp=MultiplyMSparseTranspMSparse(epsilon,fDXDtauSquared);
00731 Double_t f=0.0;
00732 if(temp->GetRowIndexArray()[1]==1) {
00733 f=temp->GetMatrixArray()[0]*one_over_epsEeps;
00734 } else if(temp->GetRowIndexArray()[1]>1) {
00735 Fatal("TUnfold::Unfold",
00736 "epsilon#fDXDtauSquared has dimension %d != 1",
00737 temp->GetRowIndexArray()[1]);
00738 }
00739 if(f!=0.0) {
00740 AddMSparse(fDXDtauSquared, -f,Eepsilon);
00741 }
00742 DeleteMatrix(&temp);
00743 }
00744 DeleteMatrix(&epsilon);
00745
00746 DeleteMatrix(&LsquaredDx);
00747
00748
00749 fDXDAM[0]=new TMatrixDSparse(*fE);
00750 fDXDAM[1]=new TMatrixDSparse(*fDXDY);
00751 fDXDAZ[0]=VyyinvDy;
00752 VyyinvDy=0;
00753 fDXDAZ[1]=new TMatrixDSparse(*fX);
00754
00755 if(fConstraint != kEConstraintNone) {
00756
00757 TMatrixDSparse *temp1=MultiplyMSparseMSparseTranspVector
00758 (Eepsilon,Eepsilon,0);
00759 AddMSparse(fDXDAM[0], -one_over_epsEeps,temp1);
00760 DeleteMatrix(&temp1);
00761
00762 Int_t *rows=new Int_t[GetNy()];
00763 Int_t *cols=new Int_t[GetNy()];
00764 Double_t *data=new Double_t[GetNy()];
00765 for(Int_t i=0;i<GetNy();i++) {
00766 rows[i]=i;
00767 cols[i]=0;
00768 data[i]=lambda_half;
00769 }
00770 TMatrixDSparse *temp2=CreateSparseMatrix
00771 (GetNy(),1,GetNy(),rows,cols,data);
00772 delete[] data;
00773 delete[] rows;
00774 delete[] cols;
00775 AddMSparse(fDXDAZ[0],1.0,temp2);
00776 DeleteMatrix(&temp2);
00777 }
00778
00779 DeleteMatrix(&Eepsilon);
00780 DeleteMatrix(&EE);
00781
00782 TMatrixD *VxxINV = InvertMSparse(fVxx);
00783 fVxxInv=new TMatrixDSparse(*VxxINV);
00784
00785
00786 const Int_t *Vxx_rows=fVxx->GetRowIndexArray();
00787 const Int_t *Vxx_cols=fVxx->GetColIndexArray();
00788 const Double_t *Vxx_data=fVxx->GetMatrixArray();
00789
00790 Double_t rho_squared_max = 0.0;
00791 Double_t rho_sum = 0.0;
00792 Int_t n_rho=0;
00793 for (int ix = 0; ix < fVxx->GetNrows(); ix++) {
00794 for(int ik=Vxx_rows[ix];ik<Vxx_rows[ix+1];ik++) {
00795 if(ix==Vxx_cols[ik]) {
00796 Double_t rho_squared =
00797 1. - 1. / ((*VxxINV) (ix, ix) * Vxx_data[ik]);
00798 if (rho_squared > rho_squared_max)
00799 rho_squared_max = rho_squared;
00800 if(rho_squared>0.0) {
00801 rho_sum += TMath::Sqrt(rho_squared);
00802 n_rho++;
00803 }
00804 break;
00805 }
00806 }
00807 }
00808 fRhoMax = TMath::Sqrt(rho_squared_max);
00809 fRhoAvg = (n_rho>0) ? (rho_sum/n_rho) : -1.0;
00810
00811 delete VxxINV;
00812
00813 return fRhoMax;
00814 }
00815
00816 TMatrixDSparse *TUnfold::CreateSparseMatrix
00817 (Int_t nrow,Int_t ncol,Int_t nel,Int_t *row,Int_t *col,Double_t *data) const {
00818 TMatrixDSparse *A=new TMatrixDSparse(nrow,ncol);
00819 if(nel>0) {
00820 A->SetMatrixArray(nel,row,col,data);
00821 }
00822 return A;
00823 }
00824
00825 TMatrixDSparse *TUnfold::MultiplyMSparseMSparse(const TMatrixDSparse *a,
00826 const TMatrixDSparse *b) const
00827 {
00828
00829
00830
00831
00832 if(a->GetNcols()!=b->GetNrows()) {
00833 Fatal("MultiplyMSparseMSparse",
00834 "inconsistent matrix col/ matrix row %d !=%d",
00835 a->GetNcols(),b->GetNrows());
00836 }
00837
00838 TMatrixDSparse *r=new TMatrixDSparse(a->GetNrows(),b->GetNcols());
00839 const Int_t *a_rows=a->GetRowIndexArray();
00840 const Int_t *a_cols=a->GetColIndexArray();
00841 const Double_t *a_data=a->GetMatrixArray();
00842 const Int_t *b_rows=b->GetRowIndexArray();
00843 const Int_t *b_cols=b->GetColIndexArray();
00844 const Double_t *b_data=b->GetMatrixArray();
00845
00846 int nMax=0;
00847 for (Int_t irow = 0; irow < a->GetNrows(); irow++) {
00848 if(a_rows[irow+1]>a_rows[irow]) nMax += b->GetNcols();
00849 }
00850 if((nMax>0)&&(a_cols)&&(b_cols)) {
00851 Int_t *r_rows=new Int_t[nMax];
00852 Int_t *r_cols=new Int_t[nMax];
00853 Double_t *r_data=new Double_t[nMax];
00854 Double_t *row_data=new Double_t[b->GetNcols()];
00855 Int_t n=0;
00856 for (Int_t irow = 0; irow < a->GetNrows(); irow++) {
00857 if(a_rows[irow+1]<=a_rows[irow]) continue;
00858
00859 for(Int_t icol=0;icol<b->GetNcols();icol++) {
00860 row_data[icol]=0.0;
00861 }
00862
00863 for(Int_t ia=a_rows[irow];ia<a_rows[irow+1];ia++) {
00864 Int_t k=a_cols[ia];
00865
00866 for(Int_t ib=b_rows[k];ib<b_rows[k+1];ib++) {
00867 row_data[b_cols[ib]] += a_data[ia]*b_data[ib];
00868 }
00869 }
00870
00871 for(Int_t icol=0;icol<b->GetNcols();icol++) {
00872 if(row_data[icol] != 0.0) {
00873 r_rows[n]=irow;
00874 r_cols[n]=icol;
00875 r_data[n]=row_data[icol];
00876 n++;
00877 }
00878 }
00879 }
00880 if(n>0) {
00881 r->SetMatrixArray(n,r_rows,r_cols,r_data);
00882 }
00883 delete[] r_rows;
00884 delete[] r_cols;
00885 delete[] r_data;
00886 delete[] row_data;
00887 }
00888
00889 return r;
00890 }
00891
00892
00893 TMatrixDSparse *TUnfold::MultiplyMSparseTranspMSparse(const TMatrixDSparse *a,
00894 const TMatrixDSparse *b) const
00895 {
00896
00897
00898
00899
00900
00901
00902 if(a->GetNrows() != b->GetNrows()) {
00903 Fatal("MultiplyMSparseTranspMSparse",
00904 "inconsistent matrix row numbers %d!=%d",
00905 a->GetNrows(),b->GetNrows());
00906 }
00907
00908 TMatrixDSparse *r=new TMatrixDSparse(a->GetNcols(),b->GetNcols());
00909 const Int_t *a_rows=a->GetRowIndexArray();
00910 const Int_t *a_cols=a->GetColIndexArray();
00911 const Double_t *a_data=a->GetMatrixArray();
00912 const Int_t *b_rows=b->GetRowIndexArray();
00913 const Int_t *b_cols=b->GetColIndexArray();
00914 const Double_t *b_data=b->GetMatrixArray();
00915
00916
00917
00918 typedef std::map<Int_t,Double_t> MMatrixRow_t;
00919 typedef std::map<Int_t, MMatrixRow_t > MMatrix_t;
00920 MMatrix_t matrix;
00921
00922 for(Int_t iRowAB=0;iRowAB<a->GetNrows();iRowAB++) {
00923 for(Int_t ia=a_rows[iRowAB];ia<a_rows[iRowAB+1];ia++) {
00924 for(Int_t ib=b_rows[iRowAB];ib<b_rows[iRowAB+1];ib++) {
00925
00926 MMatrixRow_t &row=matrix[a_cols[ia]];
00927 MMatrixRow_t::iterator icol=row.find(b_cols[ib]);
00928 if(icol!=row.end()) {
00929
00930 (*icol).second += a_data[ia]*b_data[ib];
00931 } else {
00932
00933 row[b_cols[ib]] = a_data[ia]*b_data[ib];
00934 }
00935 }
00936 }
00937 }
00938
00939 Int_t n=0;
00940 for(MMatrix_t::const_iterator irow=matrix.begin();
00941 irow!=matrix.end();irow++) {
00942 n += (*irow).second.size();
00943 }
00944 if(n>0) {
00945
00946 Int_t *r_rows=new Int_t[n];
00947 Int_t *r_cols=new Int_t[n];
00948 Double_t *r_data=new Double_t[n];
00949 n=0;
00950 for(MMatrix_t::const_iterator irow=matrix.begin();
00951 irow!=matrix.end();irow++) {
00952 for(MMatrixRow_t::const_iterator icol=(*irow).second.begin();
00953 icol!=(*irow).second.end();icol++) {
00954 r_rows[n]=(*irow).first;
00955 r_cols[n]=(*icol).first;
00956 r_data[n]=(*icol).second;
00957 n++;
00958 }
00959 }
00960
00961 if(n>0) {
00962 r->SetMatrixArray(n,r_rows,r_cols,r_data);
00963 }
00964 delete[] r_rows;
00965 delete[] r_cols;
00966 delete[] r_data;
00967 }
00968
00969 return r;
00970 }
00971
00972 TMatrixDSparse *TUnfold::MultiplyMSparseM(const TMatrixDSparse *a,
00973 const TMatrixD *b) const
00974 {
00975
00976
00977
00978
00979
00980 if(a->GetNcols()!=b->GetNrows()) {
00981 Fatal("MultiplyMSparseM","inconsistent matrix col /matrix row %d!=%d",
00982 a->GetNcols(),b->GetNrows());
00983 }
00984
00985 TMatrixDSparse *r=new TMatrixDSparse(a->GetNrows(),b->GetNcols());
00986 const Int_t *a_rows=a->GetRowIndexArray();
00987 const Int_t *a_cols=a->GetColIndexArray();
00988 const Double_t *a_data=a->GetMatrixArray();
00989
00990 int nMax=0;
00991 for (Int_t irow = 0; irow < a->GetNrows(); irow++) {
00992 if(a_rows[irow+1]-a_rows[irow]>0) nMax += b->GetNcols();
00993 }
00994 if(nMax>0) {
00995 Int_t *r_rows=new Int_t[nMax];
00996 Int_t *r_cols=new Int_t[nMax];
00997 Double_t *r_data=new Double_t[nMax];
00998
00999 Int_t n=0;
01000
01001 for (Int_t irow = 0; irow < a->GetNrows(); irow++) {
01002 if(a_rows[irow+1]-a_rows[irow]<=0) continue;
01003 for(Int_t icol=0;icol<b->GetNcols();icol++) {
01004 r_rows[n]=irow;
01005 r_cols[n]=icol;
01006 r_data[n]=0.0;
01007 for(Int_t i=a_rows[irow];i<a_rows[irow+1];i++) {
01008 Int_t j=a_cols[i];
01009 r_data[n] += a_data[i]*(*b)(j,icol);
01010 }
01011 if(r_data[n]!=0.0) n++;
01012 }
01013 }
01014 if(n>0) {
01015 r->SetMatrixArray(n,r_rows,r_cols,r_data);
01016 }
01017 delete[] r_rows;
01018 delete[] r_cols;
01019 delete[] r_data;
01020 }
01021 return r;
01022 }
01023
01024 TMatrixDSparse *TUnfold::MultiplyMSparseMSparseTranspVector
01025 (const TMatrixDSparse *m1,const TMatrixDSparse *m2,
01026 const TMatrixTBase<Double_t> *v) const {
01027
01028
01029
01030
01031 if((m1->GetNcols() != m2->GetNcols())||
01032 (v && ((m1->GetNcols()!=v->GetNrows())||(v->GetNcols()!=1)))) {
01033 if(v) {
01034 Fatal("MultiplyMSparseMSparseTranspVector",
01035 "matrix cols/vector rows %d!=%d!=%d or vector rows %d!=1\n",
01036 m1->GetNcols(),m2->GetNcols(),v->GetNrows(),v->GetNcols());
01037 } else {
01038 Fatal("MultiplyMSparseMSparseTranspVector",
01039 "matrix cols %d!=%d\n",m1->GetNcols(),m2->GetNcols());
01040 }
01041 }
01042 const Int_t *rows_m1=m1->GetRowIndexArray();
01043 const Int_t *cols_m1=m1->GetColIndexArray();
01044 const Double_t *data_m1=m1->GetMatrixArray();
01045 Int_t num_m1=0;
01046 for(Int_t i=0;i<m1->GetNrows();i++) {
01047 if(rows_m1[i]<rows_m1[i+1]) num_m1++;
01048 }
01049 const Int_t *rows_m2=m2->GetRowIndexArray();
01050 const Int_t *cols_m2=m2->GetColIndexArray();
01051 const Double_t *data_m2=m2->GetMatrixArray();
01052 Int_t num_m2=0;
01053 for(Int_t j=0;j<m2->GetNrows();j++) {
01054 if(rows_m2[j]<rows_m2[j+1]) num_m2++;
01055 }
01056 const TMatrixDSparse *v_sparse=dynamic_cast<const TMatrixDSparse *>(v);
01057 const Int_t *v_rows=0;
01058 const Double_t *v_data=0;
01059 if(v_sparse) {
01060 v_rows=v_sparse->GetRowIndexArray();
01061 v_data=v_sparse->GetMatrixArray();
01062 }
01063 Int_t num_r=num_m1*num_m2+1;
01064 Int_t *row_r=new Int_t[num_r];
01065 Int_t *col_r=new Int_t[num_r];
01066 Double_t *data_r=new Double_t[num_r];
01067 num_r=0;
01068 for(Int_t i=0;i<m1->GetNrows();i++) {
01069 for(Int_t j=0;j<m2->GetNrows();j++) {
01070 data_r[num_r]=0.0;
01071 Int_t index_m1=rows_m1[i];
01072 Int_t index_m2=rows_m2[j];
01073 while((index_m1<rows_m1[i+1])&&(index_m2<rows_m2[j+1])) {
01074 Int_t k1=cols_m1[index_m1];
01075 Int_t k2=cols_m2[index_m2];
01076 if(k1<k2) {
01077 index_m1++;
01078 } else if(k1>k2) {
01079 index_m2++;
01080 } else {
01081 if(v_sparse) {
01082 Int_t v_index=v_rows[k1];
01083 if(v_index<v_rows[k1+1]) {
01084 data_r[num_r] += data_m1[index_m1] * data_m2[index_m2]
01085 * v_data[v_index];
01086 } else {
01087 data_r[num_r] =0.0;
01088 }
01089 } else if(v) {
01090 data_r[num_r] += data_m1[index_m1] * data_m2[index_m2]
01091 * (*v)(k1,0);
01092 } else {
01093 data_r[num_r] += data_m1[index_m1] * data_m2[index_m2];
01094 }
01095 index_m1++;
01096 index_m2++;
01097 }
01098 }
01099 if(data_r[num_r] !=0.0) {
01100 row_r[num_r]=i;
01101 col_r[num_r]=j;
01102 num_r++;
01103 }
01104 }
01105 }
01106 TMatrixDSparse *r=CreateSparseMatrix(m1->GetNrows(),m2->GetNrows(),
01107 num_r,row_r,col_r,data_r);
01108 delete[] row_r;
01109 delete[] col_r;
01110 delete[] data_r;
01111 return r;
01112 }
01113
01114 void TUnfold::AddMSparse(TMatrixDSparse *dest,Double_t f,
01115 const TMatrixDSparse *src) {
01116
01117
01118 const Int_t *dest_rows=dest->GetRowIndexArray();
01119 const Int_t *dest_cols=dest->GetColIndexArray();
01120 const Double_t *dest_data=dest->GetMatrixArray();
01121 const Int_t *src_rows=src->GetRowIndexArray();
01122 const Int_t *src_cols=src->GetColIndexArray();
01123 const Double_t *src_data=src->GetMatrixArray();
01124
01125 if((dest->GetNrows()!=src->GetNrows())||
01126 (dest->GetNcols()!=src->GetNcols())) {
01127 Fatal("AddMSparse","inconsistent matrix rows %d!=%d OR cols %d!=%d",
01128 src->GetNrows(),dest->GetNrows(),src->GetNcols(),dest->GetNcols());
01129 }
01130 Int_t nmax=dest->GetNrows()*dest->GetNcols();
01131 Double_t *result_data=new Double_t[nmax];
01132 Int_t *result_rows=new Int_t[nmax];
01133 Int_t *result_cols=new Int_t[nmax];
01134 Int_t n=0;
01135 for(Int_t row=0;row<dest->GetNrows();row++) {
01136 Int_t i_dest=dest_rows[row];
01137 Int_t i_src=src_rows[row];
01138 while((i_dest<dest_rows[row+1])||(i_src<src_rows[row+1])) {
01139 Int_t col_dest=(i_dest<dest_rows[row+1]) ?
01140 dest_cols[i_dest] : dest->GetNcols();
01141 Int_t col_src =(i_src <src_rows[row+1] ) ?
01142 src_cols [i_src] : src->GetNcols();
01143 result_rows[n]=row;
01144 if(col_dest<col_src) {
01145 result_cols[n]=col_dest;
01146 result_data[n]=dest_data[i_dest++];
01147 } else if(col_dest>col_src) {
01148 result_cols[n]=col_src;
01149 result_data[n]=f*src_data[i_src++];
01150 } else {
01151 result_cols[n]=col_dest;
01152 result_data[n]=dest_data[i_dest++]+f*src_data[i_src++];
01153 }
01154 if(result_data[n] !=0.0) n++;
01155 }
01156 }
01157 if(n<=0) {
01158 n=1;
01159 result_rows[0]=0;
01160 result_cols[0]=0;
01161 result_data[0]=0.0;
01162 }
01163 dest->SetMatrixArray(n,result_rows,result_cols,result_data);
01164 delete[] result_data;
01165 delete[] result_rows;
01166 delete[] result_cols;
01167 }
01168
01169 TMatrixD *TUnfold::InvertMSparse(const TMatrixDSparse *A) const {
01170
01171
01172
01173
01174
01175
01176
01177 if(A->GetNcols()!=A->GetNrows()) {
01178 Fatal("InvertMSparse","inconsistent matrix row/col %d!=%d",
01179 A->GetNcols(),A->GetNrows());
01180 }
01181
01182 #ifdef SCHUR_COMPLEMENT_MATRIX_INVERSION
01183
01184 const Int_t *a_rows=A->GetRowIndexArray();
01185 const Int_t *a_cols=A->GetColIndexArray();
01186 const Double_t *a_data=A->GetMatrixArray();
01187
01188 Int_t nmin=0,nmax=0;
01189
01190 for(Int_t imin=0;imin<A->GetNrows();imin++) {
01191 Int_t imax=A->GetNrows();
01192 for(Int_t i2=imin;i2<imax;i2++) {
01193 for(Int_t i=a_rows[i2];i<a_rows[i2+1];i++) {
01194 if(a_data[i]==0.0) continue;
01195 Int_t icol=a_cols[i];
01196 if(icol<imin) continue;
01197 if(icol==i2) continue;
01198 if(icol<i2) {
01199
01200 imax=i2;
01201 break;
01202 } else {
01203
01204 imax=icol;
01205 break;
01206 }
01207 }
01208 }
01209 if(imax-imin>nmax-nmin) {
01210 nmin=imin;
01211 nmax=imax;
01212 }
01213 }
01214
01215 if(nmin>=nmax-1) {
01216 TMatrixD *r=new TMatrixD(*A);
01217 if(!InvertMConditioned(r)) {
01218 Fatal("InvertMSparse","InvertMConditioned(full matrix) failed");
01219 }
01220 return r;
01221 } else if((nmin==0)&&(nmax==A->GetNrows())) {
01222
01223
01224 TMatrixD *r=new TMatrixD(A->GetNrows(),A->GetNcols());
01225 Int_t error=0;
01226 for(Int_t irow=nmin;irow<nmax;irow++) {
01227 for(Int_t i=a_rows[irow];i<a_rows[irow+1];i++) {
01228 if(a_cols[i]==irow) {
01229 if(a_data[i] !=0.0) (*r)(irow,irow)=1./a_data[i];
01230 else error++;
01231 }
01232 }
01233 }
01234 if(error) {
01235 Error("InvertMSparse",
01236 "inversion failed (diagonal matrix) nerror=%d",error);
01237 }
01238 return r;
01239 }
01240
01241
01242
01243
01244
01245
01246 std::vector<Double_t> Dinv;
01247 Dinv.resize(nmax-nmin);
01248 Int_t error=0;
01249 for(Int_t irow=nmin;irow<nmax;irow++) {
01250 for(Int_t i=a_rows[irow];i<a_rows[irow+1];i++) {
01251 if(a_cols[i]==irow) {
01252 if(a_data[i]!=0.0) {
01253 Dinv[irow-nmin]=1./a_data[i];
01254 } else {
01255 Dinv[irow-nmin]=0.0;
01256 error++;
01257 }
01258 break;
01259 }
01260 }
01261 }
01262 if(error) {
01263 Error("InvertMSparse",
01264 "inversion failed (diagonal part) nerror=%d",error);
01265 }
01266
01267 Int_t nBDinv=0,nC=0;
01268 for(Int_t irow_a=0;irow_a<A->GetNrows();irow_a++) {
01269 if((irow_a<nmin)||(irow_a>=nmax)) {
01270 for(Int_t i=a_rows[irow_a];i<a_rows[irow_a+1];i++) {
01271 Int_t icol=a_cols[i];
01272 if((icol>=nmin)&&(icol<nmax)) nBDinv++;
01273 }
01274 } else {
01275 for(Int_t i=a_rows[irow_a];i<a_rows[irow_a+1];i++) {
01276 Int_t icol = a_cols[i];
01277 if((icol<nmin)||(icol>=nmax)) nC++;
01278 }
01279 }
01280 }
01281 Int_t *row_BDinv=new Int_t[nBDinv+1];
01282 Int_t *col_BDinv=new Int_t[nBDinv+1];
01283 Double_t *data_BDinv=new Double_t[nBDinv+1];
01284
01285 Int_t *row_C=new Int_t[nC+1];
01286 Int_t *col_C=new Int_t[nC+1];
01287 Double_t *data_C=new Double_t[nC+1];
01288
01289 TMatrixD Aschur(A->GetNrows()-(nmax-nmin),A->GetNcols()-(nmax-nmin));
01290
01291 nBDinv=0;
01292 nC=0;
01293 for(Int_t irow_a=0;irow_a<A->GetNrows();irow_a++) {
01294 if((irow_a<nmin)||(irow_a>=nmax)) {
01295 Int_t row=(irow_a<nmin) ? irow_a : (irow_a-(nmax-nmin));
01296 for(Int_t i=a_rows[irow_a];i<a_rows[irow_a+1];i++) {
01297 Int_t icol_a=a_cols[i];
01298 if(icol_a<nmin) {
01299 Aschur(row,icol_a)=a_data[i];
01300 } else if(icol_a>=nmax) {
01301 Aschur(row,icol_a-(nmax-nmin))=a_data[i];
01302 } else {
01303 row_BDinv[nBDinv]=row;
01304 col_BDinv[nBDinv]=icol_a-nmin;
01305 data_BDinv[nBDinv]=a_data[i]*Dinv[icol_a-nmin];
01306 nBDinv++;
01307 }
01308 }
01309 } else {
01310 for(Int_t i=a_rows[irow_a];i<a_rows[irow_a+1];i++) {
01311 Int_t icol_a=a_cols[i];
01312 if(icol_a<nmin) {
01313 row_C[nC]=irow_a-nmin;
01314 col_C[nC]=icol_a;
01315 data_C[nC]=a_data[i];
01316 nC++;
01317 } else if(icol_a>=nmax) {
01318 row_C[nC]=irow_a-nmin;
01319 col_C[nC]=icol_a-(nmax-nmin);
01320 data_C[nC]=a_data[i];
01321 nC++;
01322 }
01323 }
01324 }
01325 }
01326 TMatrixDSparse *BDinv=CreateSparseMatrix
01327 (A->GetNrows()-(nmax-nmin),nmax-nmin,
01328 nBDinv,row_BDinv,col_BDinv,data_BDinv);
01329 delete[] row_BDinv;
01330 delete[] col_BDinv;
01331 delete[] data_BDinv;
01332
01333
01334 TMatrixDSparse *C=CreateSparseMatrix(nmax-nmin,A->GetNcols()-(nmax-nmin),
01335 nC,row_C,col_C,data_C);
01336 delete[] row_C;
01337 delete[] col_C;
01338 delete[] data_C;
01339
01340 TMatrixDSparse *BDinvC=MultiplyMSparseMSparse(BDinv,C);
01341
01342 Aschur -= *BDinvC;
01343 if(!InvertMConditioned(&Aschur)) {
01344 Fatal("InvertMSparse","InvertMConditioned failed (part of matrix)");
01345 }
01346
01347 DeleteMatrix(&BDinvC);
01348
01349 TMatrixD *r=new TMatrixD(A->GetNrows(),A->GetNcols());
01350
01351 for(Int_t row_a=0;row_a<Aschur.GetNrows();row_a++) {
01352 for(Int_t col_a=0;col_a<Aschur.GetNcols();col_a++) {
01353 (*r)((row_a<nmin) ? row_a : (row_a+nmax-nmin),
01354 (col_a<nmin) ? col_a : (col_a+nmax-nmin))=Aschur(row_a,col_a);
01355 }
01356 }
01357
01358 TMatrixDSparse *CAschur=MultiplyMSparseM(C,&Aschur);
01359 TMatrixDSparse *CAschurBDinv=MultiplyMSparseMSparse(CAschur,BDinv);
01360
01361 DeleteMatrix(&C);
01362
01363 const Int_t *CAschurBDinv_row=CAschurBDinv->GetRowIndexArray();
01364 const Int_t *CAschurBDinv_col=CAschurBDinv->GetColIndexArray();
01365 const Double_t *CAschurBDinv_data=CAschurBDinv->GetMatrixArray();
01366 for(Int_t row=0;row<CAschurBDinv->GetNrows();row++) {
01367 for(Int_t i=CAschurBDinv_row[row];i<CAschurBDinv_row[row+1];i++) {
01368 Int_t col=CAschurBDinv_col[i];
01369 (*r)(row+nmin,col+nmin)=CAschurBDinv_data[i]*Dinv[row];
01370 }
01371 (*r)(row+nmin,row+nmin) += Dinv[row];
01372 }
01373
01374 DeleteMatrix(&CAschurBDinv);
01375
01376 const Int_t *CAschur_row=CAschur->GetRowIndexArray();
01377 const Int_t *CAschur_col=CAschur->GetColIndexArray();
01378 const Double_t *CAschur_data=CAschur->GetMatrixArray();
01379 for(Int_t row=0;row<CAschur->GetNrows();row++) {
01380 for(Int_t i=CAschur_row[row];i<CAschur_row[row+1];i++) {
01381 Int_t col=CAschur_col[i];
01382 (*r)(row+nmin,
01383 (col<nmin) ? col : (col+nmax-nmin))= -CAschur_data[i]*Dinv[row];
01384 }
01385 }
01386 DeleteMatrix(&CAschur);
01387
01388 const Int_t *BDinv_row=BDinv->GetRowIndexArray();
01389 const Int_t *BDinv_col=BDinv->GetColIndexArray();
01390 const Double_t *BDinv_data=BDinv->GetMatrixArray();
01391 for(Int_t row_aschur=0;row_aschur<Aschur.GetNrows();row_aschur++) {
01392 Int_t row=(row_aschur<nmin) ? row_aschur : (row_aschur+nmax-nmin);
01393 for(Int_t row_bdinv=0;row_bdinv<BDinv->GetNrows();row_bdinv++) {
01394 for(Int_t i=BDinv_row[row_bdinv];i<BDinv_row[row_bdinv+1];i++) {
01395 (*r)(row,BDinv_col[i]+nmin) -= Aschur(row_aschur,row_bdinv)*
01396 BDinv_data[i];
01397 }
01398 }
01399 }
01400 DeleteMatrix(&BDinv);
01401
01402 return r;
01403 #else
01404 TMatrixD *r=new TMatrixD(A);
01405 if(!InvertMConditioned(*r)) {
01406 Fatal("InvertMSparse","InvertMConditioned failed (full matrix)";
01407 print_backtrace();
01408 }
01409 return r;
01410 #endif
01411 }
01412
01413 Bool_t TUnfold::InvertMConditioned(TMatrixD *A) {
01414
01415
01416
01417
01418
01419
01420
01421 #ifdef INVERT_WITH_CONDITIONING
01422
01423 Double_t *A_diagonals=new Double_t[A->GetNrows()];
01424 for(Int_t i=0;i<A->GetNrows();i++) {
01425 A_diagonals[i]=TMath::Sqrt(TMath::Abs((*A)(i,i)));
01426 if(A_diagonals[i]>0.0) A_diagonals[i]=1./A_diagonals[i];
01427 else A_diagonals[i]=1.0;
01428 }
01429
01430 for(Int_t i=0;i<A->GetNrows();i++) {
01431 for(Int_t j=0;j<A->GetNcols();j++) {
01432 (*A)(i,j) *= A_diagonals[i]*A_diagonals[j];
01433 }
01434 }
01435 #endif
01436 Double_t det=0.0;
01437 A->Invert(&det);
01438 #ifdef INVERT_WITH_CONDITIONING
01439
01440 for(Int_t i=0;i<A->GetNrows();i++) {
01441 for(Int_t j=0;j<A->GetNcols();j++) {
01442 (*A)(i,j) *= A_diagonals[i]*A_diagonals[j];
01443 }
01444 }
01445 delete[] A_diagonals;
01446 #endif
01447 return (det !=0.0);
01448 }
01449
01450 TUnfold::TUnfold(const TH2 *hist_A, EHistMap histmap, ERegMode regmode,
01451 TUnfold::EConstraint constraint)
01452 {
01453
01454
01455
01456
01457
01458
01459
01460
01461
01462
01463
01464
01465
01466
01467
01468
01469
01470
01471
01472 InitTUnfold();
01473 SetConstraint(constraint);
01474 Int_t nx0, nx, ny;
01475 if (histmap == kHistMapOutputHoriz) {
01476
01477 nx0 = hist_A->GetNbinsX()+2;
01478 ny = hist_A->GetNbinsY();
01479 } else {
01480
01481 nx0 = hist_A->GetNbinsY()+2;
01482 ny = hist_A->GetNbinsX();
01483 }
01484 nx = 0;
01485
01486
01487
01488 fSumOverY.Set(nx0);
01489 fXToHist.Set(nx0);
01490 fHistToX.Set(nx0);
01491 Int_t nonzeroA=0;
01492
01493
01494
01495
01496
01497
01498
01499
01500
01501
01502 Int_t nskipped=0;
01503 for (Int_t ix = 0; ix < nx0; ix++) {
01504
01505
01506 Double_t sum = 0.0;
01507 for (Int_t iy = 0; iy < ny; iy++) {
01508 Double_t z;
01509 if (histmap == kHistMapOutputHoriz) {
01510 z = hist_A->GetBinContent(ix, iy + 1);
01511 } else {
01512 z = hist_A->GetBinContent(iy + 1, ix);
01513 }
01514 if (z > 0.0) {
01515 nonzeroA++;
01516 sum += z;
01517 }
01518 }
01519
01520 if (sum > 0.0) {
01521
01522 fXToHist[nx] = ix;
01523 fHistToX[ix] = nx;
01524
01525
01526 fSumOverY[nx] = sum;
01527 if (histmap == kHistMapOutputHoriz) {
01528 fSumOverY[nx] +=
01529 hist_A->GetBinContent(ix, 0) +
01530 hist_A->GetBinContent(ix, ny + 1);
01531 } else {
01532 fSumOverY[nx] +=
01533 hist_A->GetBinContent(0, ix) +
01534 hist_A->GetBinContent(ny + 1, ix);
01535 }
01536 nx++;
01537 } else {
01538 nskipped++;
01539
01540 fHistToX[ix] = -1;
01541 }
01542 }
01543 Int_t underflowBin=0,overflowBin=0;
01544 for (Int_t ix = 0; ix < nx; ix++) {
01545 Int_t ibinx = fXToHist[ix];
01546 if(ibinx<1) underflowBin=1;
01547 if (histmap == kHistMapOutputHoriz) {
01548 if(ibinx>hist_A->GetNbinsX()) overflowBin=1;
01549 } else {
01550 if(ibinx>hist_A->GetNbinsY()) overflowBin=1;
01551 }
01552 }
01553 if(nskipped) {
01554 Int_t nprint=0;
01555 Int_t ixfirst=-1,ixlast=-1;
01556 TString binlist;
01557 for (Int_t ix = 0; ix < nx0; ix++) {
01558 if(fHistToX[ix]<0) {
01559 nprint++;
01560 if(ixlast<0) {
01561 binlist +=" ";
01562 binlist +=ix;
01563 ixfirst=ix;
01564 }
01565 ixlast=ix;
01566 }
01567 if(((fHistToX[ix]>=0)&&(ixlast>=0))||
01568 (nprint==nskipped)) {
01569 if(ixlast>ixfirst) {
01570 binlist += "-";
01571 binlist += ixlast;
01572 }
01573 ixfirst=-1;
01574 ixlast=-1;
01575 }
01576 if(nprint==nskipped) {
01577 break;
01578 }
01579 }
01580 if(nskipped==(2-underflowBin-overflowBin)) {
01581 Info("TUnfold","the following output bins "
01582 "are not connected to the input side %s",
01583 static_cast<const char *>(binlist));
01584 } else {
01585 Warning("TUnfold","the following output bins "
01586 "are not connected to the input side %s",
01587 static_cast<const char *>(binlist));
01588 }
01589 }
01590
01591 fX0 = new TMatrixD(nx, 1, fSumOverY.GetArray());
01592
01593
01594 Int_t *rowA = new Int_t[nonzeroA];
01595 Int_t *colA = new Int_t[nonzeroA];
01596 Double_t *dataA = new Double_t[nonzeroA];
01597 Int_t index=0;
01598 for (Int_t iy = 0; iy < ny; iy++) {
01599 for (Int_t ix = 0; ix < nx; ix++) {
01600 Int_t ibinx = fXToHist[ix];
01601 Double_t z;
01602 if (histmap == kHistMapOutputHoriz) {
01603 z = hist_A->GetBinContent(ibinx, iy + 1);
01604 } else {
01605 z = hist_A->GetBinContent(iy + 1, ibinx);
01606 }
01607 if (z > 0.0) {
01608 rowA[index]=iy;
01609 colA[index] = ix;
01610 dataA[index] = z / fSumOverY[ix];
01611 index++;
01612 }
01613 }
01614 }
01615 if(underflowBin && overflowBin) {
01616 Info("TUnfold","%d input bins and %d output bins (includes 2 underflow/overflow bins)",ny,nx);
01617 } else if(underflowBin) {
01618 Info("TUnfold","%d input bins and %d output bins (includes 1 underflow bin)",ny,nx);
01619 } else if(overflowBin) {
01620 Info("TUnfold","%d input bins and %d output bins (includes 1 overflow bin)",ny,nx);
01621 } else {
01622 Info("TUnfold","%d input bins and %d output bins",ny,nx);
01623 }
01624 fA = CreateSparseMatrix(ny,nx,index,rowA,colA,dataA);
01625 if(ny<nx) {
01626 Error("TUnfold","too few (ny=%d) input bins for nx=%d output bins",ny,nx);
01627 } else if(ny==nx) {
01628 Warning("TUnfold","too few (ny=%d) input bins for nx=%d output bins",ny,nx);
01629 }
01630 delete[] rowA;
01631 delete[] colA;
01632 delete[] dataA;
01633
01634 fLsquared = new TMatrixDSparse(GetNx(), GetNx());
01635 if (regmode != kRegModeNone) {
01636 Int_t nError=RegularizeBins(0, 1, nx0, regmode);
01637 if(nError>0) {
01638 if(nError>1) {
01639 Info("TUnfold",
01640 "%d regularisation conditions have been skipped",nError);
01641 } else {
01642 Info("TUnfold",
01643 "One regularisation condition has been skipped");
01644 }
01645 }
01646 }
01647 }
01648
01649 TUnfold::~TUnfold(void)
01650 {
01651
01652
01653 DeleteMatrix(&fA);
01654 DeleteMatrix(&fLsquared);
01655 DeleteMatrix(&fVyy);
01656 DeleteMatrix(&fY);
01657 DeleteMatrix(&fX0);
01658
01659 ClearResults();
01660 }
01661
01662 void TUnfold::SetBias(const TH1 *bias)
01663 {
01664
01665
01666 DeleteMatrix(&fX0);
01667 fX0 = new TMatrixD(GetNx(), 1);
01668 for (Int_t i = 0; i < GetNx(); i++) {
01669 (*fX0) (i, 0) = bias->GetBinContent(fXToHist[i]);
01670 }
01671 }
01672
01673 Int_t TUnfold::RegularizeSize(int bin, Double_t scale)
01674 {
01675
01676
01677
01678
01679
01680
01681 if(fRegMode==kRegModeNone) fRegMode=kRegModeSize;
01682 if(fRegMode!=kRegModeSize) fRegMode=kRegModeMixed;
01683
01684 Int_t i = fHistToX[bin];
01685 if (i < 0) {
01686 return 1;
01687 }
01688 (*fLsquared) (i, i) = scale * scale;
01689 return 0;
01690 }
01691
01692 Int_t TUnfold::RegularizeDerivative(int left_bin, int right_bin,
01693 Double_t scale)
01694 {
01695
01696
01697
01698
01699
01700
01701
01702 if(fRegMode==kRegModeNone) fRegMode=kRegModeDerivative;
01703 if(fRegMode!=kRegModeDerivative) fRegMode=kRegModeMixed;
01704
01705 Int_t il = fHistToX[left_bin];
01706 Int_t ir = fHistToX[right_bin];
01707 if ((il < 0) || (ir < 0)) {
01708 return 1;
01709 }
01710 Double_t scale_squared = scale * scale;
01711 (*fLsquared) (il, il) += scale_squared;
01712 (*fLsquared) (il, ir) -= scale_squared;
01713 (*fLsquared) (ir, il) -= scale_squared;
01714 (*fLsquared) (ir, ir) += scale_squared;
01715 return 0;
01716 }
01717
01718 Int_t TUnfold::RegularizeCurvature(int left_bin, int center_bin,
01719 int right_bin,
01720 Double_t scale_left,
01721 Double_t scale_right)
01722 {
01723
01724
01725
01726
01727
01728
01729
01730
01731
01732 if(fRegMode==kRegModeNone) fRegMode=kRegModeCurvature;
01733 if(fRegMode!=kRegModeCurvature) fRegMode=kRegModeMixed;
01734
01735 Int_t il, ic, ir;
01736 il = fHistToX[left_bin];
01737 ic = fHistToX[center_bin];
01738 ir = fHistToX[right_bin];
01739 if ((il < 0) || (ic < 0) || (ir < 0)) {
01740 return 1;
01741 }
01742 Double_t r[3];
01743 r[0] = -scale_left;
01744 r[2] = -scale_right;
01745 r[1] = -r[0] - r[2];
01746
01747 (*fLsquared) (il, il) += r[0] * r[0];
01748 (*fLsquared) (il, ic) += r[0] * r[1];
01749 (*fLsquared) (il, ir) += r[0] * r[2];
01750 (*fLsquared) (ic, il) += r[1] * r[0];
01751 (*fLsquared) (ic, ic) += r[1] * r[1];
01752 (*fLsquared) (ic, ir) += r[1] * r[2];
01753 (*fLsquared) (ir, il) += r[2] * r[0];
01754 (*fLsquared) (ir, ic) += r[2] * r[1];
01755 (*fLsquared) (ir, ir) += r[2] * r[2];
01756 return 0;
01757 }
01758
01759 Int_t TUnfold::RegularizeBins(int start, int step, int nbin,
01760 ERegMode regmode)
01761 {
01762
01763
01764
01765
01766
01767
01768
01769
01770
01771 Int_t i0, i1, i2;
01772 i0 = start;
01773 i1 = i0 + step;
01774 i2 = i1 + step;
01775 Int_t nSkip = 0;
01776 Int_t nError= 0;
01777 if (regmode == kRegModeDerivative) {
01778 nSkip = 1;
01779 } else if (regmode == kRegModeCurvature) {
01780 nSkip = 2;
01781 } else if(regmode != kRegModeSize) {
01782 Error("TUnfold::RegularizeBins","regmode = %d is not valid",regmode);
01783 }
01784 for (Int_t i = nSkip; i < nbin; i++) {
01785 if (regmode == kRegModeSize) {
01786 nError += RegularizeSize(i0);
01787 } else if (regmode == kRegModeDerivative) {
01788 nError += RegularizeDerivative(i0, i1);
01789 } else if (regmode == kRegModeCurvature) {
01790 nError += RegularizeCurvature(i0, i1, i2);
01791 }
01792 i0 = i1;
01793 i1 = i2;
01794 i2 += step;
01795 }
01796 return nError;
01797 }
01798
01799 Int_t TUnfold::RegularizeBins2D(int start_bin, int step1, int nbin1,
01800 int step2, int nbin2, ERegMode regmode)
01801 {
01802
01803
01804
01805
01806
01807
01808
01809
01810
01811
01812 Int_t nError = 0;
01813 for (Int_t i1 = 0; i1 < nbin1; i1++) {
01814 nError += RegularizeBins(start_bin + step1 * i1, step2, nbin2, regmode);
01815 }
01816 for (Int_t i2 = 0; i2 < nbin2; i2++) {
01817 nError += RegularizeBins(start_bin + step2 * i2, step1, nbin1, regmode);
01818 }
01819 return nError;
01820 }
01821
01822 Double_t TUnfold::DoUnfold(Double_t tau_reg,const TH1 *input,
01823 Double_t scaleBias)
01824 {
01825
01826
01827
01828
01829
01830
01831
01832
01833
01834
01835
01836
01837
01838
01839
01840 SetInput(input,scaleBias);
01841 return DoUnfold(tau_reg);
01842 }
01843
01844 Int_t TUnfold::SetInput(const TH1 *input, Double_t scaleBias,
01845 Double_t oneOverZeroError) {
01846
01847
01848
01849
01850
01851
01852
01853
01854
01855
01856
01857
01858 fBiasScale = scaleBias;
01859
01860
01861 ClearResults();
01862
01863
01864
01865
01866 fNdf = -GetNpar();
01867 Int_t *rowColVyy=new Int_t[GetNy()];
01868 Int_t *col1Vyy=new Int_t[GetNy()];
01869 Double_t *dataVyy=new Double_t[GetNy()];
01870 Int_t nError=0;
01871 for (Int_t iy = 0; iy < GetNy(); iy++) {
01872 Double_t dy = input->GetBinError(iy + 1);
01873 rowColVyy[iy] = iy;
01874 col1Vyy[iy] = 0;
01875 if (dy <= 0.0) {
01876 nError++;
01877 if(oneOverZeroError>0.0) {
01878 dataVyy[iy] = 1./ ( oneOverZeroError*oneOverZeroError);
01879 } else {
01880 dataVyy[iy] = 0.0;
01881 }
01882 } else {
01883 dataVyy[iy] = dy * dy;
01884 }
01885 if( dataVyy[iy]>0.0) fNdf ++;
01886 }
01887 DeleteMatrix(&fVyy);
01888 fVyy = CreateSparseMatrix
01889 (GetNy(),GetNy(),GetNy(),rowColVyy,rowColVyy,dataVyy);
01890
01891 TMatrixDSparse *vecV=CreateSparseMatrix
01892 (GetNy(),1,GetNy(),rowColVyy,col1Vyy, dataVyy);
01893 delete[] rowColVyy;
01894 delete[] col1Vyy;
01895 delete[] dataVyy;
01896
01897
01898
01899 DeleteMatrix(&fY);
01900 fY = new TMatrixD(GetNy(), 1);
01901 for (Int_t i = 0; i < GetNy(); i++) {
01902 (*fY) (i, 0) = input->GetBinContent(i + 1);
01903 }
01904
01905 TMatrixDSparse *mAtV=MultiplyMSparseTranspMSparse(fA,vecV);
01906 DeleteMatrix(&vecV);
01907 Int_t nError2=0;
01908 for (Int_t i = 0; i <mAtV->GetNrows();i++) {
01909 if(mAtV->GetRowIndexArray()[i]==
01910 mAtV->GetRowIndexArray()[i+1]) {
01911 nError2 ++;
01912 }
01913 }
01914 if(nError>0) {
01915 if(oneOverZeroError !=0.0) {
01916 if(nError>1) {
01917 Warning("SetInput","%d input bins have zero error,"
01918 " 1/error set to %lf.",nError,oneOverZeroError);
01919 } else {
01920 Warning("SetInput","One input bin has zero error,"
01921 " 1/error set to %lf.",oneOverZeroError);
01922 }
01923 } else {
01924 if(nError>1) {
01925 Warning("SetInput","%d input bins have zero error,"
01926 " and are ignored.",nError);
01927 } else {
01928 Warning("SetInput","One input bin has zero error,"
01929 " and is ignored.");
01930 }
01931 }
01932 }
01933 if(nError2>0) {
01934 if(nError2>1) {
01935 Warning("SetInput","%d output bins are not constrained by any data.",
01936 nError2);
01937 } else {
01938 Warning("SetInput","One output bins is not constrained by any data.");
01939 }
01940
01941 if(oneOverZeroError<=0.0) {
01942 const Int_t *a_rows=fA->GetRowIndexArray();
01943 const Int_t *a_cols=fA->GetColIndexArray();
01944 for (Int_t col = 0; col <mAtV->GetNrows();col++) {
01945 if(mAtV->GetRowIndexArray()[col]==
01946 mAtV->GetRowIndexArray()[col+1]) {
01947 TString binlist("output bin ");
01948 binlist += fXToHist[col];
01949 binlist +=" depends on ignored input bins ";
01950 for(Int_t row=0;row<fA->GetNrows();row++) {
01951 if(input->GetBinError(row + 1)>0.0) continue;
01952 for(Int_t i=a_rows[row];i<a_rows[row+1];i++) {
01953 if(a_cols[i]!=col) continue;
01954 binlist +=" ";
01955 binlist +=row;
01956 }
01957 }
01958 Warning("SetInput","%s",binlist.Data());
01959 }
01960 }
01961 }
01962 }
01963 DeleteMatrix(&mAtV);
01964
01965 return nError+10000*nError2;
01966 }
01967
01968 Double_t TUnfold::DoUnfold(Double_t tau) {
01969
01970
01971
01972
01973
01974
01975
01976
01977
01978
01979
01980 fTauSquared=tau*tau;
01981 return DoUnfold();
01982 }
01983
01984 Int_t TUnfold::ScanLcurve(Int_t nPoint,
01985 Double_t tauMin,Double_t tauMax,
01986 TGraph **lCurve,TSpline **logTauX,
01987 TSpline **logTauY) {
01988
01989
01990
01991
01992
01993
01994
01995
01996
01997 typedef std::map<Double_t,std::pair<Double_t,Double_t> > XYtau_t;
01998 XYtau_t curve;
01999
02000
02001
02002
02003
02004
02005
02006
02007
02008
02009
02010
02011
02012
02013
02014
02015
02016
02017
02018
02019
02020
02021
02022 if((tauMin<=0)||(tauMax<=0.0)||(tauMin>=tauMax)) {
02023
02024 if(GetNdf()<=0) {
02025 Error("ScanLcurve","too few input bins, NDF<=0");
02026 }
02027
02028
02029
02030
02031
02032 DoUnfold(0.0);
02033 Double_t x0=GetLcurveX();
02034 Double_t y0=GetLcurveY();
02035 Info("ScanLcurve","logtau=-Infinity X=%lf Y=%lf",x0,y0);
02036 {
02037
02038 Double_t logTau=
02039 0.5*(TMath::Log10(fChi2A+3.*TMath::Sqrt(GetNdf()+1.0))
02040 -GetLcurveY());
02041 DoUnfold(TMath::Power(10.,logTau));
02042 curve[logTau]=std::make_pair(GetLcurveX(),GetLcurveY());
02043 Info("ScanLcurve","logtau=%lf X=%lf Y=%lf",
02044 logTau,GetLcurveX(),GetLcurveY());
02045 }
02046 if((*curve.begin()).second.first<x0) {
02047
02048
02049
02050
02051
02052
02053 do {
02054 x0=GetLcurveX();
02055 Double_t logTau=(*curve.begin()).first-0.5;
02056 DoUnfold(TMath::Power(10.,logTau));
02057 curve[logTau]=std::make_pair(GetLcurveX(),GetLcurveY());
02058 Info("ScanLcurve","logtau=%lf X=%lf Y=%lf",
02059 logTau,GetLcurveX(),GetLcurveY());
02060 }
02061 while(((int)curve.size()<(nPoint-1)/2)&&
02062 ((*curve.begin()).second.first<x0));
02063 } else {
02064
02065
02066
02067
02068
02069 while(((int)curve.size()<nPoint-1)&&
02070 (((*curve.begin()).second.first-x0>0.00432)||
02071 ((*curve.begin()).second.second-y0>0.00432)||
02072 (curve.size()<2))) {
02073 Double_t logTau=(*curve.begin()).first-0.5;
02074 DoUnfold(TMath::Power(10.,logTau));
02075 curve[logTau]=std::make_pair(GetLcurveX(),GetLcurveY());
02076 Info("ScanLcurve","logtau=%lf X=%lf Y=%lf",
02077 logTau,GetLcurveX(),GetLcurveY());
02078 }
02079 }
02080 } else {
02081 Double_t logTauMin=TMath::Log10(tauMin);
02082 Double_t logTauMax=TMath::Log10(tauMax);
02083 if(nPoint>1) {
02084
02085 DoUnfold(TMath::Power(10.,logTauMax));
02086 Info("ScanLcurve","logtau=%lf X=%lf Y=%lf",
02087 logTauMax,GetLcurveX(),GetLcurveY());
02088 curve[logTauMax]=std::make_pair(GetLcurveX(),GetLcurveY());
02089 }
02090
02091 DoUnfold(TMath::Power(10.,logTauMin));
02092 Info("ScanLcurve","logtau=%lf X=%lf Y=%lf",
02093 logTauMin,GetLcurveX(),GetLcurveY());
02094 curve[logTauMin]=std::make_pair(GetLcurveX(),GetLcurveY());
02095 }
02096
02097
02098
02099
02100
02101
02102 while(int(curve.size())<nPoint-1) {
02103
02104
02105 XYtau_t::const_iterator i0,i1;
02106 i0=curve.begin();
02107 i1=i0;
02108 Double_t logTau=(*i0).first;
02109 Double_t distMax=0.0;
02110 for(i1++;i1!=curve.end();i1++) {
02111 const std::pair<Double_t,Double_t> &xy0=(*i0).second;
02112 const std::pair<Double_t,Double_t> &xy1=(*i1).second;
02113 Double_t dx=xy1.first-xy0.first;
02114 Double_t dy=xy1.second-xy0.second;
02115 Double_t d=TMath::Sqrt(dx*dx+dy*dy);
02116 if(d>=distMax) {
02117 distMax=d;
02118 logTau=0.5*((*i0).first+(*i1).first);
02119 }
02120 i0=i1;
02121 }
02122 DoUnfold(TMath::Power(10.,logTau));
02123 Info("ScanLcurve","logtau=%lf X=%lf Y=%lf",logTau,GetLcurveX(),GetLcurveY());
02124 curve[logTau]=std::make_pair(GetLcurveX(),GetLcurveY());
02125 }
02126
02127
02128
02129
02130
02131
02132 XYtau_t::const_iterator i0,i1;
02133 i0=curve.begin();
02134 i1=i0;
02135 i1++;
02136 Double_t logTauFin=(*i0).first;
02137 if( ((int)curve.size())<nPoint) {
02138
02139 Double_t *cTi=new Double_t[curve.size()-1];
02140 Double_t *cCi=new Double_t[curve.size()-1];
02141 Int_t n=0;
02142 {
02143 Double_t *lXi=new Double_t[curve.size()];
02144 Double_t *lYi=new Double_t[curve.size()];
02145 Double_t *lTi=new Double_t[curve.size()];
02146 for( XYtau_t::const_iterator i=curve.begin();i!=curve.end();i++) {
02147 lXi[n]=(*i).second.first;
02148 lYi[n]=(*i).second.second;
02149 lTi[n]=(*i).first;
02150 n++;
02151 }
02152 TSpline3 *splineX=new TSpline3("x vs tau",lTi,lXi,n);
02153 TSpline3 *splineY=new TSpline3("y vs tau",lTi,lYi,n);
02154
02155
02156 for(Int_t i=0;i<n-1;i++) {
02157 Double_t ltau,xy,bi,ci,di;
02158 splineX->GetCoeff(i,ltau,xy,bi,ci,di);
02159 Double_t tauBar=0.5*(lTi[i]+lTi[i+1]);
02160 Double_t dTau=0.5*(lTi[i+1]-lTi[i]);
02161 Double_t dx1=bi+dTau*(2.*ci+3.*di*dTau);
02162 Double_t dx2=2.*ci+6.*di*dTau;
02163 splineY->GetCoeff(i,ltau,xy,bi,ci,di);
02164 Double_t dy1=bi+dTau*(2.*ci+3.*di*dTau);
02165 Double_t dy2=2.*ci+6.*di*dTau;
02166 cTi[i]=tauBar;
02167 cCi[i]=(dy2*dx1-dy1*dx2)/TMath::Power(dx1*dx1+dy1*dy1,1.5);
02168 }
02169 delete splineX;
02170 delete splineY;
02171 delete[] lXi;
02172 delete[] lYi;
02173 delete[] lTi;
02174 }
02175
02176 TSpline3 *splineC=new TSpline3("L curve curvature",cTi,cCi,n-1);
02177
02178
02179
02180
02181
02182 Int_t iskip=0;
02183 if(n>4) iskip=1;
02184 if(n>7) iskip=2;
02185 Double_t cCmax=cCi[iskip];
02186 Double_t cTmax=cTi[iskip];
02187 for(Int_t i=iskip;i<n-2-iskip;i++) {
02188
02189
02190 Double_t xMax=cTi[i+1];
02191 Double_t yMax=cCi[i+1];
02192 if(cCi[i]>yMax) {
02193 yMax=cCi[i];
02194 xMax=cTi[i];
02195 }
02196
02197
02198
02199 Double_t x,y,b,c,d;
02200 splineC->GetCoeff(i,x,y,b,c,d);
02201
02202 Double_t m_p_half=-c/(3.*d);
02203 Double_t q=b/(3.*d);
02204 Double_t discr=m_p_half*m_p_half-q;
02205 if(discr>=0.0) {
02206
02207 discr=TMath::Sqrt(discr);
02208 Double_t xx;
02209 if(m_p_half>0.0) {
02210 xx = m_p_half + discr;
02211 } else {
02212 xx = m_p_half - discr;
02213 }
02214 Double_t dx=cTi[i+1]-x;
02215
02216 if((xx>0.0)&&(xx<dx)) {
02217 y=splineC->Eval(x+xx);
02218 if(y>yMax) {
02219 yMax=y;
02220 xMax=x+xx;
02221 }
02222 }
02223
02224 if(xx !=0.0) {
02225 xx= q/xx;
02226 } else {
02227 xx=0.0;
02228 }
02229
02230 if((xx>0.0)&&(xx<dx)) {
02231 y=splineC->Eval(x+xx);
02232 if(y>yMax) {
02233 yMax=y;
02234 xMax=x+xx;
02235 }
02236 }
02237 }
02238
02239 if(yMax>cCmax) {
02240 cCmax=yMax;
02241 cTmax=xMax;
02242 }
02243 }
02244 #ifdef DEBUG_LCURVE
02245 {
02246 TCanvas lcc;
02247 lcc.Divide(1,1);
02248 lcc.cd(1);
02249 splineC->Draw();
02250 lcc.SaveAs("splinec.ps");
02251 }
02252 #endif
02253 delete splineC;
02254 delete[] cTi;
02255 delete[] cCi;
02256 logTauFin=cTmax;
02257 DoUnfold(TMath::Power(10.,logTauFin));
02258 Info("ScanLcurve","Result logtau=%lf X=%lf Y=%lf",
02259 logTauFin,GetLcurveX(),GetLcurveY());
02260 curve[logTauFin]=std::make_pair(GetLcurveX(),GetLcurveY());
02261 }
02262
02263
02264
02265
02266
02267
02268 Int_t bestChoice=-1;
02269 if(curve.size()>0) {
02270 Double_t *x=new Double_t[curve.size()];
02271 Double_t *y=new Double_t[curve.size()];
02272 Double_t *logT=new Double_t[curve.size()];
02273 int n=0;
02274 for( XYtau_t::const_iterator i=curve.begin();i!=curve.end();i++) {
02275 if(logTauFin==(*i).first) {
02276 bestChoice=n;
02277 }
02278 x[n]=(*i).second.first;
02279 y[n]=(*i).second.second;
02280 logT[n]=(*i).first;
02281 n++;
02282 }
02283 if(lCurve) {
02284 (*lCurve)=new TGraph(n,x,y);
02285 (*lCurve)->SetTitle("L curve");
02286 }
02287 if(logTauX) (*logTauX)=new TSpline3("log(chi**2)%log(tau)",logT,x,n);
02288 if(logTauY) (*logTauY)=new TSpline3("log(reg.cond)%log(tau)",logT,y,n);
02289 delete[] x;
02290 delete[] y;
02291 delete[] logT;
02292 }
02293
02294 return bestChoice;
02295 }
02296
02297
02298 TH1D *TUnfold::GetOutput(const char *name,const char *title,
02299 Double_t xmin, Double_t xmax) const
02300 {
02301
02302
02303
02304
02305
02306
02307 int nbins = fHistToX.GetSize() - 2;
02308 if (xmin >= xmax) {
02309 xmin = 0.0;
02310 xmax = nbins;
02311 }
02312 TH1D *out = new TH1D(name, title, nbins, xmin, xmax);
02313 GetOutput(out);
02314
02315 return out;
02316 }
02317
02318 TH1D *TUnfold::GetBias(const char *name,const char *title,
02319 Double_t xmin, Double_t xmax) const
02320 {
02321
02322
02323
02324
02325
02326
02327 int nbins = fHistToX.GetSize() - 2;
02328 if (xmin >= xmax) {
02329 xmin = 0.0;
02330 xmax = nbins;
02331 }
02332 TH1D *out = new TH1D(name, title, nbins, xmin, xmax);
02333 for (Int_t i = 0; i < GetNx(); i++) {
02334 out->SetBinContent(fXToHist[i], (*fX0) (i, 0));
02335 }
02336 return out;
02337 }
02338
02339 TH1D *TUnfold::GetFoldedOutput(const char *name,const char *title,
02340 Double_t y0, Double_t y1) const
02341 {
02342
02343
02344
02345
02346
02347
02348 if (y0 >= y1) {
02349 y0 = 0.0;
02350 y1 = GetNy();
02351 }
02352 TH1D *out = new TH1D(name, title, GetNy(), y0, y1);
02353
02354 TMatrixDSparse *AVxx=MultiplyMSparseMSparse(fA,fVxx);
02355
02356
02357 const Int_t *rows_A=fA->GetRowIndexArray();
02358 const Int_t *cols_A=fA->GetColIndexArray();
02359 const Double_t *data_A=fA->GetMatrixArray();
02360 const Int_t *rows_AVxx=AVxx->GetRowIndexArray();
02361 const Int_t *cols_AVxx=AVxx->GetColIndexArray();
02362 const Double_t *data_AVxx=AVxx->GetMatrixArray();
02363
02364 for (Int_t i = 0; i < GetNy(); i++) {
02365 out->SetBinContent(i + 1, (*fAx) (i, 0));
02366 Double_t e2=0.0;
02367 Int_t index_a=rows_A[i];
02368 Int_t index_av=rows_AVxx[i];
02369 while((index_a<rows_A[i+1])&&(index_av<rows_AVxx[i])) {
02370 Int_t j_a=cols_A[index_a];
02371 Int_t j_av=cols_AVxx[index_av];
02372 if(j_a<j_av) {
02373 index_a++;
02374 } else if(j_a>j_av) {
02375 index_av++;
02376 } else {
02377 e2 += data_AVxx[index_av] * data_A[index_a];
02378 index_a++;
02379 index_av++;
02380 }
02381 }
02382 out->SetBinError(i + 1,TMath::Sqrt(e2));
02383 }
02384 DeleteMatrix(&AVxx);
02385 return out;
02386 }
02387
02388 TH1D *TUnfold::GetInput(const char *name,const char *title,
02389 Double_t y0, Double_t y1) const
02390 {
02391
02392
02393
02394
02395
02396
02397 if (y0 >= y1) {
02398 y0 = 0.0;
02399 y1 = GetNy();
02400 }
02401 TH1D *out = new TH1D(name, title, GetNy(), y0, y1);
02402
02403 const Int_t *rows_Vyy=fVyy->GetRowIndexArray();
02404 const Int_t *cols_Vyy=fVyy->GetColIndexArray();
02405 const Double_t *data_Vyy=fVyy->GetMatrixArray();
02406
02407
02408 for (Int_t i = 0; i < GetNy(); i++) {
02409 out->SetBinContent(i + 1, (*fY) (i, 0));
02410 Double_t e=0.0;
02411 for(int index=rows_Vyy[i];index<rows_Vyy[i+1];index++) {
02412 if(cols_Vyy[index]==i) {
02413 e=TMath::Sqrt(data_Vyy[index]);
02414 }
02415 }
02416 out->SetBinError(i + 1, e);
02417 }
02418 return out;
02419 }
02420
02421 TH2D *TUnfold::GetRhoIJ(const char *name,const char *title,
02422 Double_t xmin, Double_t xmax) const
02423 {
02424
02425
02426
02427
02428
02429
02430 int nbins = fHistToX.GetSize() - 2;
02431 if (xmin >= xmax) {
02432 xmin = 0.0;
02433 xmax = nbins;
02434 }
02435 TH2D *out = new TH2D(name, title, nbins, xmin, xmax, nbins, xmin, xmax);
02436 GetRhoIJ(out);
02437 return out;
02438 }
02439
02440 TH2D *TUnfold::GetEmatrix(const char *name,const char *title,
02441 Double_t xmin, Double_t xmax) const
02442 {
02443
02444
02445
02446
02447
02448
02449 int nbins = fHistToX.GetSize() - 2;
02450 if (xmin >= xmax) {
02451 xmin = 0.0;
02452 xmax = nbins;
02453 }
02454 TH2D *out = new TH2D(name, title, nbins, xmin, xmax, nbins, xmin, xmax);
02455 GetEmatrix(out);
02456
02457 return out;
02458 }
02459
02460 TH1D *TUnfold::GetRhoI(const char *name,const char *title,
02461 Double_t xmin, Double_t xmax) const
02462 {
02463
02464
02465
02466
02467
02468
02469 int nbins = fHistToX.GetSize() - 2;
02470 if (xmin >= xmax) {
02471 xmin = 0.0;
02472 xmax = nbins;
02473 }
02474 TH1D *out = new TH1D(name, title, nbins, xmin, xmax);
02475 GetRhoI(out);
02476
02477 return out;
02478 }
02479
02480 TH2D *TUnfold::GetLsquared(const char *name,const char *title,
02481 Double_t xmin, Double_t xmax) const
02482 {
02483
02484
02485
02486
02487
02488
02489 int nbins = fHistToX.GetSize() - 2;
02490 if (xmin >= xmax) {
02491 xmin = 0.0;
02492 xmax = nbins;
02493 }
02494 TH2D *out = new TH2D(name, title, nbins, xmin, xmax, nbins, xmin, xmax);
02495 out->SetOption("BOX");
02496
02497 const Int_t *rows=fLsquared->GetRowIndexArray();
02498 const Int_t *cols=fLsquared->GetColIndexArray();
02499 const Double_t *data=fLsquared->GetMatrixArray();
02500 for (Int_t i = 0; i < GetNx(); i++) {
02501 for (Int_t cindex = rows[i]; cindex < rows[i+1]; cindex++) {
02502 Int_t j=cols[cindex];
02503 out->SetBinContent(fXToHist[i], fXToHist[j], fTauSquared * data[cindex]);
02504 }
02505 }
02506
02507 return out;
02508 }
02509
02510 void TUnfold::SetConstraint(TUnfold::EConstraint constraint) {
02511
02512 if(fConstraint !=constraint) ClearResults();
02513 fConstraint=constraint;
02514 Info("TUnfold::SetConstraint","fConstraint=%d",fConstraint);
02515 }
02516
02517 Double_t TUnfold::GetTau(void) const
02518 {
02519
02520 return TMath::Sqrt(fTauSquared);
02521 }
02522
02523 Double_t TUnfold::GetChi2L(void) const
02524 {
02525
02526 return fLXsquared*fTauSquared;
02527 }
02528
02529 Int_t TUnfold::GetNpar(void) const
02530 {
02531
02532 return GetNx();
02533 }
02534
02535 Double_t TUnfold::GetLcurveX(void) const {
02536
02537 return TMath::Log10(fChi2A);
02538 }
02539
02540 Double_t TUnfold::GetLcurveY(void) const {
02541
02542 return TMath::Log10(fLXsquared);
02543 }
02544
02545 void TUnfold::GetOutput(TH1 *output,const Int_t *binMap) const {
02546
02547
02548
02549
02550
02551
02552
02553
02554 Int_t nbin=output->GetNbinsX();
02555 Double_t *c=new Double_t[nbin+2];
02556 Double_t *e2=new Double_t[nbin+2];
02557 for(Int_t i=0;i<nbin+2;i++) {
02558 c[i]=0.0;
02559 e2[i]=0.0;
02560 }
02561
02562 const Int_t *rows_Vxx=fVxx->GetRowIndexArray();
02563 const Int_t *cols_Vxx=fVxx->GetColIndexArray();
02564 const Double_t *data_Vxx=fVxx->GetMatrixArray();
02565
02566 Int_t binMapSize = fHistToX.GetSize();
02567 for(Int_t i=0;i<binMapSize;i++) {
02568 Int_t destBinI=binMap ? binMap[i] : i;
02569 Int_t srcBinI=fHistToX[i];
02570 if((destBinI>=0)&&(destBinI<nbin+2)&&(srcBinI>=0)) {
02571 c[destBinI]+=(*fX)(srcBinI,0);
02572
02573
02574
02575
02576 Int_t j=0;
02577 Int_t index_vxx=rows_Vxx[srcBinI];
02578 while((j<binMapSize)&&(index_vxx<rows_Vxx[srcBinI+1])) {
02579 Int_t destBinJ=binMap ? binMap[j] : j;
02580 if(destBinI!=destBinJ) {
02581
02582 j++;
02583 } else {
02584 Int_t srcBinJ=fHistToX[j];
02585 if(srcBinJ<0) {
02586
02587 j++;
02588 } else {
02589 if(cols_Vxx[index_vxx]<srcBinJ) {
02590
02591 index_vxx++;
02592 } else if(cols_Vxx[index_vxx]>srcBinJ) {
02593
02594 j++;
02595 } else {
02596
02597 e2[destBinI]+= data_Vxx[index_vxx];
02598 j++;
02599 index_vxx++;
02600 }
02601 }
02602 }
02603 }
02604 }
02605 }
02606 for(Int_t i=0;i<nbin+2;i++) {
02607 output->SetBinContent(i,c[i]);
02608 output->SetBinError(i,TMath::Sqrt(e2[i]));
02609 }
02610 delete[] c;
02611 delete[] e2;
02612 }
02613
02614 void TUnfold::ErrorMatrixToHist(TH2 *ematrix,const TMatrixDSparse *emat,
02615 const Int_t *binMap,Bool_t doClear) const {
02616
02617
02618
02619
02620
02621
02622
02623
02624
02625
02626 Int_t nbin=ematrix->GetNbinsX();
02627 if(doClear) {
02628 for(Int_t i=0;i<nbin+2;i++) {
02629 for(Int_t j=0;j<nbin+2;j++) {
02630 ematrix->SetBinContent(i,j,0.0);
02631 ematrix->SetBinError(i,j,0.0);
02632 }
02633 }
02634 }
02635
02636 if(emat) {
02637 const Int_t *rows_emat=emat->GetRowIndexArray();
02638 const Int_t *cols_emat=emat->GetColIndexArray();
02639 const Double_t *data_emat=emat->GetMatrixArray();
02640
02641 Int_t binMapSize = fHistToX.GetSize();
02642 for(Int_t i=0;i<binMapSize;i++) {
02643 Int_t destBinI=binMap ? binMap[i] : i;
02644 Int_t srcBinI=fHistToX[i];
02645 if((destBinI>=0)&&(destBinI<nbin+2)&&(srcBinI>=0)) {
02646
02647
02648
02649
02650 Int_t j=0;
02651 Int_t index_vxx=rows_emat[srcBinI];
02652 while((j<binMapSize)&&(index_vxx<rows_emat[srcBinI+1])) {
02653 Int_t destBinJ=binMap ? binMap[j] : j;
02654 Int_t srcBinJ=fHistToX[j];
02655 if((destBinJ<0)||(destBinJ>=nbin+2)||(srcBinJ<0)) {
02656
02657 j++;
02658 } else {
02659 if(cols_emat[index_vxx]<srcBinJ) {
02660
02661 index_vxx++;
02662 } else if(cols_emat[index_vxx]>srcBinJ) {
02663
02664 j++;
02665 } else {
02666
02667 Double_t e2= ematrix->GetBinContent(destBinI,destBinJ)
02668 + data_emat[index_vxx];
02669 ematrix->SetBinContent(destBinI,destBinJ,e2);
02670 j++;
02671 index_vxx++;
02672 }
02673 }
02674 }
02675 }
02676 }
02677 }
02678 }
02679
02680 void TUnfold::GetEmatrix(TH2 *ematrix,const Int_t *binMap) const {
02681
02682
02683
02684
02685
02686
02687
02688
02689 ErrorMatrixToHist(ematrix,fVxx,binMap,kTRUE);
02690 }
02691
02692 Double_t TUnfold::GetRhoI(TH1 *rhoi,TH2 *ematrixinv,const Int_t *binMap) const {
02693
02694
02695
02696
02697
02698
02699
02700
02701
02702
02703
02704
02705 Int_t nbin=rhoi->GetNbinsX();
02706
02707
02708
02709
02710
02711
02712
02713
02714
02715
02716
02717
02718
02719
02720 Int_t n=0;
02721 Int_t *outputToEmat=new Int_t[nbin+2];
02722 Int_t *ematToOutput=new Int_t[nbin+2];
02723 Int_t *vxxToEmat=new Int_t[GetNx()];
02724 for(Int_t i=0;i<nbin+2;i++) {
02725 outputToEmat[i]=-1;
02726 }
02727 for(Int_t i=0;i<GetNx();i++) {
02728 Int_t matrix_bin=fXToHist[i];
02729 Int_t output_bin=binMap ? binMap[matrix_bin] : matrix_bin;
02730 if((output_bin>=0)&&(output_bin<nbin+2)) {
02731
02732 if(outputToEmat[output_bin]<0) {
02733
02734 outputToEmat[output_bin]=n;
02735 ematToOutput[n]=output_bin;
02736 n++;
02737 }
02738 vxxToEmat[i]=outputToEmat[output_bin];
02739 } else {
02740 vxxToEmat[i]=-1;
02741 }
02742 }
02743 delete[] outputToEmat;
02744
02745
02746
02747
02748
02749
02750
02751
02752
02753
02754 TMatrixD e(n,n);
02755 const Int_t *rows_Vxx=fVxx->GetRowIndexArray();
02756 const Int_t *cols_Vxx=fVxx->GetColIndexArray();
02757 const Double_t *data_Vxx=fVxx->GetMatrixArray();
02758 for(Int_t i=0;i<GetNx();i++) {
02759 Int_t ie=vxxToEmat[i];
02760 if(ie<0) continue;
02761 for(int index_vxx=rows_Vxx[i];index_vxx<rows_Vxx[i+1];index_vxx++) {
02762 Int_t je=vxxToEmat[cols_Vxx[index_vxx]];
02763 if(je<0) continue;
02764 e(ie,je) += data_Vxx[index_vxx];
02765 }
02766 }
02767 delete[] vxxToEmat;
02768
02769
02770 TMatrixD einv(e);
02771 InvertMConditioned(&einv);
02772
02773 Double_t rhoMax=0.0;
02774 for(Int_t i=0;i<n;i++) {
02775 Int_t i_out=ematToOutput[i];
02776 Double_t rho=1.-1./(einv(i,i)*e(i,i));
02777 if(rho>=0.0) rho=TMath::Sqrt(rho);
02778 else rho=-TMath::Sqrt(-rho);
02779 if(rho>rhoMax) {
02780 rhoMax = rho;
02781 }
02782 rhoi->SetBinContent(i_out,rho);
02783 if(ematrixinv) {
02784 for(Int_t j=0;j<n;j++) {
02785 ematrixinv->SetBinContent(i_out,ematToOutput[j],einv(i,j));
02786 }
02787 }
02788 }
02789 delete[] ematToOutput;
02790
02791 return rhoMax;
02792 }
02793
02794 void TUnfold::GetRhoIJ(TH2 *rhoij,const Int_t *binMap) const {
02795
02796
02797
02798
02799
02800
02801
02802
02803 GetEmatrix(rhoij,binMap);
02804 Int_t nbin=rhoij->GetNbinsX();
02805 Double_t *e=new Double_t[nbin+2];
02806 for(Int_t i=0;i<nbin+2;i++) {
02807 e[i]=TMath::Sqrt(rhoij->GetBinContent(i,i));
02808 }
02809 for(Int_t i=0;i<nbin+2;i++) {
02810 for(Int_t j=0;j<nbin+2;j++) {
02811 if((e[i]>0.0)&&(e[j]>0.0)) {
02812 rhoij->SetBinContent(i,j,rhoij->GetBinContent(i,j)/e[i]/e[j]);
02813 } else {
02814 rhoij->SetBinContent(i,j,0.0);
02815 }
02816 }
02817 }
02818 delete[] e;
02819 }