Logo ROOT  
Reference Guide
normalDist.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_math
3/// \notebook
4/// Tutorial illustrating the new statistical distributions functions (pdf, cdf and quantile)
5///
6/// \macro_image
7/// \macro_code
8///
9/// \author Anna Kreshuk
10
11#include "Math/DistFunc.h"
12#include "TF1.h"
13#include "TCanvas.h"
14#include "TSystem.h"
15#include "TLegend.h"
16#include "TAxis.h"
17
18void normalDist() {
19
20 TF1 *pdfunc = new TF1("pdf","ROOT::Math::normal_pdf(x, [0],[1])",-5,5);
21 TF1 *cdfunc = new TF1("cdf","ROOT::Math::normal_cdf(x, [0],[1])",-5,5);
22 TF1 *ccdfunc = new TF1("cdf_c","ROOT::Math::normal_cdf_c(x, [0])",-5,5);
23 TF1 *qfunc = new TF1("quantile","ROOT::Math::normal_quantile(x, [0])",0,1);
24 TF1 *cqfunc = new TF1("quantile_c","ROOT::Math::normal_quantile_c(x, [0])",0,1);
25
26 pdfunc->SetParameters(1.0,0.0); // set sigma to 1 and mean to zero
27 pdfunc->SetTitle("");
28 pdfunc->SetLineColor(kBlue);
29
30 pdfunc->GetXaxis()->SetLabelSize(0.06);
31 pdfunc->GetXaxis()->SetTitle("x");
32 pdfunc->GetXaxis()->SetTitleSize(0.07);
33 pdfunc->GetXaxis()->SetTitleOffset(0.55);
34 pdfunc->GetYaxis()->SetLabelSize(0.06);
35
36 cdfunc->SetParameters(1.0,0.0); // set sigma to 1 and mean to zero
37 cdfunc->SetTitle("");
38 cdfunc->SetLineColor(kRed);
39
40 cdfunc->GetXaxis()->SetLabelSize(0.06);
41 cdfunc->GetXaxis()->SetTitle("x");
42 cdfunc->GetXaxis()->SetTitleSize(0.07);
43 cdfunc->GetXaxis()->SetTitleOffset(0.55);
44
45 cdfunc->GetYaxis()->SetLabelSize(0.06);
46 cdfunc->GetYaxis()->SetTitle("p");
47 cdfunc->GetYaxis()->SetTitleSize(0.07);
48 cdfunc->GetYaxis()->SetTitleOffset(0.55);
49
50 ccdfunc->SetParameters(1.0,0.0); // set sigma to 1 and mean to zero
51 ccdfunc->SetTitle("");
52 ccdfunc->SetLineColor(kGreen);
53
54 qfunc->SetParameter(0, 1.0); // set sigma to 1
55 qfunc->SetTitle("");
56 qfunc->SetLineColor(kRed);
57 qfunc->SetNpx(1000); // to get more precision for p close to 0 or 1
58
59 qfunc->GetXaxis()->SetLabelSize(0.06);
60 qfunc->GetXaxis()->SetTitle("p");
61 qfunc->GetYaxis()->SetLabelSize(0.06);
62 qfunc->GetXaxis()->SetTitleSize(0.07);
63 qfunc->GetXaxis()->SetTitleOffset(0.55);
64 qfunc->GetYaxis()->SetTitle("x");
65 qfunc->GetYaxis()->SetTitleSize(0.07);
66 qfunc->GetYaxis()->SetTitleOffset(0.55);
67
68 cqfunc->SetParameter(0, 1.0); // set sigma to 1
69 cqfunc->SetTitle("");
70 cqfunc->SetLineColor(kGreen);
71 cqfunc->SetNpx(1000);
72
73 TCanvas * c1 = new TCanvas("c1","Normal Distributions",100,10,600,800);
74
75 c1->Divide(1,3);
76 c1->cd(1);
77
78 pdfunc->Draw();
79 TLegend *legend1 = new TLegend(0.583893,0.601973,0.885221,0.854151);
80 legend1->AddEntry(pdfunc,"normal_pdf","l");
81 legend1->Draw();
82
83 c1->cd(2);
84 cdfunc->Draw();
85 ccdfunc->Draw("same");
86 TLegend *legend2 = new TLegend(0.585605,0.462794,0.886933,0.710837);
87 legend2->AddEntry(cdfunc,"normal_cdf","l");
88 legend2->AddEntry(ccdfunc,"normal_cdf_c","l");
89 legend2->Draw();
90
91 c1->cd(3);
92 qfunc->Draw();
93 cqfunc->Draw("same");
94 TLegend *legend3 = new TLegend(0.315094,0.633668,0.695179,0.881711);
95 legend3->AddEntry(qfunc,"normal_quantile","l");
96 legend3->AddEntry(cqfunc,"normal_quantile_c","l");
97 legend3->Draw();
98}
@ kRed
Definition: Rtypes.h:64
@ kGreen
Definition: Rtypes.h:64
@ kBlue
Definition: Rtypes.h:64
virtual void SetTitleOffset(Float_t offset=1)
Set distance between the axis and the axis title.
Definition: TAttAxis.cxx:294
virtual void SetLabelSize(Float_t size=0.04)
Set size of axis labels.
Definition: TAttAxis.cxx:204
virtual void SetTitleSize(Float_t size=0.04)
Set size of axis title.
Definition: TAttAxis.cxx:304
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:40
The Canvas class.
Definition: TCanvas.h:27
1-Dim function class
Definition: TF1.h:210
TAxis * GetYaxis() const
Get y axis of the function.
Definition: TF1.cxx:2396
virtual void SetTitle(const char *title="")
Set function title if title has the form "fffffff;xxxx;yyyy", it is assumed that the function title i...
Definition: TF1.cxx:3551
virtual void SetNpx(Int_t npx=100)
Set the number of points used to draw the function.
Definition: TF1.cxx:3435
virtual void Draw(Option_t *option="")
Draw this function with its current attributes.
Definition: TF1.cxx:1320
virtual void SetParameters(const Double_t *params)
Definition: TF1.h:638
virtual void SetParameter(Int_t param, Double_t value)
Definition: TF1.h:628
TAxis * GetXaxis() const
Get x axis of the function.
Definition: TF1.cxx:2385
This class displays a legend box (TPaveText) containing several legend entries.
Definition: TLegend.h:23
TLegendEntry * AddEntry(const TObject *obj, const char *label="", Option_t *option="lpf")
Add a new entry to this legend.
Definition: TLegend.cxx:330
virtual void Draw(Option_t *option="")
Draw this legend with its current attributes.
Definition: TLegend.cxx:423
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition: TNamed.cxx:164
return c1
Definition: legend1.C:41