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