ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
zdemo.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_graphs
3 /// This macro is an example of graphs in log scales with annotations.
4 ///
5 /// The presented results are predictions of invariant cross-section
6 /// of Direct Photons produced at RHIC energies, based on the universality of
7 /// scaling function H(z).
8 ///
9 /// These Figures were published in JINR preprint E2-98-64, Dubna,
10 /// 1998 and submitted to CPC.
11 ///
12 /// Note that the way greek symbols, super/subscripts are obtained
13 /// illustrate the current limitations of Root in this area.
14 ///
15 /// \macro_image
16 /// \macro_code
17 ///
18 /// \authors Michael Tokarev and Elena Potrebenikova (JINR Dubna)
19 
20 #include "TCanvas.h"
21 #include "TPad.h"
22 #include "TPaveLabel.h"
23 #include "TLatex.h"
24 #include "TGraph.h"
25 #include "TFrame.h"
26 
27 const Int_t NMAX = 20;
28 Int_t NLOOP;
30 
32 
33 //__________________________________________________________________
34 void zdemo()
35 {
36 
37  Float_t energ;
38  Float_t dens;
39  Float_t tgrad;
40  Float_t ptmin;
41  Float_t ptmax;
42  Float_t delp;
43 
44  // Create a new canvas.
45  TCanvas *c1 = new TCanvas("zdemo",
46  "Monte Carlo Study of Z scaling",10,40,800,600);
47  c1->Range(0,0,25,18);
48  c1->SetFillColor(40);
49 
50  TPaveLabel *pl = new TPaveLabel(1,16.3,24,17.5,"Z-scaling of \
51  Direct Photon Productions in pp Collisions at RHIC Energies","br");
52  pl->SetFillColor(18);
53  pl->SetTextFont(32);
54  pl->SetTextColor(49);
55  pl->Draw();
56 
57  TLatex *t = new TLatex();
58  t->SetTextFont(32);
59  t->SetTextColor(1);
60  t->SetTextSize(0.03);
61  t->SetTextAlign(12);
62  t->DrawLatex(3.1,15.5,"M.Tokarev, E.Potrebenikova ");
63  t->DrawLatex(14.,15.5,"JINR preprint E2-98-64, Dubna, 1998 ");
64 
65  TPad *pad1 = new TPad("pad1","This is pad1",0.02,0.02,0.48,0.83,33);
66  TPad *pad2 = new TPad("pad2","This is pad2",0.52,0.02,0.98,0.83,33);
67 
68  pad1->Draw();
69  pad2->Draw();
70 
71 //
72 // Cross-section of direct photon production in pp collisions
73 // at 500 GeV vs Pt
74 //
75  energ = 63;
76  dens = 1.766;
77  tgrad = 90.;
78  ptmin = 4.;
79  ptmax = 24.;
80  delp = 2.;
81  hz_calc(energ, dens, tgrad, ptmin, ptmax, delp);
82  pad1->cd();
83  pad1->Range(-0.255174,-19.25,2.29657,-6.75);
84  pad1->SetLogx();
85  pad1->SetLogy();
86 
87  // create a 2-d histogram to define the range
88  pad1->DrawFrame(1,1e-18,110,1e-8);
89  pad1->GetFrame()->SetFillColor(19);
90  t = new TLatex();
91  t->SetNDC();
92  t->SetTextFont(62);
93  t->SetTextColor(36);
94  t->SetTextSize(0.08);
95  t->SetTextAlign(12);
96  t->DrawLatex(0.6,0.85,"p - p");
97 
98  t->SetTextSize(0.05);
99  t->DrawLatex(0.6,0.79,"Direct #gamma");
100  t->DrawLatex(0.6,0.75,"#theta = 90^{o}");
101 
102  t->DrawLatex(0.20,0.45,"Ed^{3}#sigma/dq^{3}");
103  t->DrawLatex(0.18,0.40,"(barn/Gev^{2})");
104 
105  t->SetTextSize(0.045);
106  t->SetTextColor(kBlue);
107  t->DrawLatex(0.22,0.260,"#sqrt{s} = 63(GeV)");
108  t->SetTextColor(kRed);
109  t->DrawLatex(0.22,0.205,"#sqrt{s} = 200(GeV)");
110  t->SetTextColor(6);
111  t->DrawLatex(0.22,0.15,"#sqrt{s} = 500(GeV)");
112 
113  t->SetTextSize(0.05);
114  t->SetTextColor(1);
115  t->DrawLatex(0.6,0.06,"q_{T} (Gev/c)");
116 
117  TGraph *gr1 = new TGraph(NLOOP,PT,INVSIG);
118 
119  gr1->SetLineColor(38);
120  gr1->SetMarkerColor(kBlue);
121  gr1->SetMarkerStyle(21);
122  gr1->SetMarkerSize(1.1);
123  gr1->Draw("LP");
124 
125 //
126 // Cross-section of direct photon production in pp collisions
127 // at 200 GeV vs Pt
128 //
129 
130  energ = 200;
131  dens = 2.25;
132  tgrad = 90.;
133  ptmin = 4.;
134  ptmax = 64.;
135  delp = 6.;
136  hz_calc(energ, dens, tgrad, ptmin, ptmax, delp);
137 
138  TGraph *gr2 = new TGraph(NLOOP,PT,INVSIG);
139  gr2->SetLineColor(38);
140  gr2->SetMarkerColor(kRed);
141  gr2->SetMarkerStyle(29);
142  gr2->SetMarkerSize(1.5);
143  gr2->Draw("LP");
144 
145 //
146 // Cross-section of direct photon production in pp collisions
147 // at 500 GeV vs Pt
148 //
149  energ = 500;
150  dens = 2.73;
151  tgrad = 90.;
152  ptmin = 4.;
153  ptmax = 104.;
154  delp = 10.;
155  hz_calc(energ, dens, tgrad, ptmin, ptmax, delp);
156 
157  TGraph *gr3 = new TGraph(NLOOP,PT,INVSIG);
158 
159  gr3->SetLineColor(38);
160  gr3->SetMarkerColor(6);
161  gr3->SetMarkerStyle(8);
162  gr3->SetMarkerSize(1.1);
163  gr3->Draw("LP");
164 
165  Float_t *dum = 0;
166  TGraph *graph = new TGraph(1,dum,dum);
167  graph->SetMarkerColor(kBlue);
168  graph->SetMarkerStyle(21);
169  graph->SetMarkerSize(1.1);
170  graph->SetPoint(0,1.7,1.e-16);
171  graph->Draw("LP");
172 
173  graph = new TGraph(1,dum,dum);
174  graph->SetMarkerColor(kRed);
175  graph->SetMarkerStyle(29);
176  graph->SetMarkerSize(1.5);
177  graph->SetPoint(0,1.7,2.e-17);
178  graph->Draw("LP");
179 
180  graph = new TGraph(1,dum,dum);
181  graph->SetMarkerColor(6);
182  graph->SetMarkerStyle(8);
183  graph->SetMarkerSize(1.1);
184  graph->SetPoint(0,1.7,4.e-18);
185  graph->Draw("LP");
186 
187  pad2->cd();
188  pad2->Range(-0.43642,-23.75,3.92778,-6.25);
189  pad2->SetLogx();
190  pad2->SetLogy();
191 
192  pad2->DrawFrame(1,1e-22,3100,1e-8);
193  pad2->GetFrame()->SetFillColor(19);
194 
195  TGraph *gr = new TGraph(NLOOP,Z,HZ);
196  gr->SetTitle("HZ vs Z");
197  gr->SetFillColor(19);
198  gr->SetLineColor(9);
199  gr->SetMarkerColor(50);
200  gr->SetMarkerStyle(29);
201  gr->SetMarkerSize(1.5);
202  gr->Draw("LP");
203 
204  t = new TLatex();
205  t->SetNDC();
206  t->SetTextFont(62);
207  t->SetTextColor(36);
208  t->SetTextSize(0.08);
209  t->SetTextAlign(12);
210  t->DrawLatex(0.6,0.85,"p - p");
211 
212  t->SetTextSize(0.05);
213  t->DrawLatex(0.6,0.79,"Direct #gamma");
214  t->DrawLatex(0.6,0.75,"#theta = 90^{o}");
215 
216  t->DrawLatex(0.70,0.55,"H(z)");
217  t->DrawLatex(0.68,0.50,"(barn)");
218 
219  t->SetTextSize(0.045);
220  t->SetTextColor(46);
221  t->DrawLatex(0.20,0.30,"#sqrt{s}, GeV");
222  t->DrawLatex(0.22,0.26,"63");
223  t->DrawLatex(0.22,0.22,"200");
224  t->DrawLatex(0.22,0.18,"500");
225 
226  t->SetTextSize(0.05);
227  t->SetTextColor(1);
228  t->DrawLatex(0.88,0.06,"z");
229 
230  c1->Modified();
231  c1->Update();
232 }
233 
234 void hz_calc(Float_t ENERG, Float_t DENS, Float_t TGRAD, Float_t PTMIN,
235  Float_t PTMAX, Float_t DELP)
236 {
237  Int_t I;
238 
239  Float_t GM1 = 0.00001;
240  Float_t GM2 = 0.00001;
241  Float_t A1 = 1.;
242  Float_t A2 = 1.;
243  Float_t ALX = 2.;
244  Float_t BETA = 1.;
245  Float_t KF1 = 8.E-7;
246  Float_t KF2 = 5.215;
247 
248  Float_t MN = 0.9383;
249  Float_t DEGRAD=0.01745329;
250 
251  Float_t EB1, EB2, PB1, PB2, MB1, MB2, M1, M2;
252  Float_t DNDETA;
253 
254  Float_t P1P2, P1P3, P2P3;
255  Float_t Y1, Y2, S, SMIN, SX1, SX2, SX1X2, DELM;
256  Float_t Y1X1, Y1X2, Y2X1, Y2X2, Y2X1X2, Y1X1X2;
257  Float_t KX1, KX2, ZX1, ZX2;
258  Float_t H1;
259 
260  Float_t PTOT, THET, ETOT, X1, X2;
261 
262  DNDETA= DENS;
263  MB1 = MN*A1;
264  MB2 = MN*A2;
265  EB1 = ENERG/2.*A1;
266  EB2 = ENERG/2.*A2;
267  M1 = GM1;
268  M2 = GM2;
269  THET = TGRAD*DEGRAD;
270  NLOOP = (PTMAX-PTMIN)/DELP;
271 
272  for (I=0; I<NLOOP;I++) {
273  PT[I]=PTMIN+I*DELP;
274  PTOT = PT[I]/sin(THET);
275 
276  ETOT = sqrt(M1*M1 + PTOT*PTOT);
277  PB1 = sqrt(EB1*EB1 - MB1*MB1);
278  PB2 = sqrt(EB2*EB2 - MB2*MB2);
279  P2P3 = EB2*ETOT+PB2*PTOT*cos(THET);
280  P1P2 = EB2*EB1+PB2*PB1;
281  P1P3 = EB1*ETOT-PB1*PTOT*cos(THET);
282 
283  X1 = P2P3/P1P2;
284  X2 = P1P3/P1P2;
285  Y1 = X1+sqrt(X1*X2*(1.-X1)/(1.-X2));
286  Y2 = X2+sqrt(X1*X2*(1.-X2)/(1.-X1));
287 
288  S = (MB1*MB1)+2.*P1P2+(MB2*MB2);
289  SMIN = 4.*((MB1*MB1)*(X1*X1) +2.*X1*X2*P1P2+(MB2*MB2)*(X2*X2));
290  SX1 = 4.*( 2*(MB1*MB1)*X1+2*X2*P1P2);
291  SX2 = 4.*( 2*(MB2*MB2)*X2+2*X1*P1P2);
292  SX1X2= 4.*(2*P1P2);
293  DELM = pow((1.-Y1)*(1.-Y2),ALX);
294 
295  Z[I] = sqrt(SMIN)/DELM/pow(DNDETA,BETA);
296 
297  Y1X1 = 1. +X2*(1-2.*X1)/(2.*(Y1-X1)*(1.-X2));
298  Y1X2 = X1*(1-X1)/(2.*(Y1-X1)*(1.-X2)*(1.-X2));
299  Y2X1 = X2*(1-X2)/(2.*(Y2-X2)*(1.-X1)*(1.-X1));
300  Y2X2 = 1. +X1*(1-2.*X2)/(2.*(Y2-X2)*(1.-X1));
301  Y2X1X2= Y2X1*( (1.-2.*X2)/(X2*(1-X2)) -( Y2X2-1.)/(Y2-X2));
302  Y1X1X2= Y1X2*( (1.-2.*X1)/(X1*(1-X1)) -( Y1X1-1.)/(Y1-X1));
303 
304  KX1=-DELM*(Y1X1*ALX/(1.-Y1) + Y2X1*ALX/(1.-Y2));
305  KX2=-DELM*(Y2X2*ALX/(1.-Y2) + Y1X2*ALX/(1.-Y1));
306  ZX1=Z[I]*(SX1/(2.*SMIN)-KX1/DELM);
307  ZX2=Z[I]*(SX2/(2.*SMIN)-KX2/DELM);
308 
309  H1=ZX1*ZX2;
310 
311  HZ[I]=KF1/pow(Z[I],KF2);
312  INVSIG[I]=(HZ[I]*H1*16.)/S;
313 
314  }
315 }
float Float_t
Definition: RtypesCore.h:53
Definition: Rtypes.h:61
#define DEGRAD
def zdemo
Definition: zdemo.py:30
TH1F * DrawFrame(Double_t xmin, Double_t ymin, Double_t xmax, Double_t ymax, const char *title="")
Draw an empty pad frame with X and Y axis.
Definition: TPad.cxx:1475
int Int_t
Definition: RtypesCore.h:41
virtual void SetTitle(const char *title="")
Set graph title.
Definition: TGraph.cxx:2153
virtual void Draw(Option_t *chopt="")
Draw this graph with its current attributes.
Definition: TGraph.cxx:740
double cos(double)
virtual void SetTextFont(Font_t tfont=62)
Definition: TAttText.h:59
double sqrt(double)
TVirtualPad * cd(Int_t subpadnumber=0)
Set Current pad.
Definition: TPad.cxx:514
TFrame * GetFrame()
Get frame.
Definition: TPad.cxx:2729
def hz_calc
Definition: zdemo.py:235
To draw Mathematical Formula.
Definition: TLatex.h:33
#define NLOOP
Definition: piRandom.C:21
double pow(double, double)
virtual void SetMarkerColor(Color_t mcolor=1)
Definition: TAttMarker.h:51
TLatex * DrawLatex(Double_t x, Double_t y, const char *text)
Make a copy of this object with the new parameters And copy object attributes.
Definition: TLatex.cxx:1901
virtual void SetNDC(Bool_t isNDC=kTRUE)
Set NDC mode on if isNDC = kTRUE, off otherwise.
Definition: TText.cxx:809
virtual void SetLogx(Int_t value=1)
Set Lin/Log scale for X.
Definition: TPad.cxx:5301
double sin(double)
virtual void Draw(Option_t *option="")
Draw Pad in Current pad (re-parent pad if necessary).
Definition: TPad.cxx:1192
virtual void SetTextAlign(Short_t align=11)
Definition: TAttText.h:55
A Pave (see TPave) with a text centered in the Pave.
Definition: TPaveLabel.h:32
tuple HZ
Definition: zdemo.py:22
virtual void SetLineColor(Color_t lcolor)
Definition: TAttLine.h:54
virtual void SetFillColor(Color_t fcolor)
Definition: TAttFill.h:50
virtual void Range(Double_t x1, Double_t y1, Double_t x2, Double_t y2)
Set world coordinate system for the pad.
Definition: TPad.cxx:4623
tuple INVSIG
Definition: zdemo.py:24
The most important graphics class in the ROOT system.
Definition: TPad.h:46
virtual void SetMarkerStyle(Style_t mstyle=1)
Definition: TAttMarker.h:53
virtual void SetMarkerSize(Size_t msize=1)
Definition: TAttMarker.h:54
The Canvas class.
Definition: TCanvas.h:48
virtual void Draw(Option_t *option="")
Draw this pavelabel with its current attributes.
Definition: TPaveLabel.cxx:77
static const float S
Definition: mandel.cpp:113
virtual void SetPoint(Int_t i, Double_t x, Double_t y)
Set x and y values for point number i.
Definition: TGraph.cxx:2127
tuple PT
Definition: zdemo.py:23
A Graph is a graphics object made of two arrays X and Y with npoints each.
Definition: TGraph.h:53
virtual void SetTextColor(Color_t tcolor=1)
Definition: TAttText.h:57
Definition: Rtypes.h:61
virtual void SetTextSize(Float_t tsize=1)
Definition: TAttText.h:60
virtual void Update()
Update canvas pad buffers.
Definition: TCanvas.cxx:2179
#define I(x, y, z)
std::complex< float_v > Z
Definition: main.cpp:120
void Modified(Bool_t flag=1)
Definition: TPad.h:407
virtual void SetLogy(Int_t value=1)
Set Lin/Log scale for Y.
Definition: TPad.cxx:5314