ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
LegendreAssoc.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_math
3 /// Example describing the usage of different kinds of Associate Legendre Polynomials
4 /// To execute the macro type in:
5 ///
6 /// ~~~ {.cpp}
7 /// root[0] .x LegendreAssoc.C
8 /// ~~~
9 ///
10 /// It draws common graphs for first 5
11 /// Associate Legendre Polynomials
12 /// and Spherical Associate Legendre Polynomials
13 /// Their integrals on the range [-1, 1] are calculated
14 ///
15 /// \macro_image
16 /// \macro_output
17 /// \macro_code
18 ///
19 /// \author Magdalena Slawinska
20 
21 #if defined(__CINT__) && !defined(__MAKECINT__)
22 {
23  gSystem->CompileMacro("LegendreAssoc.C", "k");
24  LegendreAssoc();
25 }
26 #else
27 
28 #include "TMath.h"
29 #include "TF1.h"
30 #include "TCanvas.h"
31 
32 #include <Riostream.h>
33 #include "TLegend.h"
34 #include "TLegendEntry.h"
35 
36 #include "Math/IFunction.h"
37 #include <cmath>
38 #include "TSystem.h"
39 
40 void LegendreAssoc()
41 {
42  gSystem->Load("libMathMore");
43 
44  std::cout <<"Drawing associate Legendre Polynomials.." << std::endl;
45  TCanvas *Canvas = new TCanvas("DistCanvas", "Associate Legendre polynomials", 10, 10, 800, 500);
46  Canvas->SetFillColor(17);
47  Canvas->Divide(2,1);
48  Canvas->SetFrameFillColor(19);
49  TLegend *leg1 = new TLegend(0.5, 0.7, 0.8, 0.89);
50  TLegend *leg2 = new TLegend(0.5, 0.7, 0.8, 0.89);
51 
52  //+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
53  //drawing the set of Legendre functions
54  TF1* L[5];
55 
56  L[0]= new TF1("L_0", "ROOT::Math::assoc_legendre(1, 0,x)", -1, 1);
57  L[1]= new TF1("L_1", "ROOT::Math::assoc_legendre(1, 1,x)", -1, 1);
58  L[2]= new TF1("L_2", "ROOT::Math::assoc_legendre(2, 0,x)", -1, 1);
59  L[3]= new TF1("L_3", "ROOT::Math::assoc_legendre(2, 1,x)", -1, 1);
60  L[4]= new TF1("L_4", "ROOT::Math::assoc_legendre(2, 2,x)", -1, 1);
61 
62  TF1* SL[5];
63 
64  SL[0]= new TF1("SL_0", "ROOT::Math::sph_legendre(1, 0,x)", -TMath::Pi(), TMath::Pi());
65  SL[1]= new TF1("SL_1", "ROOT::Math::sph_legendre(1, 1,x)", -TMath::Pi(), TMath::Pi());
66  SL[2]= new TF1("SL_2", "ROOT::Math::sph_legendre(2, 0,x)", -TMath::Pi(), TMath::Pi());
67  SL[3]= new TF1("SL_3", "ROOT::Math::sph_legendre(2, 1,x)", -TMath::Pi(), TMath::Pi());
68  SL[4]= new TF1("SL_4", "ROOT::Math::sph_legendre(2, 2,x)", -TMath::Pi(), TMath::Pi() );
69 
70  Canvas->cd(1);
71  gPad->SetGrid();
72  gPad->SetFillColor(kWhite);
73  L[0]->SetMaximum(3);
74  L[0]->SetMinimum(-2);
75  L[0]->SetTitle("Associate Legendre Polynomials");
76  for (int nu = 0; nu < 5; nu++) {
77  L[nu]->SetLineStyle(1);
78  L[nu]->SetLineWidth(2);
79  L[nu]->SetLineColor(nu+1);
80  }
81 
82  leg1->AddEntry(L[0]->DrawCopy(), " P^{1}_{0}(x)", "l");
83  leg1->AddEntry(L[1]->DrawCopy("same"), " P^{1}_{1}(x)", "l");
84  leg1->AddEntry(L[2]->DrawCopy("same"), " P^{2}_{0}(x)", "l");
85  leg1->AddEntry(L[3]->DrawCopy("same"), " P^{2}_{1}(x)", "l");
86  leg1->AddEntry(L[4]->DrawCopy("same"), " P^{2}_{2}(x)", "l");
87  leg1->Draw();
88 
89  Canvas->cd(2);
90  gPad->SetGrid();
91  gPad->SetFillColor(kWhite);
92  SL[0]->SetMaximum(1);
93  SL[0]->SetMinimum(-1);
94  SL[0]->SetTitle("Spherical Legendre Polynomials");
95  for (int nu = 0; nu < 5; nu++) {
96  SL[nu]->SetLineStyle(1);
97  SL[nu]->SetLineWidth(2);
98  SL[nu]->SetLineColor(nu+1);
99  }
100 
101  leg2->AddEntry(SL[0]->DrawCopy(), " P^{1}_{0}(x)", "l");
102  leg2->AddEntry(SL[1]->DrawCopy("same"), " P^{1}_{1}(x)", "l");
103  leg2->AddEntry(SL[2]->DrawCopy("same"), " P^{2}_{0}(x)", "l");
104  leg2->AddEntry(SL[3]->DrawCopy("same"), " P^{2}_{1}(x)", "l");
105  leg2->AddEntry(SL[4]->DrawCopy("same"), " P^{2}_{2}(x)", "l");
106  leg2->Draw();
107 
108 
109  //integration
110 
111  std::cout << "Calculating integrals of Associate Legendre Polynomials on [-1, 1]" << std::endl;
112  double integral[5];
113  for (int nu = 0; nu < 5; nu++) {
114  integral[nu] = L[nu]->Integral(-1.0, 1.0);
115  std::cout <<"Integral [-1,1] for Associated Legendre Polynomial of Degree " << nu << "\t = \t" << integral[nu] << std::endl;
116  }
117 }
118 
119 #endif
120 
virtual void SetLineWidth(Width_t lwidth)
Definition: TAttLine.h:57
This class displays a legend box (TPaveText) containing several legend entries.
Definition: TLegend.h:35
virtual void Draw(Option_t *option="")
Draw this legend with its current attributes.
Definition: TLegend.cxx:373
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:659
virtual int Load(const char *module, const char *entry="", Bool_t system=kFALSE)
Load a shared library.
Definition: TSystem.cxx:1766
virtual void SetMinimum(Double_t minimum=-1111)
Set the minimum value along Y for this function In case the function is already drawn, set also the minimum in the helper histogram.
Definition: TF1.cxx:3090
virtual Double_t Integral(Double_t a, Double_t b, Double_t epsrel=1.e-12)
IntegralOneDim or analytical integral.
Definition: TF1.cxx:2241
virtual int CompileMacro(const char *filename, Option_t *opt="", const char *library_name="", const char *build_dir="", UInt_t dirmode=0)
This method compiles and loads a shared library containing the code from the file "filename"...
Definition: TSystem.cxx:2736
Definition: Rtypes.h:60
virtual void SetLineColor(Color_t lcolor)
Definition: TAttLine.h:54
R__EXTERN TSystem * gSystem
Definition: TSystem.h:545
virtual void SetFillColor(Color_t fcolor)
Definition: TAttFill.h:50
Double_t Pi()
Definition: TMath.h:44
virtual void SetMaximum(Double_t maximum=-1111)
Set the maximum value along Y for this function In case the function is already drawn, set also the maximum in the helper histogram.
Definition: TF1.cxx:3077
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:3227
The Canvas class.
Definition: TCanvas.h:48
TLegendEntry * AddEntry(const TObject *obj, const char *label="", Option_t *option="lpf")
Add a new entry to this legend.
Definition: TLegend.cxx:280
virtual void SetLineStyle(Style_t lstyle)
Definition: TAttLine.h:56
virtual void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0)
Automatic pad generation by division.
Definition: TPad.cxx:1073
1-Dim function class
Definition: TF1.h:149
#define gPad
Definition: TVirtualPad.h:288
void SetFrameFillColor(Color_t color=1)
Definition: TAttPad.h:83