Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
fitExclude.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_fit
3 /// \notebook -js
4 /// Illustrates how to fit excluding points in a given range.
5 ///
6 /// \macro_image
7 /// \macro_output
8 /// \macro_code
9 ///
10 /// \author Rene Brun
11 
12 #include "TH1.h"
13 #include "TF1.h"
14 #include "TList.h"
15 
16 Bool_t reject;
17 Double_t fline(Double_t *x, Double_t *par)
18 {
19  if (reject && x[0] > 2.5 && x[0] < 3.5) {
20  TF1::RejectPoint();
21  return 0;
22  }
23  return par[0] + par[1]*x[0];
24 }
25 
26 void fitExclude() {
27  //Create a source function
28  TF1 *f1 = new TF1("f1","[0] +[1]*x +gaus(2)",0,5);
29  f1->SetParameters(6,-1,5,3,0.2);
30  // create and fill histogram according to the source function
31  TH1F *h = new TH1F("h","background + signal",100,0,5);
32  h->FillRandom("f1",2000);
33  TF1 *fl = new TF1("fl",fline,0,5,2);
34  fl->SetParameters(2,-1);
35  //fit only the linear background excluding the signal area
36  reject = kTRUE;
37  h->Fit(fl,"0");
38  reject = kFALSE;
39  //store 2 separate functions for visualization
40  TF1 *fleft = new TF1("fleft",fline,0,2.5,2);
41  fleft->SetParameters(fl->GetParameters());
42  h->GetListOfFunctions()->Add(fleft);
43  gROOT->GetListOfFunctions()->Remove(fleft);
44  TF1 *fright = new TF1("fright",fline,3.5,5,2);
45  fright->SetParameters(fl->GetParameters());
46  h->GetListOfFunctions()->Add(fright);
47  gROOT->GetListOfFunctions()->Remove(fright);
48  h->Draw();
49 }
50