Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
mathBeta.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_math
3 /// \notebook
4 /// Test the TMath::BetaDist and TMath::BetaDistI functions
5 ///
6 /// \macro_image
7 /// \macro_code
8 ///
9 /// \author Anna Kreshuk
10 
11 #include "TMath.h"
12 #include "TCanvas.h"
13 #include "TF1.h"
14 #include "TLegend.h"
15 
16 void mathBeta() {
17  TCanvas *c1=new TCanvas("c1", "TMath::BetaDist",600,800);
18  c1->Divide(1, 2);
19  TVirtualPad *pad1 = c1->cd(1);
20  pad1->SetGrid();
21  TF1 *fbeta = new TF1("fbeta", "TMath::BetaDist(x, [0], [1])", 0, 1);
22  fbeta->SetParameters(0.5, 0.5);
23  TF1 *f1 = fbeta->DrawCopy();
24  f1->SetLineColor(kRed);
25  f1->SetLineWidth(1);
26  fbeta->SetParameters(0.5, 2);
27  TF1 *f2 = fbeta->DrawCopy("same");
28  f2->SetLineColor(kGreen);
29  f2->SetLineWidth(1);
30  fbeta->SetParameters(2, 0.5);
31  TF1 *f3 = fbeta->DrawCopy("same");
32  f3->SetLineColor(kBlue);
33  f3->SetLineWidth(1);
34  fbeta->SetParameters(2, 2);
35  TF1 *f4 = fbeta->DrawCopy("same");
36  f4->SetLineColor(kMagenta);
37  f4->SetLineWidth(1);
38  TLegend *legend1 = new TLegend(.5,.7,.8,.9);
39  legend1->AddEntry(f1,"p=0.5 q=0.5","l");
40  legend1->AddEntry(f2,"p=0.5 q=2","l");
41  legend1->AddEntry(f3,"p=2 q=0.5","l");
42  legend1->AddEntry(f4,"p=2 q=2","l");
43  legend1->Draw();
44 
45  TVirtualPad *pad2 = c1->cd(2);
46  pad2->SetGrid();
47  TF1 *fbetai=new TF1("fbetai", "TMath::BetaDistI(x, [0], [1])", 0, 1);
48  fbetai->SetParameters(0.5, 0.5);
49  TF1 *g1=fbetai->DrawCopy();
50  g1->SetLineColor(kRed);
51  g1->SetLineWidth(1);
52  fbetai->SetParameters(0.5, 2);
53  TF1 *g2=fbetai->DrawCopy("same");
54  g2->SetLineColor(kGreen);
55  g2->SetLineWidth(1);
56  fbetai->SetParameters(2, 0.5);
57  TF1 *g3=fbetai->DrawCopy("same");
58  g3->SetLineColor(kBlue);
59  g3->SetLineWidth(1);
60  fbetai->SetParameters(2, 2);
61  TF1 *g4=fbetai->DrawCopy("same");
62  g4->SetLineColor(kMagenta);
63  g4->SetLineWidth(1);
64 
65  TLegend *legend2 = new TLegend(.7,.15,0.9,.35);
66  legend2->AddEntry(f1,"p=0.5 q=0.5","l");
67  legend2->AddEntry(f2,"p=0.5 q=2","l");
68  legend2->AddEntry(f3,"p=2 q=0.5","l");
69  legend2->AddEntry(f4,"p=2 q=2","l");
70  legend2->Draw();
71  c1->cd();
72 }