ROOT logo

From $ROOTSYS/tutorials/graphs/zdemo.C

// This macro is an example of graphs in log scales with annotations.
//
//  The  begin_html <a href="gif/zdemo.gif" >presented results</a> end_html
//  are predictions of invariant cross-section of Direct Photons produced
//  at RHIC energies, based on the universality of scaling function H(z).
//
//Authors: Michael Tokarev and Elena Potrebenikova (JINR Dubna)
//
//  These Figures were published in JINR preprint E2-98-64, Dubna,
//  1998 and submitted to CPC.
//
// Note that the way greek symbols, super/subscripts are obtained
// illustrate the current limitations of Root in this area.
//

#include "TCanvas.h"
#include "TPad.h"
#include "TPaveLabel.h"
#include "TLatex.h"
#include "TGraph.h"
#include "TFrame.h"

const Int_t NMAX = 20;
Int_t NLOOP;
Float_t Z[NMAX], HZ[NMAX], PT[NMAX], INVSIG[NMAX];

void hz_calc(Float_t, Float_t, Float_t, Float_t, Float_t, Float_t);

//__________________________________________________________________
void zdemo()
{

   Float_t energ;
   Float_t dens;
   Float_t tgrad;
   Float_t ptmin;
   Float_t ptmax;
   Float_t delp;

   // Create a new canvas.
   TCanvas *c1 = new TCanvas("zdemo",
      "Monte Carlo Study of Z scaling",10,40,800,600);
   c1->Range(0,0,25,18);
   c1->SetFillColor(40);

   TPaveLabel *pl = new TPaveLabel(1,16.3,24,17.5,"Z-scaling of \
      Direct Photon Productions in pp Collisions at RHIC Energies","br");
   pl->SetFillColor(18);
   pl->SetTextFont(32);
   pl->SetTextColor(49);
   pl->Draw();

   TLatex *t = new TLatex();
   t->SetTextFont(32);
   t->SetTextColor(1);
   t->SetTextSize(0.03);
   t->SetTextAlign(12);
   t->DrawLatex(3.1,15.5,"M.Tokarev, E.Potrebenikova ");
   t->DrawLatex(14.,15.5,"JINR preprint E2-98-64, Dubna, 1998 ");

   TPad *pad1 = new TPad("pad1","This is pad1",0.02,0.02,0.48,0.83,33);
   TPad *pad2 = new TPad("pad2","This is pad2",0.52,0.02,0.98,0.83,33);

   pad1->Draw();
   pad2->Draw();

//
// Cross-section of direct photon production in pp collisions 
// at 500 GeV vs Pt
//
   energ = 63;
   dens  = 1.766;
   tgrad = 90.;
   ptmin = 4.;
   ptmax = 24.;
   delp  = 2.;
   hz_calc(energ, dens, tgrad, ptmin, ptmax, delp);
   pad1->cd();
   pad1->Range(-0.255174,-19.25,2.29657,-6.75);
   pad1->SetLogx();
   pad1->SetLogy();

   // create a 2-d histogram to define the range
   pad1->DrawFrame(1,1e-18,110,1e-8);
   pad1->GetFrame()->SetFillColor(19);
   t = new TLatex();
   t->SetNDC();
   t->SetTextFont(62);
   t->SetTextColor(36);
   t->SetTextSize(0.08);
   t->SetTextAlign(12);
   t->DrawLatex(0.6,0.85,"p - p");

   t->SetTextSize(0.05);
   t->DrawLatex(0.6,0.79,"Direct #gamma");
   t->DrawLatex(0.6,0.75,"#theta = 90^{o}");

   t->DrawLatex(0.20,0.45,"Ed^{3}#sigma/dq^{3}");
   t->DrawLatex(0.18,0.40,"(barn/Gev^{2})");

   t->SetTextSize(0.045);
   t->SetTextColor(kBlue);
   t->DrawLatex(0.22,0.260,"#sqrt{s} = 63(GeV)");
   t->SetTextColor(kRed);
   t->DrawLatex(0.22,0.205,"#sqrt{s} = 200(GeV)");
   t->SetTextColor(6);
   t->DrawLatex(0.22,0.15,"#sqrt{s} = 500(GeV)");

   t->SetTextSize(0.05);
   t->SetTextColor(1);
   t->DrawLatex(0.6,0.06,"q_{T} (Gev/c)");

   TGraph *gr1 = new TGraph(NLOOP,PT,INVSIG);

   gr1->SetLineColor(38);
   gr1->SetMarkerColor(kBlue);
   gr1->SetMarkerStyle(21);
   gr1->SetMarkerSize(1.1);
   gr1->Draw("LP");

//
// Cross-section of direct photon production in pp collisions 
// at 200 GeV vs Pt
//

   energ = 200;
   dens  = 2.25;
   tgrad = 90.;
   ptmin = 4.;
   ptmax = 64.;
   delp  = 6.;
   hz_calc(energ, dens, tgrad, ptmin, ptmax, delp);

   TGraph *gr2 = new TGraph(NLOOP,PT,INVSIG);
   gr2->SetLineColor(38);
   gr2->SetMarkerColor(kRed);
   gr2->SetMarkerStyle(29);
   gr2->SetMarkerSize(1.5);
   gr2->Draw("LP");

//
// Cross-section of direct photon production in pp collisions 
// at 500 GeV vs Pt
//
   energ = 500;
   dens  = 2.73;
   tgrad = 90.;
   ptmin = 4.;
   ptmax = 104.;
   delp  = 10.;
   hz_calc(energ, dens, tgrad, ptmin, ptmax, delp);

   TGraph *gr3 = new TGraph(NLOOP,PT,INVSIG);

   gr3->SetLineColor(38);
   gr3->SetMarkerColor(6);
   gr3->SetMarkerStyle(8);
   gr3->SetMarkerSize(1.1);
   gr3->Draw("LP");

   Float_t *dum = 0;
   TGraph *graph = new TGraph(1,dum,dum);
   graph->SetMarkerColor(kBlue);
   graph->SetMarkerStyle(21);
   graph->SetMarkerSize(1.1);
   graph->SetPoint(0,1.7,1.e-16);
   graph->Draw("LP");

   graph = new TGraph(1,dum,dum);
   graph->SetMarkerColor(kRed);
   graph->SetMarkerStyle(29);
   graph->SetMarkerSize(1.5);
   graph->SetPoint(0,1.7,2.e-17);
   graph->Draw("LP");

   graph = new TGraph(1,dum,dum);
   graph->SetMarkerColor(6);
   graph->SetMarkerStyle(8);
   graph->SetMarkerSize(1.1);
   graph->SetPoint(0,1.7,4.e-18);
   graph->Draw("LP");

   pad2->cd();
   pad2->Range(-0.43642,-23.75,3.92778,-6.25);
   pad2->SetLogx();
   pad2->SetLogy();

   pad2->DrawFrame(1,1e-22,3100,1e-8);
   pad2->GetFrame()->SetFillColor(19);

   TGraph *gr = new TGraph(NLOOP,Z,HZ);
   gr->SetTitle("HZ vs Z");
   gr->SetFillColor(19);
   gr->SetLineColor(9);
   gr->SetMarkerColor(50);
   gr->SetMarkerStyle(29);
   gr->SetMarkerSize(1.5);
   gr->Draw("LP");

   t = new TLatex();
   t->SetNDC();
   t->SetTextFont(62);
   t->SetTextColor(36);
   t->SetTextSize(0.08);
   t->SetTextAlign(12);
   t->DrawLatex(0.6,0.85,"p - p");

   t->SetTextSize(0.05);
   t->DrawLatex(0.6,0.79,"Direct #gamma");
   t->DrawLatex(0.6,0.75,"#theta = 90^{o}");

   t->DrawLatex(0.70,0.55,"H(z)");
   t->DrawLatex(0.68,0.50,"(barn)");

   t->SetTextSize(0.045);
   t->SetTextColor(46);
   t->DrawLatex(0.20,0.30,"#sqrt{s}, GeV");
   t->DrawLatex(0.22,0.26,"63");
   t->DrawLatex(0.22,0.22,"200");
   t->DrawLatex(0.22,0.18,"500");

   t->SetTextSize(0.05);
   t->SetTextColor(1);
   t->DrawLatex(0.88,0.06,"z");

   c1->Modified();
   c1->Update();
}

void hz_calc(Float_t ENERG, Float_t DENS, Float_t TGRAD, Float_t PTMIN, 
   Float_t PTMAX, Float_t DELP)
{
  Int_t I;

  Float_t GM1  = 0.00001;
  Float_t GM2  = 0.00001;
  Float_t A1   = 1.;
  Float_t A2   = 1.;
  Float_t ALX  = 2.;
  Float_t BETA = 1.;
  Float_t KF1  = 8.E-7;
  Float_t KF2  = 5.215;

  Float_t MN = 0.9383;
  Float_t DEGRAD=0.01745329;

  Float_t EB1, EB2, PB1, PB2, MB1, MB2, M1, M2;
  Float_t DNDETA;

  Float_t P1P2, P1P3, P2P3;
  Float_t Y1, Y2, S, SMIN,  SX1,  SX2, SX1X2, DELM;
  Float_t Y1X1,  Y1X2,   Y2X1,   Y2X2,   Y2X1X2,   Y1X1X2;
  Float_t KX1, KX2,  ZX1, ZX2;
  Float_t H1;

  Float_t PTOT, THET, ETOT, X1, X2;

  DNDETA= DENS;
  MB1   = MN*A1;
  MB2   = MN*A2;
  EB1   = ENERG/2.*A1;
  EB2   = ENERG/2.*A2;
  M1    = GM1;
  M2    = GM2;
  THET  = TGRAD*DEGRAD;
  NLOOP = (PTMAX-PTMIN)/DELP;

  for (I=0; I<NLOOP;I++) {
     PT[I]=PTMIN+I*DELP;
     PTOT = PT[I]/sin(THET);

     ETOT = sqrt(M1*M1 + PTOT*PTOT);
     PB1  = sqrt(EB1*EB1 - MB1*MB1);
     PB2  = sqrt(EB2*EB2 - MB2*MB2);
     P2P3 = EB2*ETOT+PB2*PTOT*cos(THET);
     P1P2 = EB2*EB1+PB2*PB1;
     P1P3 = EB1*ETOT-PB1*PTOT*cos(THET);

     X1 = P2P3/P1P2;
     X2 = P1P3/P1P2;
     Y1 = X1+sqrt(X1*X2*(1.-X1)/(1.-X2));
     Y2 = X2+sqrt(X1*X2*(1.-X2)/(1.-X1));

     S    = (MB1*MB1)+2.*P1P2+(MB2*MB2);
     SMIN = 4.*((MB1*MB1)*(X1*X1) +2.*X1*X2*P1P2+(MB2*MB2)*(X2*X2));
     SX1  = 4.*( 2*(MB1*MB1)*X1+2*X2*P1P2);
     SX2  = 4.*( 2*(MB2*MB2)*X2+2*X1*P1P2);
     SX1X2= 4.*(2*P1P2);
     DELM = pow((1.-Y1)*(1.-Y2),ALX);

     Z[I] = sqrt(SMIN)/DELM/pow(DNDETA,BETA);

     Y1X1  = 1. +X2*(1-2.*X1)/(2.*(Y1-X1)*(1.-X2));
     Y1X2  =     X1*(1-X1)/(2.*(Y1-X1)*(1.-X2)*(1.-X2));
     Y2X1  =     X2*(1-X2)/(2.*(Y2-X2)*(1.-X1)*(1.-X1));
     Y2X2  = 1. +X1*(1-2.*X2)/(2.*(Y2-X2)*(1.-X1));
     Y2X1X2= Y2X1*( (1.-2.*X2)/(X2*(1-X2)) -( Y2X2-1.)/(Y2-X2));
     Y1X1X2= Y1X2*( (1.-2.*X1)/(X1*(1-X1)) -( Y1X1-1.)/(Y1-X1));

     KX1=-DELM*(Y1X1*ALX/(1.-Y1) + Y2X1*ALX/(1.-Y2));
     KX2=-DELM*(Y2X2*ALX/(1.-Y2) + Y1X2*ALX/(1.-Y1));
     ZX1=Z[I]*(SX1/(2.*SMIN)-KX1/DELM);
     ZX2=Z[I]*(SX2/(2.*SMIN)-KX2/DELM);

     H1=ZX1*ZX2;

     HZ[I]=KF1/pow(Z[I],KF2);
     INVSIG[I]=(HZ[I]*H1*16.)/S;

  }
}
 zdemo.C:1
 zdemo.C:2
 zdemo.C:3
 zdemo.C:4
 zdemo.C:5
 zdemo.C:6
 zdemo.C:7
 zdemo.C:8
 zdemo.C:9
 zdemo.C:10
 zdemo.C:11
 zdemo.C:12
 zdemo.C:13
 zdemo.C:14
 zdemo.C:15
 zdemo.C:16
 zdemo.C:17
 zdemo.C:18
 zdemo.C:19
 zdemo.C:20
 zdemo.C:21
 zdemo.C:22
 zdemo.C:23
 zdemo.C:24
 zdemo.C:25
 zdemo.C:26
 zdemo.C:27
 zdemo.C:28
 zdemo.C:29
 zdemo.C:30
 zdemo.C:31
 zdemo.C:32
 zdemo.C:33
 zdemo.C:34
 zdemo.C:35
 zdemo.C:36
 zdemo.C:37
 zdemo.C:38
 zdemo.C:39
 zdemo.C:40
 zdemo.C:41
 zdemo.C:42
 zdemo.C:43
 zdemo.C:44
 zdemo.C:45
 zdemo.C:46
 zdemo.C:47
 zdemo.C:48
 zdemo.C:49
 zdemo.C:50
 zdemo.C:51
 zdemo.C:52
 zdemo.C:53
 zdemo.C:54
 zdemo.C:55
 zdemo.C:56
 zdemo.C:57
 zdemo.C:58
 zdemo.C:59
 zdemo.C:60
 zdemo.C:61
 zdemo.C:62
 zdemo.C:63
 zdemo.C:64
 zdemo.C:65
 zdemo.C:66
 zdemo.C:67
 zdemo.C:68
 zdemo.C:69
 zdemo.C:70
 zdemo.C:71
 zdemo.C:72
 zdemo.C:73
 zdemo.C:74
 zdemo.C:75
 zdemo.C:76
 zdemo.C:77
 zdemo.C:78
 zdemo.C:79
 zdemo.C:80
 zdemo.C:81
 zdemo.C:82
 zdemo.C:83
 zdemo.C:84
 zdemo.C:85
 zdemo.C:86
 zdemo.C:87
 zdemo.C:88
 zdemo.C:89
 zdemo.C:90
 zdemo.C:91
 zdemo.C:92
 zdemo.C:93
 zdemo.C:94
 zdemo.C:95
 zdemo.C:96
 zdemo.C:97
 zdemo.C:98
 zdemo.C:99
 zdemo.C:100
 zdemo.C:101
 zdemo.C:102
 zdemo.C:103
 zdemo.C:104
 zdemo.C:105
 zdemo.C:106
 zdemo.C:107
 zdemo.C:108
 zdemo.C:109
 zdemo.C:110
 zdemo.C:111
 zdemo.C:112
 zdemo.C:113
 zdemo.C:114
 zdemo.C:115
 zdemo.C:116
 zdemo.C:117
 zdemo.C:118
 zdemo.C:119
 zdemo.C:120
 zdemo.C:121
 zdemo.C:122
 zdemo.C:123
 zdemo.C:124
 zdemo.C:125
 zdemo.C:126
 zdemo.C:127
 zdemo.C:128
 zdemo.C:129
 zdemo.C:130
 zdemo.C:131
 zdemo.C:132
 zdemo.C:133
 zdemo.C:134
 zdemo.C:135
 zdemo.C:136
 zdemo.C:137
 zdemo.C:138
 zdemo.C:139
 zdemo.C:140
 zdemo.C:141
 zdemo.C:142
 zdemo.C:143
 zdemo.C:144
 zdemo.C:145
 zdemo.C:146
 zdemo.C:147
 zdemo.C:148
 zdemo.C:149
 zdemo.C:150
 zdemo.C:151
 zdemo.C:152
 zdemo.C:153
 zdemo.C:154
 zdemo.C:155
 zdemo.C:156
 zdemo.C:157
 zdemo.C:158
 zdemo.C:159
 zdemo.C:160
 zdemo.C:161
 zdemo.C:162
 zdemo.C:163
 zdemo.C:164
 zdemo.C:165
 zdemo.C:166
 zdemo.C:167
 zdemo.C:168
 zdemo.C:169
 zdemo.C:170
 zdemo.C:171
 zdemo.C:172
 zdemo.C:173
 zdemo.C:174
 zdemo.C:175
 zdemo.C:176
 zdemo.C:177
 zdemo.C:178
 zdemo.C:179
 zdemo.C:180
 zdemo.C:181
 zdemo.C:182
 zdemo.C:183
 zdemo.C:184
 zdemo.C:185
 zdemo.C:186
 zdemo.C:187
 zdemo.C:188
 zdemo.C:189
 zdemo.C:190
 zdemo.C:191
 zdemo.C:192
 zdemo.C:193
 zdemo.C:194
 zdemo.C:195
 zdemo.C:196
 zdemo.C:197
 zdemo.C:198
 zdemo.C:199
 zdemo.C:200
 zdemo.C:201
 zdemo.C:202
 zdemo.C:203
 zdemo.C:204
 zdemo.C:205
 zdemo.C:206
 zdemo.C:207
 zdemo.C:208
 zdemo.C:209
 zdemo.C:210
 zdemo.C:211
 zdemo.C:212
 zdemo.C:213
 zdemo.C:214
 zdemo.C:215
 zdemo.C:216
 zdemo.C:217
 zdemo.C:218
 zdemo.C:219
 zdemo.C:220
 zdemo.C:221
 zdemo.C:222
 zdemo.C:223
 zdemo.C:224
 zdemo.C:225
 zdemo.C:226
 zdemo.C:227
 zdemo.C:228
 zdemo.C:229
 zdemo.C:230
 zdemo.C:231
 zdemo.C:232
 zdemo.C:233
 zdemo.C:234
 zdemo.C:235
 zdemo.C:236
 zdemo.C:237
 zdemo.C:238
 zdemo.C:239
 zdemo.C:240
 zdemo.C:241
 zdemo.C:242
 zdemo.C:243
 zdemo.C:244
 zdemo.C:245
 zdemo.C:246
 zdemo.C:247
 zdemo.C:248
 zdemo.C:249
 zdemo.C:250
 zdemo.C:251
 zdemo.C:252
 zdemo.C:253
 zdemo.C:254
 zdemo.C:255
 zdemo.C:256
 zdemo.C:257
 zdemo.C:258
 zdemo.C:259
 zdemo.C:260
 zdemo.C:261
 zdemo.C:262
 zdemo.C:263
 zdemo.C:264
 zdemo.C:265
 zdemo.C:266
 zdemo.C:267
 zdemo.C:268
 zdemo.C:269
 zdemo.C:270
 zdemo.C:271
 zdemo.C:272
 zdemo.C:273
 zdemo.C:274
 zdemo.C:275
 zdemo.C:276
 zdemo.C:277
 zdemo.C:278
 zdemo.C:279
 zdemo.C:280
 zdemo.C:281
 zdemo.C:282
 zdemo.C:283
 zdemo.C:284
 zdemo.C:285
 zdemo.C:286
 zdemo.C:287
 zdemo.C:288
 zdemo.C:289
 zdemo.C:290
 zdemo.C:291
 zdemo.C:292
 zdemo.C:293
 zdemo.C:294
 zdemo.C:295
 zdemo.C:296
 zdemo.C:297
 zdemo.C:298
 zdemo.C:299
 zdemo.C:300
 zdemo.C:301
 zdemo.C:302
 zdemo.C:303
 zdemo.C:304
 zdemo.C:305
 zdemo.C:306
 zdemo.C:307
 zdemo.C:308
 zdemo.C:309
 zdemo.C:310
 zdemo.C:311
 zdemo.C:312