ROOT logo

From $ROOTSYS/tutorials/fit/combinedFit.C

//+ Combined (simultaneous) fit of two histogram with separate functions 
//  and some common parameters
//
// See http://root.cern.ch/phpBB3//viewtopic.php?f=3&t=11740#p50908
// for a modified version working with Fumili or GSLMultiFit 
//
// N.B. this macro must be compiled with ACliC 
//
//Author: L. Moneta - Dec 2010

#include "Fit/Fitter.h"
#include "Fit/BinData.h"
#include "Fit/Chi2FCN.h"
#include "TH1.h"
#include "TList.h"
#include "Math/WrappedMultiTF1.h"
#include "HFitInterface.h"
#include "TCanvas.h"
#include "TStyle.h"


// definition of shared parameter
// background function 
int iparB[2] = { 0,      // exp amplitude in B histo
                 2    // exp common parameter 
};

// signal + background function 
int iparSB[5] = { 1, // exp amplitude in S+B histo
                  2, // exp common parameter
                  3, // gaussian amplitude
                  4, // gaussian mean
                  5  // gaussian sigma
};

struct GlobalChi2 { 
   GlobalChi2(  ROOT::Math::IMultiGenFunction & f1,  
                ROOT::Math::IMultiGenFunction & f2) : 
      fChi2_1(&f1), fChi2_2(&f2) {}

   // parameter vector is first background (in common 1 and 2) 
   // and then is signal (only in 2)
   double operator() (const double *par) const {
      double p1[2];
      for (int i = 0; i < 2; ++i) p1[i] = par[iparB[i] ];

      double p2[5]; 
      for (int i = 0; i < 5; ++i) p2[i] = par[iparSB[i] ];

      return (*fChi2_1)(p1) + (*fChi2_2)(p2);
   } 

   const  ROOT::Math::IMultiGenFunction * fChi2_1;
   const  ROOT::Math::IMultiGenFunction * fChi2_2;
};

void combinedFit() { 

#if defined(__CINT__) && !defined(__MAKECINT__) 
   cout << "ERROR: This tutorial can run only using ACliC, you must run it by doing: " << endl;
   cout << "\t .x $ROOTSYS/tutorials/fit/combinedFit.C+" << endl;
   return;
#endif 

  TH1D * hB = new TH1D("hB","histo B",100,0,100);
  TH1D * hSB = new TH1D("hSB","histo S+B",100, 0,100);

  TF1 * fB = new TF1("fB","expo",0,100);
  fB->SetParameters(1,-0.05);
  hB->FillRandom("fB");

  TF1 * fS = new TF1("fS","gaus",0,100);
  fS->SetParameters(1,30,5);

  hSB->FillRandom("fB",2000);
  hSB->FillRandom("fS",1000);

  // perform now global fit 

  TF1 * fSB = new TF1("fSB","expo + gaus(2)",0,100);

  ROOT::Math::WrappedMultiTF1 wfB(*fB,1);
  ROOT::Math::WrappedMultiTF1 wfSB(*fSB,1);

  ROOT::Fit::DataOptions opt; 
  ROOT::Fit::DataRange rangeB; 
  // set the data range
  rangeB.SetRange(10,90);
  ROOT::Fit::BinData dataB(opt,rangeB); 
  ROOT::Fit::FillData(dataB, hB);

  ROOT::Fit::DataRange rangeSB; 
  rangeSB.SetRange(10,50);
  ROOT::Fit::BinData dataSB(opt,rangeSB); 
  ROOT::Fit::FillData(dataSB, hSB);

  ROOT::Fit::Chi2Function chi2_B(dataB, wfB);
  ROOT::Fit::Chi2Function chi2_SB(dataSB, wfSB);

  GlobalChi2 globalChi2(chi2_B, chi2_SB);

  ROOT::Fit::Fitter fitter;

  const int Npar = 6; 
  double par0[Npar] = { 5,5,-0.1,100, 30,10};

  // create before the parameter settings in order to fix or set range on them
  fitter.Config().SetParamsSettings(6,par0);
  // fix 5-th parameter  
  fitter.Config().ParSettings(4).Fix();
  // set limits on the third and 4-th parameter
  fitter.Config().ParSettings(2).SetLimits(-10,-1.E-4);
  fitter.Config().ParSettings(3).SetLimits(0,10000);
  fitter.Config().ParSettings(3).SetStepSize(5);

  fitter.Config().MinimizerOptions().SetPrintLevel(0);
  fitter.Config().SetMinimizer("Minuit2","Migrad"); 

  // fit FCN function directly 
  // (specify optionally data size and flag to indicate that is a chi2 fit)
  fitter.FitFCN(6,globalChi2,0,dataB.Size()+dataSB.Size(),true);
  ROOT::Fit::FitResult result = fitter.Result();
  result.Print(std::cout);

  TCanvas * c1 = new TCanvas("Simfit","Simultaneous fit of two histograms",
                             10,10,700,700);
  c1->Divide(1,2);
  c1->cd(1);
  gStyle->SetOptFit(1111);

  fB->SetFitResult( result, iparB);
  fB->SetRange(rangeB().first, rangeB().second);   
  fB->SetLineColor(kBlue);
  hB->GetListOfFunctions()->Add(fB);
  hB->Draw(); 

  c1->cd(2);
  fSB->SetFitResult( result, iparSB);
  fSB->SetRange(rangeSB().first, rangeSB().second);   
  fSB->SetLineColor(kRed);
  hSB->GetListOfFunctions()->Add(fSB);
  hSB->Draw(); 
    

}
 combinedFit.C:1
 combinedFit.C:2
 combinedFit.C:3
 combinedFit.C:4
 combinedFit.C:5
 combinedFit.C:6
 combinedFit.C:7
 combinedFit.C:8
 combinedFit.C:9
 combinedFit.C:10
 combinedFit.C:11
 combinedFit.C:12
 combinedFit.C:13
 combinedFit.C:14
 combinedFit.C:15
 combinedFit.C:16
 combinedFit.C:17
 combinedFit.C:18
 combinedFit.C:19
 combinedFit.C:20
 combinedFit.C:21
 combinedFit.C:22
 combinedFit.C:23
 combinedFit.C:24
 combinedFit.C:25
 combinedFit.C:26
 combinedFit.C:27
 combinedFit.C:28
 combinedFit.C:29
 combinedFit.C:30
 combinedFit.C:31
 combinedFit.C:32
 combinedFit.C:33
 combinedFit.C:34
 combinedFit.C:35
 combinedFit.C:36
 combinedFit.C:37
 combinedFit.C:38
 combinedFit.C:39
 combinedFit.C:40
 combinedFit.C:41
 combinedFit.C:42
 combinedFit.C:43
 combinedFit.C:44
 combinedFit.C:45
 combinedFit.C:46
 combinedFit.C:47
 combinedFit.C:48
 combinedFit.C:49
 combinedFit.C:50
 combinedFit.C:51
 combinedFit.C:52
 combinedFit.C:53
 combinedFit.C:54
 combinedFit.C:55
 combinedFit.C:56
 combinedFit.C:57
 combinedFit.C:58
 combinedFit.C:59
 combinedFit.C:60
 combinedFit.C:61
 combinedFit.C:62
 combinedFit.C:63
 combinedFit.C:64
 combinedFit.C:65
 combinedFit.C:66
 combinedFit.C:67
 combinedFit.C:68
 combinedFit.C:69
 combinedFit.C:70
 combinedFit.C:71
 combinedFit.C:72
 combinedFit.C:73
 combinedFit.C:74
 combinedFit.C:75
 combinedFit.C:76
 combinedFit.C:77
 combinedFit.C:78
 combinedFit.C:79
 combinedFit.C:80
 combinedFit.C:81
 combinedFit.C:82
 combinedFit.C:83
 combinedFit.C:84
 combinedFit.C:85
 combinedFit.C:86
 combinedFit.C:87
 combinedFit.C:88
 combinedFit.C:89
 combinedFit.C:90
 combinedFit.C:91
 combinedFit.C:92
 combinedFit.C:93
 combinedFit.C:94
 combinedFit.C:95
 combinedFit.C:96
 combinedFit.C:97
 combinedFit.C:98
 combinedFit.C:99
 combinedFit.C:100
 combinedFit.C:101
 combinedFit.C:102
 combinedFit.C:103
 combinedFit.C:104
 combinedFit.C:105
 combinedFit.C:106
 combinedFit.C:107
 combinedFit.C:108
 combinedFit.C:109
 combinedFit.C:110
 combinedFit.C:111
 combinedFit.C:112
 combinedFit.C:113
 combinedFit.C:114
 combinedFit.C:115
 combinedFit.C:116
 combinedFit.C:117
 combinedFit.C:118
 combinedFit.C:119
 combinedFit.C:120
 combinedFit.C:121
 combinedFit.C:122
 combinedFit.C:123
 combinedFit.C:124
 combinedFit.C:125
 combinedFit.C:126
 combinedFit.C:127
 combinedFit.C:128
 combinedFit.C:129
 combinedFit.C:130
 combinedFit.C:131
 combinedFit.C:132
 combinedFit.C:133
 combinedFit.C:134
 combinedFit.C:135
 combinedFit.C:136
 combinedFit.C:137
 combinedFit.C:138
 combinedFit.C:139
 combinedFit.C:140
 combinedFit.C:141
 combinedFit.C:142
 combinedFit.C:143
 combinedFit.C:144
 combinedFit.C:145
 combinedFit.C:146