tStudent.C

Go to the documentation of this file.
00001 // Example macro describing the student t distribution
00002 //
00003 // root[0]: .x tStudent.C 
00004 //
00005 // It draws the pdf, the cdf and then 10 quantiles of the t Student distribution
00006 //
00007 // Author: Magdalena Slawinska
00008 
00009 #include "TH1.h"
00010 #include "TF1.h"
00011 #include "TCanvas.h"
00012 #include "TSystem.h"
00013 #include "TLegend.h"
00014 #include "TLegendEntry.h"
00015 #ifndef __CINT__
00016 #include "Math/DistFunc.h"
00017 #endif
00018 
00019 
00020 void tStudent()
00021 {
00022    gSystem->Load("libMathCore");
00023    gSystem->Load("libMathMore");
00024 
00025   int n=100;
00026   double a=-5.;
00027   double b=5.;
00028   //double r  = 3; 
00029   TF1* pdf = new TF1("pdf", "ROOT::Math::tdistribution_pdf(x,3.0)", a,b);
00030   TF1* cum = new TF1("cum", "ROOT::Math::tdistribution_cdf(x,3.0)", a,b);
00031 
00032 
00033   TH1D* quant = new TH1D("quant", "", 9, 0, 0.9);
00034   
00035 
00036   for(int i=1; i < 10; i++)
00037      quant->Fill((i-0.5)/10.0, ROOT::Math::tdistribution_quantile((1.0*i)/10, 3.0 ) );
00038 
00039   double xx[10];
00040   xx[0] = -1.5;
00041   for(int i=1; i<9; i++)
00042     xx[i]= quant->GetBinContent(i);
00043   xx[9] = 1.5;
00044   TH1D* pdfq[10];
00045   //int nbin = n/10.0;
00046   for(int i=0; i < 10; i++)
00047   {
00048     int nbin = n * (xx[i+1]-xx[i])/3.0+1.0; 
00049     TString name = "pdf" + i; 
00050     pdfq[i]= new TH1D(name, "", nbin,xx[i],xx[i+1] );
00051     for(int j=1; j<nbin; j++)
00052     {
00053        double x= j*(xx[i+1]-xx[i])/nbin + xx[i];
00054        pdfq[i]->SetBinContent(j, ROOT::Math::tdistribution_pdf(x,3));
00055 
00056     }
00057  
00058   }
00059 
00060  TCanvas *Canvas = new TCanvas("DistCanvas", "Student Distribution graphs", 10, 10, 1000, 800); 
00061  pdf->SetTitle("Student t distribution function");
00062  cum->SetTitle("Cumulative for Student t");
00063  quant->SetTitle("10-quantiles  for Student t");
00064  Canvas->Divide(2, 2);
00065  Canvas->cd(1);
00066  pdf->SetLineWidth(2);
00067  pdf->DrawCopy();
00068  Canvas->cd(2);
00069  cum->SetLineWidth(2);
00070  cum->SetLineColor(kRed);
00071  cum->Draw();
00072  Canvas->cd(3);
00073  quant->Draw();
00074  quant->SetLineWidth(2);
00075  quant->SetLineColor(kBlue);
00076  quant->SetStats(0);
00077  Canvas->cd(4);
00078  pdfq[0]->SetTitle("Student t & its quantiles");
00079  pdf->SetTitle("");
00080  pdf->Draw();
00081  //pdfq[0]->SetAxisRange(-1.5, 0, 1.5,1.0);
00082  pdfq[0]->SetTitle("Student t & its quantiles");
00083  for(int i=0; i < 10; i++)
00084  {
00085    pdfq[i]->SetStats(0);
00086    pdfq[i]->SetFillColor(i+1);
00087    pdfq[i]->Draw("same");
00088  }
00089  Canvas->Modified();
00090  Canvas->cd();
00091 }

Generated on Tue Jul 5 15:44:51 2011 for ROOT_528-00b_version by  doxygen 1.5.1