00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013 #if defined(__CINT__) && !defined(__MAKECINT__)
00014 {
00015 gSystem->CompileMacro("LegendreAssoc.C", "k");
00016 LegendreAssoc();
00017 }
00018 #else
00019
00020 #include "TMath.h"
00021 #include "TF1.h"
00022 #include "TCanvas.h"
00023
00024 #include <Riostream.h>
00025 #include "TLegend.h"
00026 #include "TLegendEntry.h"
00027
00028 #include "Math/IFunction.h"
00029 #include <cmath>
00030 #include "TSystem.h"
00031
00032 void LegendreAssoc()
00033 {
00034
00035
00036 gSystem->Load("libMathMore");
00037
00038 std::cout <<"Drawing associate Legendre Polynomials.." << std::endl;
00039 TCanvas *Canvas = new TCanvas("DistCanvas", "Associate Legendre polynomials", 10, 10, 1000, 600);
00040 Canvas->SetFillColor(17);
00041 Canvas->Divide(2,1);
00042 Canvas->SetFrameFillColor(19);
00043 TLegend *leg1 = new TLegend(0.5, 0.7, 0.8, 0.89);
00044 TLegend *leg2 = new TLegend(0.5, 0.7, 0.8, 0.89);
00045
00046
00047
00048 TF1* L[5];
00049
00050 L[0]= new TF1("L_0", "ROOT::Math::assoc_legendre(1, 0,x)", -1, 1);
00051 L[1]= new TF1("L_1", "ROOT::Math::assoc_legendre(1, 1,x)", -1, 1);
00052 L[2]= new TF1("L_2", "ROOT::Math::assoc_legendre(2, 0,x)", -1, 1);
00053 L[3]= new TF1("L_3", "ROOT::Math::assoc_legendre(2, 1,x)", -1, 1);
00054 L[4]= new TF1("L_4", "ROOT::Math::assoc_legendre(2, 2,x)", -1, 1);
00055
00056
00057 TF1* SL[5];
00058
00059 SL[0]= new TF1("SL_0", "ROOT::Math::sph_legendre(1, 0,x)", -TMath::Pi(), TMath::Pi());
00060 SL[1]= new TF1("SL_1", "ROOT::Math::sph_legendre(1, 1,x)", -TMath::Pi(), TMath::Pi());
00061 SL[2]= new TF1("SL_2", "ROOT::Math::sph_legendre(2, 0,x)", -TMath::Pi(), TMath::Pi());
00062 SL[3]= new TF1("SL_3", "ROOT::Math::sph_legendre(2, 1,x)", -TMath::Pi(), TMath::Pi());
00063 SL[4]= new TF1("SL_4", "ROOT::Math::sph_legendre(2, 2,x)", -TMath::Pi(), TMath::Pi() );
00064
00065
00066 Canvas->cd(1);
00067 gPad->SetGrid();
00068 gPad->SetFillColor(kWhite);
00069 L[0]->SetMaximum(3);
00070 L[0]->SetMinimum(-2);
00071 L[0]->SetTitle("Associate Legendre Polynomials");
00072 for(int nu = 0; nu < 5; nu++)
00073 {
00074
00075
00076 L[nu]->SetLineStyle(1);
00077 L[nu]->SetLineWidth(2);
00078 L[nu]->SetLineColor(nu+1);
00079 }
00080
00081 leg1->AddEntry(L[0]->DrawCopy(), " P^{1}_{0}(x)", "l");
00082 leg1->AddEntry(L[1]->DrawCopy("same"), " P^{1}_{1}(x)", "l");
00083 leg1->AddEntry(L[2]->DrawCopy("same"), " P^{2}_{0}(x)", "l");
00084 leg1->AddEntry(L[3]->DrawCopy("same"), " P^{2}_{1}(x)", "l");
00085 leg1->AddEntry(L[4]->DrawCopy("same"), " P^{2}_{2}(x)", "l");
00086 leg1->Draw();
00087
00088 Canvas->cd(2);
00089 gPad->SetGrid();
00090 gPad->SetFillColor(kWhite);
00091 SL[0]->SetMaximum(1);
00092 SL[0]->SetMinimum(-1);
00093 SL[0]->SetTitle("Spherical Legendre Polynomials");
00094 for(int nu = 0; nu < 5; nu++)
00095 {
00096
00097
00098 SL[nu]->SetLineStyle(1);
00099 SL[nu]->SetLineWidth(2);
00100 SL[nu]->SetLineColor(nu+1);
00101 }
00102
00103 leg2->AddEntry(SL[0]->DrawCopy(), " P^{1}_{0}(x)", "l");
00104 leg2->AddEntry(SL[1]->DrawCopy("same"), " P^{1}_{1}(x)", "l");
00105 leg2->AddEntry(SL[2]->DrawCopy("same"), " P^{2}_{0}(x)", "l");
00106 leg2->AddEntry(SL[3]->DrawCopy("same"), " P^{2}_{1}(x)", "l");
00107 leg2->AddEntry(SL[4]->DrawCopy("same"), " P^{2}_{2}(x)", "l");
00108 leg2->Draw();
00109
00110
00111
00112
00113 std::cout << "Calculating integrals of Associate Legendre Polynomials on [-1, 1]" << std::endl;
00114 double integral[5];
00115 for(int nu = 0; nu < 5; nu++)
00116 {
00117 integral[nu] = L[nu]->Integral(-1.0, 1.0);
00118 std::cout <<"Integral [-1,1] for Associated Legendre Polynomial of Degree " << nu << "\t = \t" << integral[nu] << std::endl;
00119 }
00120
00121
00122
00123 }
00124
00125 #endif
00126