TFractionFitter.cxx

Go to the documentation of this file.
00001 // @(#)root/hist:$Id: TFractionFitter.cxx 35510 2010-09-21 08:47:21Z moneta $
00002 // Author: Frank Filthaut filthaut@hef.kun.nl  20/05/2002
00003 // with additions by Bran Wijngaarden <dwijngaa@hef.kun.nl>
00004 
00005 ///////////////////////////////////////////////////////////////////////////////
00006 //
00007 // Fits MC fractions to data histogram (a la HMCMLL, see R. Barlow and C. Beeston,
00008 // Comp. Phys. Comm. 77 (1993) 219-228, and http://www.hep.man.ac.uk/~roger/hfrac.f).
00009 //
00010 // The virtue of this fit is that it takes into account both data and Monte Carlo
00011 // statistical uncertainties. The way in which this is done is through a standard
00012 // likelihood fit using Poisson statistics; however, the template (MC) predictions
00013 // are also varied within statistics, leading to additional contributions to the
00014 // overall likelihood. This leads to many more fit parameters (one per bin per
00015 // template), but the minimisation with respect to these additional parameters is
00016 // done analytically rather than introducing them as formal fit parameters. Some
00017 // special care needs to be taken in the case of bins with zero content. For more
00018 // details please see the original publication cited above.
00019 //
00020 // An example application of this fit is given below. For a TH1* histogram
00021 // ("data") fitted as the sum of three Monte Carlo sources ("mc"):
00022 //
00023 // {
00024 //   TH1F *data;                              //data histogram
00025 //   TH1F *mc0;                               // first MC histogram
00026 //   TH1F *mc1;                               // second MC histogram
00027 //   TH1F *mc2;                               // third MC histogram
00028 //   ....                                     // retrieve histograms
00029 //   TObjArray *mc = new TObjArray(3);        // MC histograms are put in this array
00030 //   mc->Add(mc0);
00031 //   mc->Add(mc1);
00032 //   mc->Add(mc2);
00033 //   TFractionFitter* fit = new TFractionFitter(data, mc); // initialise
00034 //   fit->Constrain(1,0.0,1.0);               // constrain fraction 1 to be between 0 and 1
00035 //   fit->SetRangeX(1,15);                    // use only the first 15 bins in the fit
00036 //   Int_t status = fit->Fit();               // perform the fit
00037 //   cout << "fit status: " << status << endl;
00038 //   if (status == 0) {                       // check on fit status
00039 //     TH1F* result = (TH1F*) fit->GetPlot();
00040 //     data->Draw("Ep");
00041 //     result->Draw("same");
00042 //   }
00043 // }
00044 //
00045 //
00046 // Assumptions
00047 // ===========
00048 // A few assumptions need to be made for the fit procedure to be carried out:
00049 //
00050 // (1) The total number of events in each template is not too small
00051 //     (so that its Poisson uncertainty can be neglected).
00052 // (2) The number of events in each bin is much smaller than the total
00053 //     number of events in each template (so that multinomial
00054 //     uncertainties can be replaced with Poisson uncertainties).
00055 //
00056 // Biased fit uncertainties may result if these conditions are not fulfilled
00057 // (see e.g. arXiv:0803.2711).
00058 //
00059 // Instantiation
00060 // =============
00061 // A fit object is instantiated through
00062 //     TFractionFitter* fit = new TFractionFitter(data, mc);
00063 // A number of basic checks (intended to ensure that the template histograms
00064 // represent the same "kind" of distribution as the data one) are carried out.
00065 // The TVirtualFitter object is then addressed and all fit parameters (the
00066 // template fractions) declared (initially unbounded).
00067 //
00068 // Applying constraints
00069 // ====================
00070 // Fit parameters can be constrained through
00071 //     fit->Constrain(parameter #, lower bound, upper bound);
00072 // Setting lower bound = upper bound = 0 removes the constraint (a la Minuit);
00073 // however, a function
00074 //     fit->Unconstrain(parameter #)
00075 // is also provided to simplify this.
00076 //
00077 // Setting parameter values
00078 // ========================
00079 // The function
00080 //     TVirtualFitter* vFit = fit->GetFitter();
00081 // is provided for direct access to the TVirtualFitter object. This allows to
00082 // set and fix parameter values, and set step sizes directly.
00083 //
00084 // Restricting the fit range
00085 // =========================
00086 // The fit range can be restricted through
00087 //     fit->SetRangeX(first bin #, last bin #);
00088 // and freed using
00089 //     fit->ReleaseRangeX();
00090 // For 2D histograms the Y range can be similarly restricted using
00091 //     fit->SetRangeY(first bin #, last bin #);
00092 //     fit->ReleaseRangeY();
00093 // and for 3D histograms also
00094 //     fit->SetRangeZ(first bin #, last bin #);
00095 //     fit->ReleaseRangeZ();
00096 //
00097 // Weights histograms
00098 // ==================
00099 // Weights histograms (for a motivation see the above publication) can be specified
00100 // for the individual MC sources through
00101 //     fit->SetWeight(parameter #, pointer to weights histogram);
00102 // and unset by specifying a null pointer.
00103 //
00104 // Obtaining fit results
00105 // =====================
00106 // The fit is carried out through
00107 //     Int_t status = fit->Fit();
00108 // where  status  is the code returned from the "MINIMIZE" command. For fits
00109 // that converged, parameter values and errors can be obtained through
00110 //     fit->GetResult(parameter #, value, error);
00111 // and the histogram corresponding to the total Monte Carlo prediction (which
00112 // is not the same as a simple weighted sum of the input Monte Carlo distributions)
00113 // can be obtained by
00114 //     TH1* result = fit->GetPlot();
00115 //
00116 // Using different histograms
00117 // ==========================
00118 // It is possible to change the histogram being fitted through
00119 //     fit->SetData(TH1* data);
00120 // and to change the template histogram for a given parameter number through
00121 //     fit->SetMC(parameter #, TH1* MC);
00122 // This can speed up code in case of multiple data or template histograms;
00123 // however, it should be done with care as any settings are taken over from
00124 // the previous fit. In addition, neither the dimensionality nor the numbers of
00125 // bins of the histograms should change (in that case it is better to instantiate
00126 // a new TFractionFitter object).
00127 //
00128 // Errors
00129 // ======
00130 // Any serious inconsistency results in an error.
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    // TFractionFitter default constructor.
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    // TFractionFitter constructor. Does a complete initialisation (including
00170    // consistency checks, default fit range as the whole histogram but without
00171    // under- and overflows, and declaration of the fit parameters). Note that
00172    // the histograms are not copied, only references are used.
00173    // Arguments:
00174    //     data: histogram to be fitted
00175    //     MCs:  array of TH1* corresponding template distributions
00176 
00177    fData = data;
00178    // Default: include all of the histogram (but without under- and overflows)
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       // Histogram containing template prediction
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   // TFractionFitter default destructor
00221 
00222    delete fractionFitter;
00223    delete[] fIntegralMCs;
00224    delete[] fFractions;
00225 }
00226 
00227 //______________________________________________________________________________
00228 void TFractionFitter::SetData(TH1* data) {
00229    // Change the histogram to be fitted to. Notes:
00230    // - Parameter constraints and settings are retained from a possible previous fit.
00231    // - Modifying the dimension or number of bins results in an error (in this case
00232    //   rather instantiate a new TFractionFitter object)
00233 
00234    fData = data;
00235    fFitDone = kFALSE;
00236    CheckConsistency();
00237 }
00238 
00239 //______________________________________________________________________________
00240 void TFractionFitter::SetMC(Int_t parm, TH1* MC) {
00241    // Change the histogram for template number <parm>. Notes:
00242    // - Parameter constraints and settings are retained from a possible previous fit.
00243    // - Modifying the dimension or number of bins results in an error (in this case
00244    //   rather instantiate a new TFractionFitter object)
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    // Set bin by bin weights for template number <parm> (the parameter numbering
00256    // follows that of the input template vector).
00257    // Weights can be "unset" by passing a null pointer.
00258    // Consistency of the weights histogram with the data histogram is checked at
00259    // this point, and an error in case of problems.
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    // Give direct access to the underlying minimisation class. This can be
00280    // used e.g. to modify parameter values or step sizes.
00281 
00282    return fractionFitter;
00283 }
00284 
00285 //______________________________________________________________________________
00286 void TFractionFitter::CheckParNo(Int_t parm) const {
00287    // Function for internal use, checking parameter validity
00288    // An invalid parameter results in an error.
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    // Set the X range of the histogram to be used in the fit.
00298    // Use ReleaseRangeX() to go back to fitting the full histogram.
00299    // The consistency check ensures that no empty fit range occurs (and also
00300    // recomputes the bin content integrals).
00301    // Arguments:
00302    //     low:  lower X bin number
00303    //     high: upper X bin number
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    // Release restrictions on the X range of the histogram to be used in the fit.
00313 
00314    fLowLimitX = 1;
00315    fHighLimitX = fData->GetNbinsX();
00316    CheckConsistency();
00317 }
00318 
00319 //______________________________________________________________________________
00320 void TFractionFitter::SetRangeY(Int_t low, Int_t high) {
00321    // Set the Y range of the histogram to be used in the fit (2D or 3D histograms only).
00322    // Use ReleaseRangeY() to go back to fitting the full histogram.
00323    // The consistency check ensures that no empty fit range occurs (and also
00324    // recomputes the bin content integrals).
00325    // Arguments:
00326    //     low:  lower Y bin number
00327    //     high: upper Y bin number
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    // Release restrictions on the Y range of the histogram to be used in the fit.
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    // Set the Z range of the histogram to be used in the fit (3D histograms only).
00352    // Use ReleaseRangeY() to go back to fitting the full histogram.
00353    // The consistency check ensures that no empty fit range occurs (and also
00354    // recomputes the bin content integrals).
00355    // Arguments:
00356    //     low:  lower Y bin number
00357    //     high: upper Y bin number
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    // Release restrictions on the Z range of the histogram to be used in the fit.
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    // Constrain the values of parameter number <parm> (the parameter numbering
00382    // follows that of the input template vector).
00383    // Use UnConstrain() to remove this constraint.
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    // Remove the constraints on the possible values of parameter <parm>.
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    // Function used internally to check the consistency between the
00408    // various histograms. Checks are performed on nonexistent or empty
00409    // histograms, the precise histogram class, and the number of bins.
00410    // In addition, integrals over the "allowed" bin ranges are computed.
00411    // Any inconsistency results in a error.
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    // Perform the fit with the default UP value.
00472    // The value returned is the minimisation status.
00473 
00474    Double_t plist[1];
00475    plist[0] = 0.5;
00476    // set the UP value to 0.5
00477    fractionFitter->ExecuteCommand("SET ERRDEF",plist,1);
00478    // remove any existing output histogram
00479    if (fPlot) {
00480       delete fPlot; fPlot = 0;
00481    }
00482 
00483    // Make sure the correct likelihood computation is used
00484    fractionFitter->SetObjectFit(this);
00485 
00486    // fit
00487    Int_t status = fractionFitter->ExecuteCommand("MINIMIZE",0,0);
00488    if (status == 0) fFitDone = kTRUE;
00489 
00490    // determine goodness of fit
00491    ComputeChisquareLambda();
00492 
00493    return status;
00494 }
00495 
00496 //______________________________________________________________________________
00497 void TFractionFitter::ErrorAnalysis(Double_t UP) {
00498    // Set UP to the given value (see class TMinuit), and perform a MINOS minimisation.
00499 
00500    if (! fFitDone) {
00501       Error("ErrorAnalysis","Fit not yet performed");
00502       return;
00503    }
00504 
00505    // Make sure the correct likelihood computation is used
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    // Obtain the fit result for parameter <parm> (the parameter numbering
00520    // follows that of the input template vector).
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    // Return the "template prediction" corresponding to the fit result (this is not
00536    // the same as the weighted sum of template distributions, as template statistical
00537    // uncertainties are taken into account).
00538    // Note that the name of this histogram will simply be the same as that of the
00539    // "data" histogram, prefixed with the string "Fraction fit to hist: ".
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    // Used internally to obtain the bin ranges according to the dimensionality of
00557    // the histogram and the limits set by hand.
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& /*npar*/, Double_t* /*gin*/,
00581                                  Double_t& f, Double_t* xx, Int_t flag)
00582 {
00583    // Used internally to compute the likelihood value.
00584 
00585    // normalise the fit parameters
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    // likelihood computation
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             // Solve for the "predictions"
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    // Function used internally to obtain the template prediction in the individual bins
00670 
00671    // Sanity check: none of the fractions should be =0
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    // Case data = 0
00685    if (TMath::Nint(fData->GetBinContent(bin)) == 0) {
00686       ti = 1;
00687       k0 = -1;
00688       aki=0;
00689       return;
00690    }
00691 
00692    // Case one or more of the MC bin contents = 0: find largest fraction
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    // Case of nonzero histogram contents: solve for Ti using Newton's method
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    // Function called by the minimisation package. The actual functionality is passed
00772    // on to the TFractionFitter::ComputeFCN member function.
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    // Return the likelihood ratio Chi-squared (chi2) for the fit.
00787    // The value is computed when the fit is executed successfully.
00788    // Chi2 calculation is based on the "likelihood ratio" lambda,
00789    // lambda = L(y;n) / L(m;n),
00790    // where L(y;n) is the likelihood of the fit result <y> describing the data <n>
00791    // and L(m;n) is the likelihood of an unknown "true" underlying distribution
00792    // <m> describing the data <n>. Since <m> is unknown, the data distribution is
00793    // used instead,
00794    // lambda = L(y;n) / L(n;n).
00795    // Note that this ratio is 1 if the fit is perfect. The chi2 value is then
00796    // computed according to
00797    // chi2 = -2*ln(lambda).
00798    // This parameter can be shown to follow a Chi-square distribution. See for
00799    // example S. Baker and R. Cousins, "Clarification of the use of chi-square
00800    // and likelihood functions in fits to histograms", Nucl. Instr. Meth. A221, 
00801    // pp. 437-442 (1984)
00802 
00803    return fChisquare;
00804 }
00805 
00806 //______________________________________________________________________________
00807 Int_t TFractionFitter::GetNDF() const
00808 {
00809    // return the number of degrees of freedom in the fit
00810    // the fNDF parameter has been previously computed during a fit.
00811    // The number of degrees of freedom corresponds to the number of points
00812    // used in the fit minus the number of templates.
00813 
00814    if (fNDF == 0) return fNpfits-fNpar;
00815    return fNDF;
00816 }
00817 
00818 //______________________________________________________________________________
00819 Double_t TFractionFitter::GetProb() const
00820 {
00821    // return the fit probability
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    // Method used internally to compute the likelihood ratio chi2
00832    // See the function GetChisquare() for details
00833 
00834    if ( !fFitDone ) {
00835       Error("ComputeChisquareLambda","Fit not yet (successfully) performed");
00836       fChisquare = 0;
00837       return;
00838    }
00839 
00840    // fPlot must be initialized and filled. Leave this to the GetPlot() method.
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; // likelihood of prediction
00848    Double_t logLmn = 0; // likelihood of data ("true" distribution)
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    // Return the adjusted MC template (Aji) for template (parm).
00875    // Note that the (Aji) times fractions only sum to the total prediction
00876    // of the fit if all weights are 1.
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 }

Generated on Tue Jul 5 14:24:17 2011 for ROOT_528-00b_version by  doxygen 1.5.1