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