Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
pdf009_Bessel.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_pdf
3## \notebook
4## Show the different kinds of Bessel functions available in ROOT
5## To execute the macro type in:
6##
7## ~~~{.cpp}
8## root[0] .x Bessel.C
9## ~~~
10##
11## It will create one canvas with the representation
12## of the cylindrical and spherical Bessel functions
13## regular and modified
14##
15## Based on Bessel.C by Magdalena Slawinska
16##
17## \macro_image
18## \macro_code
19##
20## \author Juan Fernando Jaramillo Botero
21
22import ROOT
23from ROOT import TCanvas, TF1, gSystem, gPad, TLegend, TPaveLabel, kBlack
24
25
26gSystem.Load("libMathMore")
27
28DistCanvas = TCanvas("DistCanvas", "Bessel functions example", 10, 10, 800, 600)
34leg = TLegend(0.75, 0.7, 0.89, 0.89)
35
36# Drawing the set of Bessel J functions
37#
38# n is the number of functions in each pad
39n = 5
40JBessel = []
41for nu in range(n):
42 jbessel = TF1("J_0", "ROOT::Math::cyl_bessel_j([0],x)", 0, 10)
48 JBessel.append(jbessel)
49
50# Setting x axis for JBessel
51xaxis = JBessel[0].GetXaxis()
55
56# setting the title in a label style
57p1 = TPaveLabel(.0, .90, .0 + .50, .90 + .10, "Bessel J functions", "NDC")
60p1.SetTextColor(kBlack)
61
62# setting the legend
63leg.AddEntry(JBessel[0].DrawCopy(), " J_0(x)", "l")
64leg.AddEntry(JBessel[1].DrawCopy("same"), " J_1(x)", "l")
65leg.AddEntry(JBessel[2].DrawCopy("same"), " J_2(x)", "l")
66leg.AddEntry(JBessel[3].DrawCopy("same"), " J_3(x)", "l")
67leg.AddEntry(JBessel[4].DrawCopy("same"), " J_4(x)", "l")
68
70p1.Draw()
71
72# Set canvas 2
76leg2 = TLegend(0.75, 0.7, 0.89, 0.89)
77
78# Drawing Bessel k
79KBessel = []
80for nu in range(n):
81 kbessel = TF1("J_0", "ROOT::Math::cyl_bessel_k([0],x)", 0, 10)
83 kbessel.SetTitle("Bessel K functions")
87 KBessel.append(kbessel)
88kxaxis = KBessel[0].GetXaxis()
92
93# setting title
94p2 = TPaveLabel(.0, .90, .0 + .50, .90 + .10, "Bessel K functions", "NDC")
97p2.SetTextColor(kBlack)
98
99# setting legend
100leg2.AddEntry(KBessel[0].DrawCopy(), " K_0(x)", "l")
101leg2.AddEntry(KBessel[1].DrawCopy("same"), " K_1(x)", "l")
102leg2.AddEntry(KBessel[2].DrawCopy("same"), " K_2(x)", "l")
103leg2.AddEntry(KBessel[3].DrawCopy("same"), " K_3(x)", "l")
104leg2.AddEntry(KBessel[4].DrawCopy("same"), " K_4(x)", "l")
105leg2.Draw()
106p2.Draw()
107
108# Set canvas 3
112leg3 = TLegend(0.75, 0.7, 0.89, 0.89)
113
114# Drawing Bessel i
115iBessel = []
116for nu in range(5):
117 ibessel = TF1("J_0", "ROOT::Math::cyl_bessel_i([0],x)", 0, 10)
118 ibessel.SetParameters(nu, 0.0)
119 ibessel.SetTitle("Bessel I functions")
123 iBessel.append(ibessel)
124
125iaxis = iBessel[0].GetXaxis()
129
130# setting title
131p3 = TPaveLabel(.0, .90, .0 + .50, .90 + .10, "Bessel I functions", "NDC")
134p3.SetTextColor(kBlack)
135
136# setting legend
137leg3.AddEntry(iBessel[0].DrawCopy(), " I_0", "l")
138leg3.AddEntry(iBessel[1].DrawCopy("same"), " I_1(x)", "l")
139leg3.AddEntry(iBessel[2].DrawCopy("same"), " I_2(x)", "l")
140leg3.AddEntry(iBessel[3].DrawCopy("same"), " I_3(x)", "l")
141leg3.AddEntry(iBessel[4].DrawCopy("same"), " I_4(x)", "l")
142leg3.Draw()
143p3.Draw()
144
145# Set canvas 4
149leg4 = TLegend(0.75, 0.7, 0.89, 0.89)
150
151# Drawing sph_bessel
152jBessel = []
153for nu in range(5):
154 jbessel = TF1("J_0", "ROOT::Math::sph_bessel([0],x)", 0, 10)
155 jbessel.SetParameters(nu, 0.0)
156 jbessel.SetTitle("Bessel j functions")
160 jBessel.append(jbessel)
161jaxis = jBessel[0].GetXaxis()
165
166# setting title
167p4 = TPaveLabel(.0, .90, .0 + .50, .90 + .10, "Bessel j functions", "NDC")
170p4.SetTextColor(kBlack)
171
172# setting legend
173leg4.AddEntry(jBessel[0].DrawCopy(), " j_0(x)", "l")
174leg4.AddEntry(jBessel[1].DrawCopy("same"), " j_1(x)", "l")
175leg4.AddEntry(jBessel[2].DrawCopy("same"), " j_2(x)", "l")
176leg4.AddEntry(jBessel[3].DrawCopy("same"), " j_3(x)", "l")
177leg4.AddEntry(jBessel[4].DrawCopy("same"), " j_4(x)", "l")
178
179leg4.Draw()
180p4.Draw()
181
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
The Canvas class.
Definition TCanvas.h:23
1-Dim function class
Definition TF1.h:234
This class displays a legend box (TPaveText) containing several legend entries.
Definition TLegend.h:23
A Pave (see TPave) with a text centered in the Pave.
Definition TPaveLabel.h:20