ROOT logo
// Example describing the usage of different kinds of Associate Legendre Polynomials
// To execute the macro type in:
//
// root[0]: .x LegendreAssoc.C 
//
// It draws common graphs for first 5
// Associate Legendre Polynomials
// and Spherical Associate Legendre Polynomials
// Their integrals on the range [-1, 1] are calculated
//
// Author: Magdalena Slawinska

#if defined(__CINT__) && !defined(__MAKECINT__)
{
   gSystem->CompileMacro("LegendreAssoc.C", "k");
   LegendreAssoc();
}
#else

#include "TMath.h"
#include "TF1.h"
#include "TCanvas.h"

#include <Riostream.h>
#include "TLegend.h"
#include "TLegendEntry.h"

#include "Math/IFunction.h"
#include <cmath>
#include "TSystem.h"

void LegendreAssoc()
{
  
  //const int n=5;
  gSystem->Load("libMathMore");

  std::cout <<"Drawing associate Legendre Polynomials.." << std::endl;
  TCanvas *Canvas = new TCanvas("DistCanvas", "Associate Legendre polynomials", 10, 10, 1000, 600);  
  Canvas->SetFillColor(17);
  Canvas->Divide(2,1);
  Canvas->SetFrameFillColor(19);
  TLegend *leg1 = new TLegend(0.5, 0.7, 0.8, 0.89); 
  TLegend *leg2 = new TLegend(0.5, 0.7, 0.8, 0.89); 
  //leg->TLegend::SetNDC();
//+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
//drawing the set of Legendre functions   
  TF1* L[5];

    L[0]= new TF1("L_0", "ROOT::Math::assoc_legendre(1, 0,x)", -1, 1);
    L[1]= new TF1("L_1", "ROOT::Math::assoc_legendre(1, 1,x)", -1, 1);
    L[2]= new TF1("L_2", "ROOT::Math::assoc_legendre(2, 0,x)", -1, 1);
    L[3]= new TF1("L_3", "ROOT::Math::assoc_legendre(2, 1,x)", -1, 1);
    L[4]= new TF1("L_4", "ROOT::Math::assoc_legendre(2, 2,x)", -1, 1);


  TF1* SL[5];

    SL[0]= new TF1("SL_0", "ROOT::Math::sph_legendre(1, 0,x)", -TMath::Pi(), TMath::Pi());
    SL[1]= new TF1("SL_1", "ROOT::Math::sph_legendre(1, 1,x)", -TMath::Pi(), TMath::Pi());
    SL[2]= new TF1("SL_2", "ROOT::Math::sph_legendre(2, 0,x)", -TMath::Pi(), TMath::Pi());
    SL[3]= new TF1("SL_3", "ROOT::Math::sph_legendre(2, 1,x)", -TMath::Pi(), TMath::Pi());
    SL[4]= new TF1("SL_4", "ROOT::Math::sph_legendre(2, 2,x)", -TMath::Pi(), TMath::Pi() );


    Canvas->cd(1);
    gPad->SetGrid();
    gPad->SetFillColor(kWhite);
    L[0]->SetMaximum(3);
    L[0]->SetMinimum(-2);
    L[0]->SetTitle("Associate Legendre Polynomials");  
    for(int nu = 0; nu < 5; nu++)
    {
      //L[nu]->SetTitle("Legendre polynomials");  

      L[nu]->SetLineStyle(1);
      L[nu]->SetLineWidth(2);
      L[nu]->SetLineColor(nu+1);
    }

    leg1->AddEntry(L[0]->DrawCopy(), " P^{1}_{0}(x)", "l");
    leg1->AddEntry(L[1]->DrawCopy("same"), " P^{1}_{1}(x)", "l");
    leg1->AddEntry(L[2]->DrawCopy("same"), " P^{2}_{0}(x)", "l");
    leg1->AddEntry(L[3]->DrawCopy("same"), " P^{2}_{1}(x)", "l");
    leg1->AddEntry(L[4]->DrawCopy("same"), " P^{2}_{2}(x)", "l");
    leg1->Draw();

    Canvas->cd(2);
    gPad->SetGrid();
    gPad->SetFillColor(kWhite);
    SL[0]->SetMaximum(1);
    SL[0]->SetMinimum(-1);
    SL[0]->SetTitle("Spherical Legendre Polynomials");  
    for(int nu = 0; nu < 5; nu++)
    {
      //L[nu]->SetTitle("Legendre polynomials");  

      SL[nu]->SetLineStyle(1);
      SL[nu]->SetLineWidth(2);
      SL[nu]->SetLineColor(nu+1);
    }

    leg2->AddEntry(SL[0]->DrawCopy(), " P^{1}_{0}(x)", "l");
    leg2->AddEntry(SL[1]->DrawCopy("same"), " P^{1}_{1}(x)", "l");
    leg2->AddEntry(SL[2]->DrawCopy("same"), " P^{2}_{0}(x)", "l");
    leg2->AddEntry(SL[3]->DrawCopy("same"), " P^{2}_{1}(x)", "l");
    leg2->AddEntry(SL[4]->DrawCopy("same"), " P^{2}_{2}(x)", "l");
    leg2->Draw();


    //integration

    std::cout << "Calculating integrals of Associate Legendre Polynomials on [-1, 1]" << std::endl;
    double integral[5];
    for(int nu = 0; nu < 5; nu++)
    {
      integral[nu] = L[nu]->Integral(-1.0, 1.0);
      std::cout <<"Integral [-1,1] for Associated Legendre Polynomial of Degree " << nu << "\t = \t" << integral[nu] <<  std::endl;
    }


  
}

#endif

 LegendreAssoc.C:1
 LegendreAssoc.C:2
 LegendreAssoc.C:3
 LegendreAssoc.C:4
 LegendreAssoc.C:5
 LegendreAssoc.C:6
 LegendreAssoc.C:7
 LegendreAssoc.C:8
 LegendreAssoc.C:9
 LegendreAssoc.C:10
 LegendreAssoc.C:11
 LegendreAssoc.C:12
 LegendreAssoc.C:13
 LegendreAssoc.C:14
 LegendreAssoc.C:15
 LegendreAssoc.C:16
 LegendreAssoc.C:17
 LegendreAssoc.C:18
 LegendreAssoc.C:19
 LegendreAssoc.C:20
 LegendreAssoc.C:21
 LegendreAssoc.C:22
 LegendreAssoc.C:23
 LegendreAssoc.C:24
 LegendreAssoc.C:25
 LegendreAssoc.C:26
 LegendreAssoc.C:27
 LegendreAssoc.C:28
 LegendreAssoc.C:29
 LegendreAssoc.C:30
 LegendreAssoc.C:31
 LegendreAssoc.C:32
 LegendreAssoc.C:33
 LegendreAssoc.C:34
 LegendreAssoc.C:35
 LegendreAssoc.C:36
 LegendreAssoc.C:37
 LegendreAssoc.C:38
 LegendreAssoc.C:39
 LegendreAssoc.C:40
 LegendreAssoc.C:41
 LegendreAssoc.C:42
 LegendreAssoc.C:43
 LegendreAssoc.C:44
 LegendreAssoc.C:45
 LegendreAssoc.C:46
 LegendreAssoc.C:47
 LegendreAssoc.C:48
 LegendreAssoc.C:49
 LegendreAssoc.C:50
 LegendreAssoc.C:51
 LegendreAssoc.C:52
 LegendreAssoc.C:53
 LegendreAssoc.C:54
 LegendreAssoc.C:55
 LegendreAssoc.C:56
 LegendreAssoc.C:57
 LegendreAssoc.C:58
 LegendreAssoc.C:59
 LegendreAssoc.C:60
 LegendreAssoc.C:61
 LegendreAssoc.C:62
 LegendreAssoc.C:63
 LegendreAssoc.C:64
 LegendreAssoc.C:65
 LegendreAssoc.C:66
 LegendreAssoc.C:67
 LegendreAssoc.C:68
 LegendreAssoc.C:69
 LegendreAssoc.C:70
 LegendreAssoc.C:71
 LegendreAssoc.C:72
 LegendreAssoc.C:73
 LegendreAssoc.C:74
 LegendreAssoc.C:75
 LegendreAssoc.C:76
 LegendreAssoc.C:77
 LegendreAssoc.C:78
 LegendreAssoc.C:79
 LegendreAssoc.C:80
 LegendreAssoc.C:81
 LegendreAssoc.C:82
 LegendreAssoc.C:83
 LegendreAssoc.C:84
 LegendreAssoc.C:85
 LegendreAssoc.C:86
 LegendreAssoc.C:87
 LegendreAssoc.C:88
 LegendreAssoc.C:89
 LegendreAssoc.C:90
 LegendreAssoc.C:91
 LegendreAssoc.C:92
 LegendreAssoc.C:93
 LegendreAssoc.C:94
 LegendreAssoc.C:95
 LegendreAssoc.C:96
 LegendreAssoc.C:97
 LegendreAssoc.C:98
 LegendreAssoc.C:99
 LegendreAssoc.C:100
 LegendreAssoc.C:101
 LegendreAssoc.C:102
 LegendreAssoc.C:103
 LegendreAssoc.C:104
 LegendreAssoc.C:105
 LegendreAssoc.C:106
 LegendreAssoc.C:107
 LegendreAssoc.C:108
 LegendreAssoc.C:109
 LegendreAssoc.C:110
 LegendreAssoc.C:111
 LegendreAssoc.C:112
 LegendreAssoc.C:113
 LegendreAssoc.C:114
 LegendreAssoc.C:115
 LegendreAssoc.C:116
 LegendreAssoc.C:117
 LegendreAssoc.C:118
 LegendreAssoc.C:119
 LegendreAssoc.C:120
 LegendreAssoc.C:121
 LegendreAssoc.C:122
 LegendreAssoc.C:123
 LegendreAssoc.C:124
 LegendreAssoc.C:125
 LegendreAssoc.C:126
 LegendreAssoc.C:127