{
"cells": [
{
"cell_type": "markdown",
"id": "a1b7403b",
"metadata": {},
"source": [
"# FittingDemo\n",
"Example for fitting signal/background.\n",
"This example can be executed with:\n",
"\n",
"```cpp\n",
"root > .x FittingDemo.C (using the cling interpreter)\n",
"root > .x FittingDemo.C+ (using the native complier via ACLIC)\n",
"```\n",
"\n",
"\n",
"\n",
"\n",
"**Author:** Rene Brun \n",
"This notebook tutorial was automatically generated with ROOTBOOK-izer from the macro found in the ROOT repository on Tuesday, May 19, 2026 at 08:24 PM."
]
},
{
"cell_type": "code",
"execution_count": 1,
"id": "06a1226b",
"metadata": {
"collapsed": false,
"execution": {
"iopub.execute_input": "2026-05-19T20:25:02.361444Z",
"iopub.status.busy": "2026-05-19T20:25:02.361333Z",
"iopub.status.idle": "2026-05-19T20:25:02.368335Z",
"shell.execute_reply": "2026-05-19T20:25:02.367869Z"
}
},
"outputs": [],
"source": [
"%%cpp -d\n",
"#include \"TH1.h\"\n",
"#include \"TMath.h\"\n",
"#include \"TF1.h\"\n",
"#include \"TLegend.h\"\n",
"#include \"TCanvas.h\""
]
},
{
"cell_type": "markdown",
"id": "4cdd52ec",
"metadata": {},
"source": [
" Quadratic background function\n",
" "
]
},
{
"cell_type": "code",
"execution_count": 2,
"id": "48ecf513",
"metadata": {
"collapsed": false,
"execution": {
"iopub.execute_input": "2026-05-19T20:25:02.370182Z",
"iopub.status.busy": "2026-05-19T20:25:02.370067Z",
"iopub.status.idle": "2026-05-19T20:25:02.372808Z",
"shell.execute_reply": "2026-05-19T20:25:02.372290Z"
}
},
"outputs": [],
"source": [
"%%cpp -d\n",
"double background(double *x, double *par) {\n",
" return par[0] + par[1]*x[0] + par[2]*x[0]*x[0];\n",
"}"
]
},
{
"cell_type": "markdown",
"id": "e59af450",
"metadata": {},
"source": [
" Lorenzian Peak function\n",
" "
]
},
{
"cell_type": "code",
"execution_count": 3,
"id": "910f7d1e",
"metadata": {
"collapsed": false,
"execution": {
"iopub.execute_input": "2026-05-19T20:25:02.374385Z",
"iopub.status.busy": "2026-05-19T20:25:02.374269Z",
"iopub.status.idle": "2026-05-19T20:25:02.376911Z",
"shell.execute_reply": "2026-05-19T20:25:02.376416Z"
}
},
"outputs": [],
"source": [
"%%cpp -d\n",
"double lorentzianPeak(double *x, double *par) {\n",
" return (0.5*par[0]*par[1]/TMath::Pi()) /\n",
" TMath::Max( 1.e-10,(x[0]-par[2])*(x[0]-par[2])\n",
" + .25*par[1]*par[1]);\n",
"}"
]
},
{
"cell_type": "markdown",
"id": "e513dcca",
"metadata": {},
"source": [
" Sum of background and peak function\n",
" "
]
},
{
"cell_type": "code",
"execution_count": 4,
"id": "b05ae6df",
"metadata": {
"collapsed": false,
"execution": {
"iopub.execute_input": "2026-05-19T20:25:02.378350Z",
"iopub.status.busy": "2026-05-19T20:25:02.378239Z",
"iopub.status.idle": "2026-05-19T20:25:02.380427Z",
"shell.execute_reply": "2026-05-19T20:25:02.379972Z"
}
},
"outputs": [],
"source": [
"%%cpp -d\n",
"double fitFunction(double *x, double *par) {\n",
" return background(x,par) + lorentzianPeak(x,&par[3]);\n",
"}"
]
},
{
"cell_type": "markdown",
"id": "0448db75",
"metadata": {},
"source": [
"Bevington Exercise by Peter Malzacher, modified by Rene Brun"
]
},
{
"cell_type": "code",
"execution_count": 5,
"id": "4ebc3dff",
"metadata": {
"collapsed": false,
"execution": {
"iopub.execute_input": "2026-05-19T20:25:02.382000Z",
"iopub.status.busy": "2026-05-19T20:25:02.381889Z",
"iopub.status.idle": "2026-05-19T20:25:02.730322Z",
"shell.execute_reply": "2026-05-19T20:25:02.729628Z"
}
},
"outputs": [
{
"name": "stderr",
"output_type": "stream",
"text": [
"input_line_48:4:3: warning: 'data' shadows a declaration with the same name in the 'std' namespace; use '::data' to reference this declaration\n",
" double data[nBins] = { 6, 1,10,12, 6,13,23,22,15,21,\n",
" ^\n"
]
}
],
"source": [
" const int nBins = 60;\n",
"\n",
" double data[nBins] = { 6, 1,10,12, 6,13,23,22,15,21,\n",
" 23,26,36,25,27,35,40,44,66,81,\n",
" 75,57,48,45,46,41,35,36,53,32,\n",
" 40,37,38,31,36,44,42,37,32,32,\n",
" 43,44,35,33,33,39,29,41,32,44,\n",
" 26,39,29,35,32,21,21,15,25,15};\n",
" TCanvas *c1 = new TCanvas(\"c1\",\"Fitting Demo\",10,10,700,500);\n",
" c1->SetFillColor(33);\n",
" c1->SetFrameFillColor(41);\n",
" c1->SetGrid();\n",
"\n",
" TH1F *histo = new TH1F(\"histo\",\n",
" \"Lorentzian Peak on Quadratic Background\",60,0,3);\n",
" histo->SetMarkerStyle(21);\n",
" histo->SetMarkerSize(0.8);\n",
" histo->SetStats(false);\n",
"\n",
" for(int i=0; i < nBins; i++) histo->SetBinContent(i+1,data[i]);"
]
},
{
"cell_type": "markdown",
"id": "3c6a64fb",
"metadata": {},
"source": [
"create a TF1 with the range from 0 to 3 and 6 parameters"
]
},
{
"cell_type": "code",
"execution_count": 6,
"id": "1bc2b90a",
"metadata": {
"collapsed": false,
"execution": {
"iopub.execute_input": "2026-05-19T20:25:02.732174Z",
"iopub.status.busy": "2026-05-19T20:25:02.732051Z",
"iopub.status.idle": "2026-05-19T20:25:02.938697Z",
"shell.execute_reply": "2026-05-19T20:25:02.938015Z"
}
},
"outputs": [],
"source": [
" TF1 *fitFcn = new TF1(\"fitFcn\",fitFunction,0,3,6);\n",
" fitFcn->SetNpx(500);\n",
" fitFcn->SetLineWidth(4);\n",
" fitFcn->SetLineColor(kMagenta);"
]
},
{
"cell_type": "markdown",
"id": "69fb9f99",
"metadata": {},
"source": [
"first try without starting values for the parameters\n",
"This defaults to 1 for each param.\n",
"this results in an ok fit for the polynomial function\n",
"however the non-linear part (lorenzian) does not\n",
"respond well."
]
},
{
"cell_type": "code",
"execution_count": 7,
"id": "1be41deb",
"metadata": {
"collapsed": false,
"execution": {
"iopub.execute_input": "2026-05-19T20:25:02.941129Z",
"iopub.status.busy": "2026-05-19T20:25:02.941009Z",
"iopub.status.idle": "2026-05-19T20:25:03.151416Z",
"shell.execute_reply": "2026-05-19T20:25:03.150700Z"
}
},
"outputs": [
{
"name": "stdout",
"output_type": "stream",
"text": [
"****************************************\n",
"Minimizer is Minuit2 / Migrad\n",
"Chi2 = 58.9284\n",
"NDf = 54\n",
"Edm = 9.25608e-07\n",
"NCalls = 605\n",
"p0 = -0.864558 +/- 0.891619 \n",
"p1 = 45.8427 +/- 2.65587 \n",
"p2 = -13.3213 +/- 0.980096 \n",
"p3 = 13.8077 +/- 2.23262 \n",
"p4 = 0.172312 +/- 0.0364912 \n",
"p5 = 0.987278 +/- 0.0112098 \n"
]
}
],
"source": [
" fitFcn->SetParameters(1,1,1,1,1,1);\n",
" histo->Fit(\"fitFcn\",\"0\");"
]
},
{
"cell_type": "markdown",
"id": "d52060c0",
"metadata": {},
"source": [
"second try: set start values for some parameters"
]
},
{
"cell_type": "code",
"execution_count": 8,
"id": "fa6d5f28",
"metadata": {
"collapsed": false,
"execution": {
"iopub.execute_input": "2026-05-19T20:25:03.153240Z",
"iopub.status.busy": "2026-05-19T20:25:03.153121Z",
"iopub.status.idle": "2026-05-19T20:25:03.362865Z",
"shell.execute_reply": "2026-05-19T20:25:03.361765Z"
}
},
"outputs": [
{
"name": "stdout",
"output_type": "stream",
"text": [
"Minuit2Minimizer: Minimize with max-calls 1780 convergence for edm < 0.01 strategy 1\n",
"Minuit2Minimizer : Valid minimum - status = 0\n",
"FVAL = 58.9283860619491549\n",
"Edm = 4.48185715140414416e-13\n",
"Nfcn = 163\n",
"p0\t = -0.864713\t +/- 0.891795\n",
"p1\t = 45.8434\t +/- 2.64223\n",
"p2\t = -13.3214\t +/- 0.976969\n",
"p3\t = 13.8074\t +/- 2.1776\n",
"p4\t = 0.172309\t +/- 0.0358281\n",
"p5\t = 0.987281\t +/- 0.0112683\n",
"\n",
"Covariance Matrix:\n",
"\n",
" \t p0 p1 p2 p3 p4 p5\n",
"p0 \t 0.7953 -1.2054 0.34842 -0.15945 -0.0037284 0.00042265\n",
"p1 \t -1.2054 6.9814 -2.5255 -3.0272 -0.037043 -0.0018129\n",
"p2 \t 0.34842 -2.5255 0.95447 1.1587 0.014336 0.00058758\n",
"p3 \t -0.15945 -3.0272 1.1587 4.7419 0.05537 0.0017371\n",
"p4 \t -0.0037284 -0.037043 0.014336 0.05537 0.0012837 2.8395e-05\n",
"p5 \t 0.00042265 -0.0018129 0.00058758 0.0017371 2.8395e-05 0.00012697\n",
"\n",
"Correlation Matrix:\n",
"\n",
" \t p0 p1 p2 p3 p4 p5\n",
"p0 \t 1 -0.51155 0.39991 -0.08211 -0.11669 0.042058\n",
"p1 \t -0.51155 1 -0.97834 -0.52614 -0.3913 -0.06089\n",
"p2 \t 0.39991 -0.97834 1 0.54464 0.40958 0.053373\n",
"p3 \t -0.08211 -0.52614 0.54464 1 0.7097 0.070793\n",
"p4 \t -0.11669 -0.3913 0.40958 0.7097 1 0.070334\n",
"p5 \t 0.042058 -0.06089 0.053373 0.070793 0.070334 1\n",
"****************************************\n",
"Minimizer is Minuit2 / Migrad\n",
"Chi2 = 58.9284\n",
"NDf = 54\n",
"Edm = 4.48186e-13\n",
"NCalls = 163\n",
"p0 = -0.864713 +/- 0.891795 \n",
"p1 = 45.8434 +/- 2.64223 \n",
"p2 = -13.3214 +/- 0.976969 \n",
"p3 = 13.8074 +/- 2.1776 \n",
"p4 = 0.172309 +/- 0.0358281 \n",
"p5 = 0.987281 +/- 0.0112683 \n"
]
},
{
"name": "stderr",
"output_type": "stream",
"text": [
"Info in : MnSeedGenerator Computing seed using NumericalGradient calculator\n",
"Info in : MnSeedGenerator Evaluated function and gradient in 34.499 μs\n",
"Info in : MnSeedGenerator Initial state: FCN = 60.85761877 Edm = 2.042788608 NCalls = 23\n",
"Info in : MnSeedGenerator Initial state \n",
" Minimum value : 60.85761877\n",
" Edm : 2.042788608\n",
" Internal parameters:\t[ -0.8645583494 45.84274925 -13.32134505 13.80774266 0.2 1]\t\n",
" Internal gradient :\t[ 0.4275838297 0.4487206938 0.8008642186 -0.6674368476 53.51638699 153.2778371]\t\n",
" Internal covariance matrix:\n",
"[[ 0.30404743 0 0 0 0 0]\n",
" [ 0 0.16270541 0 0 0 0]\n",
" [ 0 0 0.027215553 0 0 0]\n",
" [ 0 0 0 1.9959967 0 0]\n",
" [ 0 0 0 0 0.0011646905 0]\n",
" [ 0 0 0 0 0 0.00016346709]]]\n",
"Info in : VariableMetricBuilder Start iterating until Edm is < 2e-05 with call limit = 1780\n",
"Info in : VariableMetricBuilder 0 - FCN = 60.85761877 Edm = 2.042788608 NCalls = 23\n",
"Info in : VariableMetricBuilder 1 - FCN = 59.09847664 Edm = 0.1924298892 NCalls = 40\n",
"Info in : VariableMetricBuilder 2 - FCN = 58.95990622 Edm = 0.01887091831 NCalls = 54\n",
"Info in : VariableMetricBuilder 3 - FCN = 58.9375889 Edm = 0.009530429569 NCalls = 68\n",
"Info in : VariableMetricBuilder 4 - FCN = 58.93183994 Edm = 0.001565170778 NCalls = 82\n",
"Info in : VariableMetricBuilder 5 - FCN = 58.92877036 Edm = 0.0004383416974 NCalls = 96\n",
"Info in : VariableMetricBuilder 6 - FCN = 58.92841786 Edm = 9.082253755e-07 NCalls = 110\n",
"Info in : MnHesse Done after 86.301 μs\n",
"Info in : VariableMetricBuilder After Hessian\n",
"Info in : VariableMetricBuilder 7 - FCN = 58.92841786 Edm = 3.179425109e-05 NCalls = 150\n",
"Info in : VariableMetricBuilder Tolerance not sufficient, continue minimization; Edm 3.17943e-05 Required 2e-05\n",
"Info in : VariableMetricBuilder 8 - FCN = 58.92838606 Edm = 4.481857151e-13 NCalls = 163\n",
"Info in : VariableMetricBuilder Stop iterating after 288.121 μs\n"
]
}
],
"source": [
" fitFcn->SetParameter(4,0.2); // width\n",
" fitFcn->SetParameter(5,1); // peak\n",
"\n",
" histo->Fit(\"fitFcn\",\"V+\",\"ep\");"
]
},
{
"cell_type": "markdown",
"id": "f30baf4e",
"metadata": {},
"source": [
"improve the picture:"
]
},
{
"cell_type": "code",
"execution_count": 9,
"id": "13cb5ac9",
"metadata": {
"collapsed": false,
"execution": {
"iopub.execute_input": "2026-05-19T20:25:03.364493Z",
"iopub.status.busy": "2026-05-19T20:25:03.364373Z",
"iopub.status.idle": "2026-05-19T20:25:03.570961Z",
"shell.execute_reply": "2026-05-19T20:25:03.570273Z"
}
},
"outputs": [],
"source": [
" TF1 *backFcn = new TF1(\"backFcn\",background,0,3,3);\n",
" backFcn->SetLineColor(kRed);\n",
" TF1 *signalFcn = new TF1(\"signalFcn\",lorentzianPeak,0,3,3);\n",
" signalFcn->SetLineColor(kBlue);\n",
" signalFcn->SetNpx(500);\n",
" double par[6];"
]
},
{
"cell_type": "markdown",
"id": "bb981b30",
"metadata": {},
"source": [
"writes the fit results into the par array"
]
},
{
"cell_type": "code",
"execution_count": 10,
"id": "21b928f5",
"metadata": {
"collapsed": false,
"execution": {
"iopub.execute_input": "2026-05-19T20:25:03.573242Z",
"iopub.status.busy": "2026-05-19T20:25:03.573119Z",
"iopub.status.idle": "2026-05-19T20:25:03.791824Z",
"shell.execute_reply": "2026-05-19T20:25:03.791281Z"
}
},
"outputs": [],
"source": [
" fitFcn->GetParameters(par);\n",
"\n",
" backFcn->SetParameters(par);\n",
" backFcn->Draw(\"same\");\n",
"\n",
" signalFcn->SetParameters(&par[3]);\n",
" signalFcn->Draw(\"same\");"
]
},
{
"cell_type": "markdown",
"id": "70beece6",
"metadata": {},
"source": [
"draw the legend"
]
},
{
"cell_type": "code",
"execution_count": 11,
"id": "44a001f8",
"metadata": {
"collapsed": false,
"execution": {
"iopub.execute_input": "2026-05-19T20:25:03.800348Z",
"iopub.status.busy": "2026-05-19T20:25:03.800213Z",
"iopub.status.idle": "2026-05-19T20:25:04.015818Z",
"shell.execute_reply": "2026-05-19T20:25:04.015416Z"
}
},
"outputs": [],
"source": [
" TLegend *legend=new TLegend(0.6,0.65,0.88,0.85);\n",
" legend->SetTextFont(72);\n",
" legend->SetTextSize(0.04);\n",
" legend->AddEntry(histo,\"Data\",\"lpe\");\n",
" legend->AddEntry(backFcn,\"Background fit\",\"l\");\n",
" legend->AddEntry(signalFcn,\"Signal fit\",\"l\");\n",
" legend->AddEntry(fitFcn,\"Global Fit\",\"l\");\n",
" legend->Draw();"
]
},
{
"cell_type": "markdown",
"id": "586cad7e",
"metadata": {},
"source": [
"Draw all canvases "
]
},
{
"cell_type": "code",
"execution_count": 12,
"id": "5ea78b19",
"metadata": {
"collapsed": false,
"execution": {
"iopub.execute_input": "2026-05-19T20:25:04.031231Z",
"iopub.status.busy": "2026-05-19T20:25:04.031089Z",
"iopub.status.idle": "2026-05-19T20:25:04.268032Z",
"shell.execute_reply": "2026-05-19T20:25:04.261326Z"
}
},
"outputs": [
{
"data": {
"text/html": [
"\n",
"\n",
"\n",
"
\n",
"\n",
"\n",
"\n"
],
"text/plain": [
""
]
},
"metadata": {},
"output_type": "display_data"
}
],
"source": [
"%jsroot on\n",
"gROOT->GetListOfCanvases()->Draw()"
]
}
],
"metadata": {
"kernelspec": {
"display_name": "ROOT C++",
"language": "c++",
"name": "root"
},
"language_info": {
"codemirror_mode": "text/x-c++src",
"file_extension": ".C",
"mimetype": " text/x-c++src",
"name": "c++"
}
},
"nbformat": 4,
"nbformat_minor": 5
}