Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
LegendreAssoc.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_math
3/// \notebook
4/// Example describing the usage of different kinds of Associate Legendre Polynomials
5/// To execute the macro type in:
6///
7/// ~~~{.cpp}
8/// root[0] .x LegendreAssoc.C
9/// ~~~
10///
11/// It draws common graphs for first 5
12/// Associate Legendre Polynomials
13/// and Spherical Associate Legendre Polynomials
14/// Their integrals on the range [-1, 1] are calculated
15///
16/// \macro_image
17/// \macro_output
18/// \macro_code
19///
20/// \author Magdalena Slawinska
21
22
23#include "TMath.h"
24#include "TF1.h"
25#include "TCanvas.h"
26
27#include <Riostream.h>
28#include "TLegend.h"
29#include "TLegendEntry.h"
30
31#include "Math/IFunction.h"
32#include <cmath>
33#include "TSystem.h"
34
35void LegendreAssoc()
36{
37 std::cout <<"Drawing associate Legendre Polynomials.." << std::endl;
38 TCanvas *Canvas = new TCanvas("DistCanvas", "Associate Legendre polynomials", 10, 10, 800, 500);
39 Canvas->Divide(2,1);
40 TLegend *leg1 = new TLegend(0.5, 0.7, 0.8, 0.89);
41 TLegend *leg2 = new TLegend(0.5, 0.7, 0.8, 0.89);
42
43 //-------------------------------------------
44 //drawing the set of Legendre functions
45 TF1* L[5];
46
47 L[0]= new TF1("L_0", "ROOT::Math::assoc_legendre(1, 0,x)", -1, 1);
48 L[1]= new TF1("L_1", "ROOT::Math::assoc_legendre(1, 1,x)", -1, 1);
49 L[2]= new TF1("L_2", "ROOT::Math::assoc_legendre(2, 0,x)", -1, 1);
50 L[3]= new TF1("L_3", "ROOT::Math::assoc_legendre(2, 1,x)", -1, 1);
51 L[4]= new TF1("L_4", "ROOT::Math::assoc_legendre(2, 2,x)", -1, 1);
52
53 TF1* SL[5];
54
55 SL[0]= new TF1("SL_0", "ROOT::Math::sph_legendre(1, 0,x)", -TMath::Pi(), TMath::Pi());
56 SL[1]= new TF1("SL_1", "ROOT::Math::sph_legendre(1, 1,x)", -TMath::Pi(), TMath::Pi());
57 SL[2]= new TF1("SL_2", "ROOT::Math::sph_legendre(2, 0,x)", -TMath::Pi(), TMath::Pi());
58 SL[3]= new TF1("SL_3", "ROOT::Math::sph_legendre(2, 1,x)", -TMath::Pi(), TMath::Pi());
59 SL[4]= new TF1("SL_4", "ROOT::Math::sph_legendre(2, 2,x)", -TMath::Pi(), TMath::Pi() );
60
61 Canvas->cd(1);
62 gPad->SetGrid();
63 gPad->SetFillColor(kWhite);
64 L[0]->SetMaximum(3);
65 L[0]->SetMinimum(-2);
66 L[0]->SetTitle("Associate Legendre Polynomials");
67 for (int nu = 0; nu < 5; nu++) {
68 L[nu]->SetLineStyle(1);
69 L[nu]->SetLineWidth(2);
70 L[nu]->SetLineColor(nu+1);
71 }
72
73 leg1->AddEntry(L[0]->DrawCopy(), " P^{1}_{0}(x)", "l");
74 leg1->AddEntry(L[1]->DrawCopy("same"), " P^{1}_{1}(x)", "l");
75 leg1->AddEntry(L[2]->DrawCopy("same"), " P^{2}_{0}(x)", "l");
76 leg1->AddEntry(L[3]->DrawCopy("same"), " P^{2}_{1}(x)", "l");
77 leg1->AddEntry(L[4]->DrawCopy("same"), " P^{2}_{2}(x)", "l");
78 leg1->Draw();
79
80 Canvas->cd(2);
81 gPad->SetGrid();
82 gPad->SetFillColor(kWhite);
83 SL[0]->SetMaximum(1);
84 SL[0]->SetMinimum(-1);
85 SL[0]->SetTitle("Spherical Legendre Polynomials");
86 for (int nu = 0; nu < 5; nu++) {
87 SL[nu]->SetLineStyle(1);
88 SL[nu]->SetLineWidth(2);
89 SL[nu]->SetLineColor(nu+1);
90 }
91
92 leg2->AddEntry(SL[0]->DrawCopy(), " P^{1}_{0}(x)", "l");
93 leg2->AddEntry(SL[1]->DrawCopy("same"), " P^{1}_{1}(x)", "l");
94 leg2->AddEntry(SL[2]->DrawCopy("same"), " P^{2}_{0}(x)", "l");
95 leg2->AddEntry(SL[3]->DrawCopy("same"), " P^{2}_{1}(x)", "l");
96 leg2->AddEntry(SL[4]->DrawCopy("same"), " P^{2}_{2}(x)", "l");
97 leg2->Draw();
98
99
100 //integration
101
102 std::cout << "Calculating integrals of Associate Legendre Polynomials on [-1, 1]" << std::endl;
103 double integral[5];
104 for (int nu = 0; nu < 5; nu++) {
105 integral[nu] = L[nu]->Integral(-1.0, 1.0);
106 std::cout <<"Integral [-1,1] for Associated Legendre Polynomial of Degree " << nu << "\t = \t" << integral[nu] << std::endl;
107 }
108}
@ kWhite
Definition Rtypes.h:65
#define gPad
virtual void SetLineStyle(Style_t lstyle)
Set the line style.
Definition TAttLine.h:42
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition TAttLine.h:43
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:40
The Canvas class.
Definition TCanvas.h:23
TVirtualPad * cd(Int_t subpadnumber=0) override
Set current canvas & pad.
Definition TCanvas.cxx:716
1-Dim function class
Definition TF1.h:233
virtual void SetMaximum(Double_t maximum=-1111)
Set the maximum value along Y for this function In case the function is already drawn,...
Definition TF1.cxx:3394
void SetTitle(const char *title="") override
Set function title if title has the form "fffffff;xxxx;yyyy", it is assumed that the function title i...
Definition TF1.cxx:3558
virtual void SetMinimum(Double_t minimum=-1111)
Set the minimum value along Y for this function In case the function is already drawn,...
Definition TF1.cxx:3407
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:317
void Draw(Option_t *option="") override
Draw this legend with its current attributes.
Definition TLegend.cxx:422
void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0) override
Automatic pad generation by division.
Definition TPad.cxx:1196
virtual void SetGrid(Int_t valuex=1, Int_t valuey=1)=0
RooArgList L(Args_t &&... args)
Definition RooArgList.h:156
constexpr Double_t Pi()
Definition TMath.h:37