ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
langaus.C
Go to the documentation of this file.
1 //-----------------------------------------------------------------------
2 //
3 // Convoluted Landau and Gaussian Fitting Function
4 // (using ROOT's Landau and Gauss functions)
5 //
6 // Based on a Fortran code by R.Fruehwirth (fruhwirth@hephy.oeaw.ac.at)
7 // Adapted for C++/ROOT by H.Pernegger (Heinz.Pernegger@cern.ch) and
8 // Markus Friedl (Markus.Friedl@cern.ch)
9 //
10 // to execute this example, do:
11 // root > .x langaus.C
12 // or
13 // root > .x langaus.C++
14 //
15 //-----------------------------------------------------------------------
16 
17 #include "TH1.h"
18 #include "TF1.h"
19 #include "TROOT.h"
20 #include "TStyle.h"
21 #include "TMath.h"
22 
24 
25  //Fit parameters:
26  //par[0]=Width (scale) parameter of Landau density
27  //par[1]=Most Probable (MP, location) parameter of Landau density
28  //par[2]=Total area (integral -inf to inf, normalization constant)
29  //par[3]=Width (sigma) of convoluted Gaussian function
30  //
31  //In the Landau distribution (represented by the CERNLIB approximation),
32  //the maximum is located at x=-0.22278298 with the location parameter=0.
33  //This shift is corrected within this function, so that the actual
34  //maximum is identical to the MP parameter.
35 
36  // Numeric constants
37  Double_t invsq2pi = 0.3989422804014; // (2 pi)^(-1/2)
38  Double_t mpshift = -0.22278298; // Landau maximum location
39 
40  // Control constants
41  Double_t np = 100.0; // number of convolution steps
42  Double_t sc = 5.0; // convolution extends to +-sc Gaussian sigmas
43 
44  // Variables
45  Double_t xx;
46  Double_t mpc;
47  Double_t fland;
48  Double_t sum = 0.0;
49  Double_t xlow,xupp;
50  Double_t step;
51  Double_t i;
52 
53 
54  // MP shift correction
55  mpc = par[1] - mpshift * par[0];
56 
57  // Range of convolution integral
58  xlow = x[0] - sc * par[3];
59  xupp = x[0] + sc * par[3];
60 
61  step = (xupp-xlow) / np;
62 
63  // Convolution integral of Landau and Gaussian by sum
64  for(i=1.0; i<=np/2; i++) {
65  xx = xlow + (i-.5) * step;
66  fland = TMath::Landau(xx,mpc,par[0]) / par[0];
67  sum += fland * TMath::Gaus(x[0],xx,par[3]);
68 
69  xx = xupp - (i-.5) * step;
70  fland = TMath::Landau(xx,mpc,par[0]) / par[0];
71  sum += fland * TMath::Gaus(x[0],xx,par[3]);
72  }
73 
74  return (par[2] * step * sum * invsq2pi / par[3]);
75 }
76 
77 
78 
79 TF1 *langaufit(TH1F *his, Double_t *fitrange, Double_t *startvalues, Double_t *parlimitslo, Double_t *parlimitshi, Double_t *fitparams, Double_t *fiterrors, Double_t *ChiSqr, Int_t *NDF)
80 {
81  // Once again, here are the Landau * Gaussian parameters:
82  // par[0]=Width (scale) parameter of Landau density
83  // par[1]=Most Probable (MP, location) parameter of Landau density
84  // par[2]=Total area (integral -inf to inf, normalization constant)
85  // par[3]=Width (sigma) of convoluted Gaussian function
86  //
87  // Variables for langaufit call:
88  // his histogram to fit
89  // fitrange[2] lo and hi boundaries of fit range
90  // startvalues[4] reasonable start values for the fit
91  // parlimitslo[4] lower parameter limits
92  // parlimitshi[4] upper parameter limits
93  // fitparams[4] returns the final fit parameters
94  // fiterrors[4] returns the final fit errors
95  // ChiSqr returns the chi square
96  // NDF returns ndf
97 
98  Int_t i;
99  Char_t FunName[100];
100 
101  sprintf(FunName,"Fitfcn_%s",his->GetName());
102 
103  TF1 *ffitold = (TF1*)gROOT->GetListOfFunctions()->FindObject(FunName);
104  if (ffitold) delete ffitold;
105 
106  TF1 *ffit = new TF1(FunName,langaufun,fitrange[0],fitrange[1],4);
107  ffit->SetParameters(startvalues);
108  ffit->SetParNames("Width","MP","Area","GSigma");
109 
110  for (i=0; i<4; i++) {
111  ffit->SetParLimits(i, parlimitslo[i], parlimitshi[i]);
112  }
113 
114  his->Fit(FunName,"RB0"); // fit within specified range, use ParLimits, do not plot
115 
116  ffit->GetParameters(fitparams); // obtain fit parameters
117  for (i=0; i<4; i++) {
118  fiterrors[i] = ffit->GetParError(i); // obtain fit parameter errors
119  }
120  ChiSqr[0] = ffit->GetChisquare(); // obtain chi^2
121  NDF[0] = ffit->GetNDF(); // obtain ndf
122 
123  return (ffit); // return fit function
124 
125 }
126 
127 
128 Int_t langaupro(Double_t *params, Double_t &maxx, Double_t &FWHM) {
129 
130  // Seaches for the location (x value) at the maximum of the
131  // Landau-Gaussian convolute and its full width at half-maximum.
132  //
133  // The search is probably not very efficient, but it's a first try.
134 
135  Double_t p,x,fy,fxr,fxl;
136  Double_t step;
137  Double_t l,lold;
138  Int_t i = 0;
139  Int_t MAXCALLS = 10000;
140 
141 
142  // Search for maximum
143 
144  p = params[1] - 0.1 * params[0];
145  step = 0.05 * params[0];
146  lold = -2.0;
147  l = -1.0;
148 
149 
150  while ( (l != lold) && (i < MAXCALLS) ) {
151  i++;
152 
153  lold = l;
154  x = p + step;
155  l = langaufun(&x,params);
156 
157  if (l < lold)
158  step = -step/10;
159 
160  p += step;
161  }
162 
163  if (i == MAXCALLS)
164  return (-1);
165 
166  maxx = x;
167 
168  fy = l/2;
169 
170 
171  // Search for right x location of fy
172 
173  p = maxx + params[0];
174  step = params[0];
175  lold = -2.0;
176  l = -1e300;
177  i = 0;
178 
179 
180  while ( (l != lold) && (i < MAXCALLS) ) {
181  i++;
182 
183  lold = l;
184  x = p + step;
185  l = TMath::Abs(langaufun(&x,params) - fy);
186 
187  if (l > lold)
188  step = -step/10;
189 
190  p += step;
191  }
192 
193  if (i == MAXCALLS)
194  return (-2);
195 
196  fxr = x;
197 
198 
199  // Search for left x location of fy
200 
201  p = maxx - 0.5 * params[0];
202  step = -params[0];
203  lold = -2.0;
204  l = -1e300;
205  i = 0;
206 
207  while ( (l != lold) && (i < MAXCALLS) ) {
208  i++;
209 
210  lold = l;
211  x = p + step;
212  l = TMath::Abs(langaufun(&x,params) - fy);
213 
214  if (l > lold)
215  step = -step/10;
216 
217  p += step;
218  }
219 
220  if (i == MAXCALLS)
221  return (-3);
222 
223 
224  fxl = x;
225 
226  FWHM = fxr - fxl;
227  return (0);
228 }
229 
230 void langaus() {
231  // Fill Histogram
232  Int_t data[100] = {0,0,0,0,0,0,2,6,11,18,18,55,90,141,255,323,454,563,681,
233  737,821,796,832,720,637,558,519,460,357,291,279,241,212,
234  153,164,139,106,95,91,76,80,80,59,58,51,30,49,23,35,28,23,
235  22,27,27,24,20,16,17,14,20,12,12,13,10,17,7,6,12,6,12,4,
236  9,9,10,3,4,5,2,4,1,5,5,1,7,1,6,3,3,3,4,5,4,4,2,2,7,2,4};
237  TH1F *hSNR = new TH1F("snr","Signal-to-noise",400,0,400);
238 
239  for (Int_t i=0; i<100; i++) hSNR->Fill(i,data[i]);
240 
241  // Fitting SNR histo
242  printf("Fitting...\n");
243 
244  // Setting fit range and start values
245  Double_t fr[2];
246  Double_t sv[4], pllo[4], plhi[4], fp[4], fpe[4];
247  fr[0]=0.3*hSNR->GetMean();
248  fr[1]=3.0*hSNR->GetMean();
249 
250  pllo[0]=0.5; pllo[1]=5.0; pllo[2]=1.0; pllo[3]=0.4;
251  plhi[0]=5.0; plhi[1]=50.0; plhi[2]=1000000.0; plhi[3]=5.0;
252  sv[0]=1.8; sv[1]=20.0; sv[2]=50000.0; sv[3]=3.0;
253 
254  Double_t chisqr;
255  Int_t ndf;
256  TF1 *fitsnr = langaufit(hSNR,fr,sv,pllo,plhi,fp,fpe,&chisqr,&ndf);
257 
258  Double_t SNRPeak, SNRFWHM;
259  langaupro(fp,SNRPeak,SNRFWHM);
260 
261  printf("Fitting done\nPlotting results...\n");
262 
263  // Global style settings
264  gStyle->SetOptStat(1111);
265  gStyle->SetOptFit(111);
266  gStyle->SetLabelSize(0.03,"x");
267  gStyle->SetLabelSize(0.03,"y");
268 
269  hSNR->GetXaxis()->SetRange(0,70);
270  hSNR->Draw();
271  fitsnr->Draw("lsame");
272 }
273 
double par[1]
Definition: unuranDistr.cxx:38
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3159
Double_t Landau(Double_t x, Double_t mpv=0, Double_t sigma=1, Bool_t norm=kFALSE)
The LANDAU function.
Definition: TMath.cxx:473
virtual void SetParameters(const Double_t *params)
Definition: TF1.h:432
virtual void SetParNames(const char *name0="p0", const char *name1="p1", const char *name2="p2", const char *name3="p3", const char *name4="p4", const char *name5="p5", const char *name6="p6", const char *name7="p7", const char *name8="p8", const char *name9="p9", const char *name10="p10")
Set up to 10 parameter names.
Definition: TF1.cxx:3147
R__EXTERN TStyle * gStyle
Definition: TStyle.h:423
TF1 * langaufit(TH1F *his, Double_t *fitrange, Double_t *startvalues, Double_t *parlimitslo, Double_t *parlimitshi, Double_t *fitparams, Double_t *fiterrors, Double_t *ChiSqr, Int_t *NDF)
Definition: langaus.C:79
#define gROOT
Definition: TROOT.h:344
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:570
void langaus()
Definition: langaus.C:230
int Int_t
Definition: RtypesCore.h:41
virtual Double_t GetParError(Int_t ipar) const
Return value of parameter number ipar.
Definition: TF1.cxx:1608
virtual void Draw(Option_t *option="")
Draw this function with its current attributes.
Definition: TF1.cxx:1052
Int_t langaupro(Double_t *params, Double_t &maxx, Double_t &FWHM)
Definition: langaus.C:128
virtual Int_t GetNDF() const
Return the number of degrees of freedom in the fit the fNDF parameter has been previously computed du...
Definition: TF1.cxx:1568
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
Double_t x[n]
Definition: legend1.C:17
Double_t GetChisquare() const
Definition: TF1.h:329
tuple np
Definition: multifit.py:30
virtual void SetRange(Int_t first=0, Int_t last=0)
Set the viewing range for the axis from bin first to last.
Definition: TAxis.cxx:831
virtual void SetParLimits(Int_t ipar, Double_t parmin, Double_t parmax)
Set limits for parameter ipar.
Definition: TF1.cxx:3183
virtual Double_t GetMean(Int_t axis=1) const
For axis = 1,2 or 3 returns the mean value of the histogram along X,Y or Z axis.
Definition: TH1.cxx:7014
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2878
void SetOptFit(Int_t fit=1)
The type of information about fit parameters printed in the histogram statistics box can be selected ...
Definition: TStyle.cxx:1204
TLine * l
Definition: textangle.C:4
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
Double_t Gaus(Double_t x, Double_t mean=0, Double_t sigma=1, Bool_t norm=kFALSE)
Calculate a gaussian function with mean and sigma.
Definition: TMath.cxx:453
double Double_t
Definition: RtypesCore.h:55
ClassImp(TMCParticle) void TMCParticle printf(": p=(%7.3f,%7.3f,%9.3f) ;", fPx, fPy, fPz)
virtual Double_t * GetParameters() const
Definition: TF1.h:358
void SetLabelSize(Float_t size=0.04, Option_t *axis="X")
Set size of axis labels.
Definition: TStyle.cxx:1063
char Char_t
Definition: RtypesCore.h:29
1-Dim function class
Definition: TF1.h:149
void SetOptStat(Int_t stat=1)
The type of information printed in the histogram statistics box can be selected via the parameter mod...
Definition: TStyle.cxx:1252
Double_t langaufun(Double_t *x, Double_t *par)
Definition: langaus.C:23
virtual TFitResultPtr Fit(const char *formula, Option_t *option="", Option_t *goption="", Double_t xmin=0, Double_t xmax=0)
Fit histogram with function fname.
Definition: TH1.cxx:3607
TAxis * GetXaxis()
Definition: TH1.h:319