TwoHistoFit2D.C: Example to fit two histograms at the same time | Fitting tutorials | fit2.C: Fitting a 2-D histogram |
void fit1() { //Simple fitting example (1-d histogram with an interpreted function) //To see the output of this macro, click begin_html <a href="gif/fit1.gif">here</a>. end_html //Author: Rene Brun TCanvas *c1 = new TCanvas("c1_fit1","The Fit Canvas",200,10,700,500); c1->SetGridx(); c1->SetGridy(); c1->GetFrame()->SetFillColor(21); c1->GetFrame()->SetBorderMode(-1); c1->GetFrame()->SetBorderSize(5); gBenchmark->Start("fit1"); // // We connect the ROOT file generated in a previous tutorial // (see begin_html <a href="fillrandom.C.html">Filling histograms with random numbers from a function</a>) end_html // TString dir = gSystem->UnixPathName(gInterpreter->GetCurrentMacroName()); dir.ReplaceAll("fit1.C",""); dir.ReplaceAll("/./","/"); TFile *fill = TFile::Open("fillrandom.root"); if (!fill) { gROOT->ProcessLine(Form(".x %s../hist/fillrandom.C",dir.Data())); fill = TFile::Open("fillrandom.root"); if (!fill) return; } // // The function "ls()" lists the directory contents of this file // fill->ls(); // // Get object "sqroot" from the file. Undefined objects are searched // for using gROOT->FindObject("xxx"), e.g.: // TF1 *sqroot = (TF1*) gROOT.FindObject("sqroot") // sqroot->Print(); // // Now fit histogram h1f with the function sqroot // h1f->SetFillColor(45); h1f->Fit("sqroot"); // We now annotate the picture by creating a PaveText object // and displaying the list of commands in this macro // fitlabel = new TPaveText(0.6,0.3,0.9,0.80,"NDC"); fitlabel->SetTextAlign(12); fitlabel->SetFillColor(42); fitlabel->ReadFile(Form("%sfit1_C.C",dir.Data())); fitlabel->Draw(); c1->Update(); gBenchmark->Show("fit1"); } fit1.C:1 fit1.C:2 fit1.C:3 fit1.C:4 fit1.C:5 fit1.C:6 fit1.C:7 fit1.C:8 fit1.C:9 fit1.C:10 fit1.C:11 fit1.C:12 fit1.C:13 fit1.C:14 fit1.C:15 fit1.C:16 fit1.C:17 fit1.C:18 fit1.C:19 fit1.C:20 fit1.C:21 fit1.C:22 fit1.C:23 fit1.C:24 fit1.C:25 fit1.C:26 fit1.C:27 fit1.C:28 fit1.C:29 fit1.C:30 fit1.C:31 fit1.C:32 fit1.C:33 fit1.C:34 fit1.C:35 fit1.C:36 fit1.C:37 fit1.C:38 fit1.C:39 fit1.C:40 fit1.C:41 fit1.C:42 fit1.C:43 fit1.C:44 fit1.C:45 fit1.C:46 fit1.C:47 fit1.C:48 fit1.C:49 fit1.C:50 fit1.C:51 fit1.C:52 fit1.C:53 fit1.C:54 fit1.C:55 fit1.C:56 fit1.C:57 |
|