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 #include "Riostream.h"
00135 #include "TH1.h"
00136 #include "TMath.h"
00137 #include "TClass.h"
00138
00139 #include "TFractionFitter.h"
00140
00141 TVirtualFitter *fractionFitter=0;
00142
00143 ClassImp(TFractionFitter)
00144
00145
00146 TFractionFitter::TFractionFitter() :
00147 fFitDone(kFALSE),
00148 fLowLimitX(0), fHighLimitX(0),
00149 fLowLimitY(0), fHighLimitY(0),
00150 fLowLimitZ(0), fHighLimitZ(0),
00151 fData(0), fIntegralData(0),
00152 fPlot(0)
00153 {
00154
00155
00156 fractionFitter = 0;
00157 fIntegralMCs = 0;
00158 fFractions = 0;
00159
00160 fNpfits = 0;
00161 fNDF = 0;
00162 fChisquare = 0;
00163 fNpar = 0;
00164 }
00165
00166
00167 TFractionFitter::TFractionFitter(TH1* data, TObjArray *MCs) :
00168 fFitDone(kFALSE), fChisquare(0), fPlot(0) {
00169
00170
00171
00172
00173
00174
00175
00176
00177 fData = data;
00178
00179 fLowLimitX = 1;
00180 fHighLimitX = fData->GetNbinsX();
00181 if (fData->GetDimension() > 1) {
00182 fLowLimitY = 1;
00183 fHighLimitY = fData->GetNbinsY();
00184 if (fData->GetDimension() > 2) {
00185 fLowLimitZ = 1;
00186 fHighLimitZ = fData->GetNbinsZ();
00187 }
00188 }
00189 fNpar = MCs->GetEntries();
00190 Int_t par;
00191 for (par = 0; par < fNpar; ++par) {
00192 fMCs.Add(MCs->At(par));
00193
00194 TString s = Form("Prediction for MC sample %i",par);
00195 TH1* pred = (TH1*) ((TH1*)MCs->At(par))->Clone(s);
00196 pred->SetTitle(s);
00197 fAji.Add(pred);
00198 }
00199 fIntegralMCs = new Double_t[fNpar];
00200 fFractions = new Double_t[fNpar];
00201
00202 CheckConsistency();
00203 fWeights.Expand(fNpar);
00204
00205 fractionFitter = TVirtualFitter::Fitter(this, fNpar);
00206 fractionFitter->Clear();
00207 fractionFitter->SetObjectFit(this);
00208 fractionFitter->SetFCN(TFractionFitFCN);
00209
00210 Double_t defaultFraction = 1.0/((Double_t)fNpar);
00211 Double_t defaultStep = 0.01;
00212 for (par = 0; par < fNpar; ++par) {
00213 TString name("frac"); name += par;
00214 fractionFitter->SetParameter(par, name.Data(), defaultFraction, defaultStep, 0, 0);
00215 }
00216 }
00217
00218
00219 TFractionFitter::~TFractionFitter() {
00220
00221
00222 delete fractionFitter;
00223 delete[] fIntegralMCs;
00224 delete[] fFractions;
00225 }
00226
00227
00228 void TFractionFitter::SetData(TH1* data) {
00229
00230
00231
00232
00233
00234 fData = data;
00235 fFitDone = kFALSE;
00236 CheckConsistency();
00237 }
00238
00239
00240 void TFractionFitter::SetMC(Int_t parm, TH1* MC) {
00241
00242
00243
00244
00245
00246 CheckParNo(parm);
00247 fMCs.RemoveAt(parm);
00248 fMCs.AddAt(MC,parm);
00249 fFitDone = kFALSE;
00250 CheckConsistency();
00251 }
00252
00253
00254 void TFractionFitter::SetWeight(Int_t parm, TH1* weight) {
00255
00256
00257
00258
00259
00260
00261 CheckParNo(parm);
00262 if (fWeights[parm]) {
00263 fWeights.RemoveAt(parm);
00264 }
00265 if (weight) {
00266 if (weight->GetNbinsX() != fData->GetNbinsX() ||
00267 (fData->GetDimension() > 1 && weight->GetNbinsY() != fData->GetNbinsY()) ||
00268 (fData->GetDimension() > 2 && weight->GetNbinsZ() != fData->GetNbinsZ())) {
00269 Error("SetWeight","Inconsistent weights histogram for source %d", parm);
00270 return;
00271 }
00272 TString ts = "weight hist: "; ts += weight->GetName();
00273 fWeights.AddAt(weight,parm);
00274 }
00275 }
00276
00277
00278 TVirtualFitter* TFractionFitter::GetFitter() const {
00279
00280
00281
00282 return fractionFitter;
00283 }
00284
00285
00286 void TFractionFitter::CheckParNo(Int_t parm) const {
00287
00288
00289
00290 if (parm < 0 || parm > fNpar) {
00291 Error("CheckParNo","Invalid parameter number %d",parm);
00292 }
00293 }
00294
00295
00296 void TFractionFitter::SetRangeX(Int_t low, Int_t high) {
00297
00298
00299
00300
00301
00302
00303
00304
00305 fLowLimitX = (low > 0 ) ? low : 1;
00306 fHighLimitX = ( high > 0 && high <= fData->GetNbinsX()) ? high : fData->GetNbinsX();
00307 CheckConsistency();
00308 }
00309
00310
00311 void TFractionFitter::ReleaseRangeX() {
00312
00313
00314 fLowLimitX = 1;
00315 fHighLimitX = fData->GetNbinsX();
00316 CheckConsistency();
00317 }
00318
00319
00320 void TFractionFitter::SetRangeY(Int_t low, Int_t high) {
00321
00322
00323
00324
00325
00326
00327
00328
00329 if (fData->GetDimension() < 2) {
00330 Error("SetRangeY","Y range cannot be set for 1D histogram");
00331 return;
00332 }
00333
00334 fLowLimitY = (low > 0) ? low : 1;
00335 fHighLimitY = (high > 0 && high <= fData->GetNbinsY()) ? high : fData->GetNbinsY();
00336 CheckConsistency();
00337 }
00338
00339
00340 void TFractionFitter::ReleaseRangeY() {
00341
00342
00343 fLowLimitY = 1;
00344 fHighLimitY = fData->GetNbinsY();
00345 CheckConsistency();
00346 }
00347
00348
00349
00350 void TFractionFitter::SetRangeZ(Int_t low, Int_t high) {
00351
00352
00353
00354
00355
00356
00357
00358
00359 if (fData->GetDimension() < 3) {
00360 Error("SetRangeZ","Z range cannot be set for 1D or 2D histogram");
00361 return;
00362 }
00363
00364
00365 fLowLimitZ = (low > 0) ? low : 1;
00366 fHighLimitZ = (high > 0 && high <= fData->GetNbinsZ()) ? high : fData->GetNbinsZ();
00367 CheckConsistency();
00368 }
00369
00370
00371 void TFractionFitter::ReleaseRangeZ() {
00372
00373
00374 fLowLimitZ = 1;
00375 fHighLimitZ = fData->GetNbinsZ();
00376 CheckConsistency();
00377 }
00378
00379
00380 void TFractionFitter::Constrain(Int_t parm, Double_t low, Double_t high) {
00381
00382
00383
00384
00385 CheckParNo(parm);
00386 Double_t plist[3];
00387 plist[0] = (Double_t) parm;
00388 plist[1] = low;
00389 plist[2] = high;
00390 fractionFitter->ExecuteCommand("SET LIMIT", plist, 3);
00391 }
00392
00393
00394 void TFractionFitter::UnConstrain(Int_t parm) {
00395
00396
00397 CheckParNo(parm);
00398 Double_t plist[3];
00399 plist[0] = (Double_t) parm;
00400 plist[1] = 0;
00401 plist[2] = 0;
00402 fractionFitter->ExecuteCommand("SET LIMIT", plist, 3);
00403 }
00404
00405
00406 void TFractionFitter::CheckConsistency() {
00407
00408
00409
00410
00411
00412
00413 if (! fData) {
00414 Error("CheckConsistency","Nonexistent data histogram");
00415 return;
00416 }
00417 Int_t minX, maxX, minY, maxY, minZ, maxZ;
00418 Int_t x,y,z,par;
00419 GetRanges(minX, maxX, minY, maxY, minZ, maxZ);
00420 fIntegralData = 0;
00421 fNpfits = 0;
00422 for (z = minZ; z <= maxZ; ++z) {
00423 for (y = minY; y <= maxY; ++y) {
00424 for (x = minX; x <= maxX; ++x) {
00425 fNpfits++;
00426 fIntegralData += fData->GetBinContent(x, y, z);
00427 }
00428 }
00429 }
00430 if (fIntegralData <= 0) {
00431 Error("CheckConsistency","Empty data histogram");
00432 return;
00433 }
00434 TClass* cl = fData->Class();
00435
00436 fNDF = fNpfits - fNpar;
00437
00438 if (fNpar < 2) {
00439 Error("CheckConsistency","Need at least two MC histograms");
00440 return;
00441 }
00442
00443 for (par = 0; par < fNpar; ++par) {
00444 TH1 *h = (TH1*)fMCs.At(par);
00445 if (! h) {
00446 Error("CheckConsistency","Nonexistent MC histogram for source #%d",par);
00447 return;
00448 }
00449 if ((! h->Class()->InheritsFrom(cl)) || h->GetNbinsX() != fData->GetNbinsX() ||
00450 (fData->GetDimension() > 1 && h->GetNbinsY() != fData->GetNbinsY()) ||
00451 (fData->GetDimension() > 2 && h->GetNbinsZ() != fData->GetNbinsZ())) {
00452 Error("CheckConsistency","Histogram inconsistency for source #%d",par);
00453 return;
00454 }
00455 fIntegralMCs[par] = 0;
00456 for (z = minZ; z <= maxZ; ++z) {
00457 for (y = minY; y <= maxY; ++y) {
00458 for (x = minX; x <= maxX; ++x) {
00459 fIntegralMCs[par] += h->GetBinContent(x, y, z);
00460 }
00461 }
00462 }
00463 if (fIntegralMCs[par] <= 0) {
00464 Error("CheckConsistency","Empty MC histogram #%d",par);
00465 }
00466 }
00467 }
00468
00469
00470 Int_t TFractionFitter::Fit() {
00471
00472
00473
00474 Double_t plist[1];
00475 plist[0] = 0.5;
00476
00477 fractionFitter->ExecuteCommand("SET ERRDEF",plist,1);
00478
00479 if (fPlot) {
00480 delete fPlot; fPlot = 0;
00481 }
00482
00483
00484 fractionFitter->SetObjectFit(this);
00485
00486
00487 Int_t status = fractionFitter->ExecuteCommand("MINIMIZE",0,0);
00488 if (status == 0) fFitDone = kTRUE;
00489
00490
00491 ComputeChisquareLambda();
00492
00493 return status;
00494 }
00495
00496
00497 void TFractionFitter::ErrorAnalysis(Double_t UP) {
00498
00499
00500 if (! fFitDone) {
00501 Error("ErrorAnalysis","Fit not yet performed");
00502 return;
00503 }
00504
00505
00506 fractionFitter->SetObjectFit(this);
00507
00508 Double_t plist[1];
00509 plist[0] = UP > 0 ? UP : 0.5;
00510 fractionFitter->ExecuteCommand("SET ERRDEF",plist,1);
00511 Int_t status = fractionFitter->ExecuteCommand("MINOS",0,0);
00512 if (status != 0) {
00513 Error("ErrorAnalysis","Error return from MINOS: %d",status);
00514 }
00515 }
00516
00517
00518 void TFractionFitter::GetResult(Int_t parm, Double_t& value, Double_t& error) const {
00519
00520
00521
00522 CheckParNo(parm);
00523 if (! fFitDone) {
00524 Error("GetResult","Fit not yet performed");
00525 return;
00526 }
00527 char parname[100];
00528 Double_t vlow, vhigh;
00529
00530 fractionFitter->GetParameter(parm, parname, value, error, vlow, vhigh);
00531 }
00532
00533
00534 TH1* TFractionFitter::GetPlot() {
00535
00536
00537
00538
00539
00540
00541 if (! fFitDone) {
00542 Error("GetPlot","Fit not yet performed");
00543 return 0;
00544 }
00545 if (! fPlot) {
00546 Double_t plist[1];
00547 plist[0] = 3;
00548 fractionFitter->ExecuteCommand("CALL FCN", plist, 1);
00549 }
00550 return fPlot;
00551 }
00552
00553
00554 void TFractionFitter::GetRanges(Int_t& minX, Int_t& maxX, Int_t& minY, Int_t& maxY,
00555 Int_t& minZ, Int_t& maxZ) const {
00556
00557
00558
00559 if (fData->GetDimension() < 2) {
00560 minY = maxY = minZ = maxZ = 0;
00561 minX = fLowLimitX;
00562 maxX = fHighLimitX;
00563 } else if (fData->GetDimension() < 3) {
00564 minZ = maxZ = 0;
00565 minX = fLowLimitX;
00566 maxX = fHighLimitX;
00567 minY = fLowLimitY;
00568 maxY = fHighLimitY;
00569 } else {
00570 minX = fLowLimitX;
00571 maxX = fHighLimitX;
00572 minY = fLowLimitY;
00573 maxY = fHighLimitY;
00574 minZ = fLowLimitZ;
00575 maxZ = fHighLimitZ;
00576 }
00577 }
00578
00579
00580 void TFractionFitter::ComputeFCN(Int_t& , Double_t* ,
00581 Double_t& f, Double_t* xx, Int_t flag)
00582 {
00583
00584
00585
00586 Int_t bin, mc;
00587 Int_t minX, maxX, minY, maxY, minZ, maxZ;
00588 Int_t x,y,z;
00589 GetRanges(minX, maxX, minY, maxY, minZ, maxZ);
00590 for (mc = 0; mc < fNpar; ++mc) {
00591 Double_t tot;
00592 TH1 *h = (TH1*)fMCs[mc];
00593 TH1 *hw = (TH1*)fWeights[mc];
00594 if (hw) {
00595 tot = 0;
00596 for (z = minZ; z <= maxZ; ++z) {
00597 for (y = minY; y <= maxY; ++y) {
00598 for (x = minX; x <= maxX; ++x) {
00599 Double_t weight = hw->GetBinContent(x, y, z);
00600 if (weight <= 0) {
00601 Error("ComputeFCN","Invalid weight encountered for MC source %d",mc);
00602 return;
00603 }
00604 tot += weight * h->GetBinContent(x, y, z);
00605 }
00606 }
00607 }
00608 } else tot = fIntegralMCs[mc];
00609 fFractions[mc] = xx[mc] * fIntegralData / tot;
00610 }
00611
00612 if (flag == 3) {
00613 TString ts = "Fraction fit to hist: "; ts += fData->GetName();
00614 fPlot = (TH1*) fData->Clone(ts.Data());
00615 fPlot->Reset();
00616 }
00617
00618 Double_t result = 0;
00619 for (z = minZ; z <= maxZ; ++z) {
00620 for (y = minY; y <= maxY; ++y) {
00621 for (x = minX; x <= maxX; ++x) {
00622 bin = fData->GetBin(x, y, z);
00623
00624
00625 int k0;
00626 Double_t ti; Double_t aki;
00627 FindPrediction(bin, fFractions, ti, k0, aki);
00628
00629 Double_t prediction = 0;
00630 for (mc = 0; mc < fNpar; ++mc) {
00631 TH1 *h = (TH1*)fMCs[mc];
00632 TH1 *hw = (TH1*)fWeights[mc];
00633 Double_t binPrediction;
00634 Double_t binContent = h->GetBinContent(bin);
00635 Double_t weight = hw ? hw->GetBinContent(bin) : 1;
00636 if (k0 >= 0 && fFractions[mc] == fFractions[k0]) {
00637 binPrediction = aki;
00638 } else {
00639 binPrediction = binContent > 0 ? binContent / (1+weight*fFractions[mc]*ti) : 0;
00640 }
00641
00642 prediction += fFractions[mc]*weight*binPrediction;
00643 result -= binPrediction;
00644 if (binContent > 0 && binPrediction > 0)
00645 result += binContent*TMath::Log(binPrediction);
00646
00647 if (flag == 3) {
00648 ((TH1*)fAji.At(mc))->SetBinContent(bin, binPrediction);
00649 }
00650 }
00651
00652 if (flag == 3) {
00653 fPlot->SetBinContent(bin, prediction);
00654 }
00655
00656 result -= prediction;
00657 Double_t found = fData->GetBinContent(bin);
00658 if (found > 0 && prediction > 0)
00659 result += found*TMath::Log(prediction);
00660 }
00661 }
00662 }
00663
00664 f = -result;
00665 }
00666
00667
00668 void TFractionFitter::FindPrediction(int bin, Double_t *fractions, Double_t &ti, int& k0, Double_t &aki) const {
00669
00670
00671
00672 Int_t par;
00673 TH1 *hw;
00674 for (par = 0; par < fNpar; ++par) {
00675 hw = (TH1*)fWeights.At(par);
00676 Double_t weightedFraction = hw ?
00677 hw->GetBinContent(bin) * fractions[par] : fractions[par];
00678 if (weightedFraction == 0) {
00679 Error("FindPrediction","Fraction[%d] = 0!",par);
00680 return;
00681 }
00682 }
00683
00684
00685 if (TMath::Nint(fData->GetBinContent(bin)) == 0) {
00686 ti = 1;
00687 k0 = -1;
00688 aki=0;
00689 return;
00690 }
00691
00692
00693 k0 = 0;
00694 TH1 *hw0 = (TH1*)fWeights.At(0);
00695 Double_t refWeightedFraction = hw0 ?
00696 hw0->GetBinContent(bin) * fractions[0] : fractions[0];
00697 for (par = 1; par < fNpar; ++par) {
00698 hw = (TH1*)fWeights.At(par);
00699 Double_t weightedFraction = hw ?
00700 hw->GetBinContent(bin) * fractions[par] : fractions[par];
00701 if (weightedFraction > refWeightedFraction) {
00702 k0 = par;
00703 refWeightedFraction = weightedFraction;
00704 }
00705 }
00706 int nMax = 1;
00707 Double_t contentsMax = ((TH1*)fMCs.At(k0))->GetBinContent(bin);
00708 for (par = 0; par < fNpar; ++par) {
00709 if (par == k0) continue;
00710 hw = (TH1*)fWeights.At(par);
00711 Double_t weightedFraction = hw ?
00712 hw->GetBinContent(bin) * fractions[par] : fractions[par];
00713 if (weightedFraction == refWeightedFraction) {
00714 nMax++;
00715 contentsMax += ((TH1*)fMCs.At(par))->GetBinContent(bin);
00716 }
00717 }
00718 Double_t tmin = -1/refWeightedFraction;
00719
00720 if (contentsMax == 0) {
00721 aki = fData->GetBinContent(bin)/(1+refWeightedFraction);
00722 for (par = 0; par < fNpar; ++par) {
00723 hw = (TH1*)fWeights.At(par);
00724 if (par != k0) {
00725 Double_t weightedFraction = hw ?
00726 hw->GetBinContent(bin) * fractions[par] : fractions[par];
00727 if (weightedFraction != refWeightedFraction)
00728 aki -= ((TH1*)fMCs.At(par))->GetBinContent(bin)*weightedFraction/ (refWeightedFraction - weightedFraction);
00729 }
00730 }
00731 if (aki > 0) {
00732 if (nMax) aki /= nMax;
00733 ti = tmin;
00734 return;
00735 }
00736 }
00737 k0 = -1;
00738
00739
00740 ti = 0;
00741 for (Double_t step = 0.2;;) {
00742 if (ti >= 1 || ti < tmin) {
00743 step /= 10;
00744 ti = 0;
00745 }
00746 Double_t aFunction = - fData->GetBinContent(bin)/(1-ti);
00747 Double_t aDerivative = aFunction / (1-ti);
00748 for (par = 0; par < fNpar; ++par) {
00749 TH1 *h = (TH1*)fMCs.At(par);
00750 hw = (TH1*)fWeights.At(par);
00751 Double_t weightedFraction = hw ?
00752 hw->GetBinContent(bin) * fractions[par] : fractions[par];
00753 Double_t d = 1/(ti+1/weightedFraction);
00754 aFunction += h->GetBinContent(bin)*d;
00755 aDerivative -= h->GetBinContent(bin)*d*d;
00756 }
00757 if (TMath::Abs(aFunction) < 1e-12) break;
00758 Double_t delta = -aFunction/aDerivative;
00759 if (TMath::Abs(delta) > step)
00760 delta = (delta > 0) ? step : -step;
00761 ti += delta;
00762 if (TMath::Abs(delta) < 1e-13) break;
00763 }
00764
00765 return;
00766 }
00767
00768
00769
00770 void TFractionFitFCN(Int_t& npar, Double_t* gin, Double_t& f, Double_t* par, Int_t flag) {
00771
00772
00773
00774 TFractionFitter* fitter = dynamic_cast<TFractionFitter*>(fractionFitter->GetObjectFit());
00775 if (!fitter) {
00776 Error("TFractionFitFCN","Invalid fit object encountered!");
00777 return;
00778 }
00779 fitter->ComputeFCN(npar, gin, f, par, flag);
00780 }
00781
00782
00783
00784 Double_t TFractionFitter::GetChisquare() const
00785 {
00786
00787
00788
00789
00790
00791
00792
00793
00794
00795
00796
00797
00798
00799
00800
00801
00802
00803 return fChisquare;
00804 }
00805
00806
00807 Int_t TFractionFitter::GetNDF() const
00808 {
00809
00810
00811
00812
00813
00814 if (fNDF == 0) return fNpfits-fNpar;
00815 return fNDF;
00816 }
00817
00818
00819 Double_t TFractionFitter::GetProb() const
00820 {
00821
00822
00823 Int_t ndf = fNpfits - fNpar;
00824 if (ndf <= 0) return 0;
00825 return TMath::Prob(fChisquare,ndf);
00826 }
00827
00828
00829 void TFractionFitter::ComputeChisquareLambda()
00830 {
00831
00832
00833
00834 if ( !fFitDone ) {
00835 Error("ComputeChisquareLambda","Fit not yet (successfully) performed");
00836 fChisquare = 0;
00837 return;
00838 }
00839
00840
00841 if (! fPlot)
00842 GetPlot();
00843
00844 Int_t minX, maxX, minY, maxY, minZ, maxZ;
00845 GetRanges(minX, maxX, minY, maxY, minZ, maxZ);
00846
00847 Double_t logLyn = 0;
00848 Double_t logLmn = 0;
00849 for(Int_t x = minX; x <= maxX; x++) {
00850 for(Int_t y = minY; y <= maxY; y++) {
00851 for(Int_t z = minZ; z <= maxZ; z++) {
00852 Double_t di = fData->GetBinContent(x, y, z);
00853 Double_t fi = fPlot->GetBinContent(x, y, z);
00854 if(fi != 0) logLyn += di * TMath::Log(fi) - fi;
00855 if(di != 0) logLmn += di * TMath::Log(di) - di;
00856 for(Int_t j = 0; j < fNpar; j++) {
00857 Double_t aji = ((TH1*)fMCs.At(j))->GetBinContent(x, y, z);
00858 Double_t bji = ((TH1*)fAji.At(j))->GetBinContent(x, y, z);
00859 if(bji != 0) logLyn += aji * TMath::Log(bji) - bji;
00860 if(aji != 0) logLmn += aji * TMath::Log(aji) - aji;
00861 }
00862 }
00863 }
00864 }
00865
00866 fChisquare = -2*logLyn + 2*logLmn;
00867
00868 return;
00869 }
00870
00871
00872 TH1* TFractionFitter::GetMCPrediction(Int_t parm) const
00873 {
00874
00875
00876
00877
00878 CheckParNo(parm);
00879 if ( !fFitDone ) {
00880 Error("GetMCPrediction","Fit not yet performed");
00881 return 0;
00882 }
00883 return (TH1*) fAji.At(parm);
00884 }