From $ROOTSYS/tutorials/fit/fit2a.C

#include "TF2.h"
#include "TH2.h"
#include "TCutG.h"
#include "TMath.h"
#include "TCanvas.h"
#include "TStyle.h"


//+ Fitting a 2-D histogram (a variant)
// This tutorial illustrates :
//  - how to create a 2-d function
//  - fill a 2-d histogram randomly from this function
//  - fit the histogram
//  - display the fitted function on top of the histogram (lego-plot)
//    using a surface plot in a sub-range of the histogram.
//
// This example can be executed via the interpreter or/and the compiler
//   root > .x fit2a.C
//   root > .x fit2a.C++
//Author: Rene Brun
         
Double_t g2(Double_t *x, Double_t *par) {
   Double_t r1 = Double_t((x[0]-par[1])/par[2]);
   Double_t r2 = Double_t((x[1]-par[3])/par[4]);
   return par[0]*TMath::Exp(-0.5*(r1*r1+r2*r2));
}   
Double_t fun2(Double_t *x, Double_t *par) {
   Double_t *p1 = &par[0];
   Double_t *p2 = &par[5];
   Double_t *p3 = &par[10];
   Double_t result = g2(x,p1) + g2(x,p2) + g2(x,p3);
   return result;
}

TCanvas *fit2a() {
   TCanvas *c = new TCanvas();
   gStyle->SetOptStat(kTRUE); 
   gStyle->SetPalette(1); 
   const Int_t npar = 15;
   Double_t f2params[npar] = {100,-3,3,-3,3,160,0,0.8,0,0.9,40,4,0.7,4,0.7};
   TF2 *f2 = new TF2("f2",fun2,-10,10,-10,10, npar);
   f2->SetParameters(f2params);

   //Create an histogram and fill it randomly with f2
   TH2F *h2 = new TH2F("h2","From f2",40,-10,10,40,-10,10);
   Int_t nentries = 100000;
   h2->FillRandom("f2",nentries);
   //Fit h2 with original function f2
   Float_t ratio = 4*nentries/100000;
   f2params[ 0] *= ratio;
   f2params[ 5] *= ratio;
   f2params[10] *= ratio;
   f2->SetParameters(f2params);
   h2->Fit("f2","N");
   TCutG *cutg = new TCutG("cutg",5);
   cutg->SetPoint(0,-7,-7);
   cutg->SetPoint(1, 2,-7);
   cutg->SetPoint(2, 2, 2);
   cutg->SetPoint(3,-7, 2);
   cutg->SetPoint(4,-7,-7);
   h2->Draw("lego2 0");
   h2->SetFillColor(38);
   f2->SetNpx(80);
   f2->SetNpy(80);
   f2->Draw("surf1 same bb [cutg]");
   return c;
}
 fit2a.C:1
 fit2a.C:2
 fit2a.C:3
 fit2a.C:4
 fit2a.C:5
 fit2a.C:6
 fit2a.C:7
 fit2a.C:8
 fit2a.C:9
 fit2a.C:10
 fit2a.C:11
 fit2a.C:12
 fit2a.C:13
 fit2a.C:14
 fit2a.C:15
 fit2a.C:16
 fit2a.C:17
 fit2a.C:18
 fit2a.C:19
 fit2a.C:20
 fit2a.C:21
 fit2a.C:22
 fit2a.C:23
 fit2a.C:24
 fit2a.C:25
 fit2a.C:26
 fit2a.C:27
 fit2a.C:28
 fit2a.C:29
 fit2a.C:30
 fit2a.C:31
 fit2a.C:32
 fit2a.C:33
 fit2a.C:34
 fit2a.C:35
 fit2a.C:36
 fit2a.C:37
 fit2a.C:38
 fit2a.C:39
 fit2a.C:40
 fit2a.C:41
 fit2a.C:42
 fit2a.C:43
 fit2a.C:44
 fit2a.C:45
 fit2a.C:46
 fit2a.C:47
 fit2a.C:48
 fit2a.C:49
 fit2a.C:50
 fit2a.C:51
 fit2a.C:52
 fit2a.C:53
 fit2a.C:54
 fit2a.C:55
 fit2a.C:56
 fit2a.C:57
 fit2a.C:58
 fit2a.C:59
 fit2a.C:60
 fit2a.C:61
 fit2a.C:62
 fit2a.C:63
 fit2a.C:64
 fit2a.C:65
 fit2a.C:66
 fit2a.C:67
 fit2a.C:68