1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
18 // Utility class to make simple Glauber type calculations
19 // for SYMMETRIC collision geometries (AA):
20 // Impact parameter, production points, reaction plane dependence
22 // The SimulateTrigger method can be used for simple MB and hard-process
23 // (binary scaling) trigger studies.
25 // Some basic quantities can be visualized directly.
27 // The default set-up for PbPb or AUAu collisions can be read from a file
28 // calling Init(1) or Init(2) if you want to read Almonds too.
30 // ***** If you change settings dont forget to call init afterwards, *****
31 // ***** in order to update the formulas with the new parameters. *****
33 // Author: andreas.morsch@cern.ch
34 //=================== Added by A. Dainese 11/02/04 ===========================
35 // Calculate path length for a parton with production point (x0,y0)
36 // and propagation direction (ux=cos(phi0),uy=sin(phi0))
37 // in a collision with impact parameter b and functions that make use
39 //=================== Added by A. Dainese 05/03/04 ===========================
40 // Calculation of line integrals I0 and I1
41 // integral0 = \int_0^ellCut dl*(T_A*T_B)(x0+l*ux,y0+l*uy)
42 // integral1 = \int_0^ellCut dl*l*(T_A*T_B)(x0+l*ux,y0+l*uy)
43 // mostly for use in the Quenching class
44 //=================== Added by C. Loizdes 27/03/04 ===========================
45 // Handling of AuAu collisions
46 // More get/set functions
47 // Comments, units and clearing of code
51 #include "AliFastGlauber.h"
53 #include <Riostream.h>
66 ClassImp(AliFastGlauber)
68 Float_t AliFastGlauber::fgBMax = 0.;
69 TF1* AliFastGlauber::fgWSb = NULL;
70 TF2* AliFastGlauber::fgWSbz = NULL;
71 TF1* AliFastGlauber::fgWSz = NULL;
72 TF1* AliFastGlauber::fgWSta = NULL;
73 TF2* AliFastGlauber::fgWStarfi = NULL;
74 TF2* AliFastGlauber::fgWAlmond = NULL;
75 TF1* AliFastGlauber::fgWStaa = NULL;
76 TF1* AliFastGlauber::fgWSgeo = NULL;
77 TF1* AliFastGlauber::fgWSbinary = NULL;
78 TF1* AliFastGlauber::fgWSN = NULL;
79 TF1* AliFastGlauber::fgWPathLength0 = NULL;
80 TF1* AliFastGlauber::fgWPathLength = NULL;
81 TF1* AliFastGlauber::fgWEnergyDensity = NULL;
82 TF1* AliFastGlauber::fgWIntRadius = NULL;
83 TF2* AliFastGlauber::fgWKParticipants = NULL;
84 TF1* AliFastGlauber::fgWParticipants = NULL;
85 TF2* AliFastGlauber::fgWAlmondCurrent = NULL;
86 TF2* AliFastGlauber::fgWAlmondFixedB[40];
87 const Int_t AliFastGlauber::fgkMCInts = 100000;
88 AliFastGlauber* AliFastGlauber::fgGlauber = NULL;
91 AliFastGlauber::AliFastGlauber():
104 // Default Constructor
107 SetLengthDefinition();
110 fI0I1[0] = fI0I1[1] = 0;
113 AliFastGlauber::AliFastGlauber(const AliFastGlauber & gl)
130 fI0I1[0] = fI0I1[1] = 0;
133 AliFastGlauber* AliFastGlauber::Instance()
135 // Set random number generator
139 fgGlauber = new AliFastGlauber();
144 AliFastGlauber::~AliFastGlauber()
147 for(Int_t k=0; k<40; k++) delete fgWAlmondFixedB[k];
150 void AliFastGlauber::SetAuAuRhic()
152 //Set all parameters for RHIC
153 SetWoodSaxonParametersAu();
154 SetHardCrossSection();
155 SetNNCrossSection(42);
157 SetFileName("$(ALICE_ROOT)/FASTSIM/data/glauberAuAu.root");
160 void AliFastGlauber::SetPbPbLHC()
162 //Set all parameters for LHC
163 SetWoodSaxonParametersPb();
164 SetHardCrossSection();
170 void AliFastGlauber::Init(Int_t mode)
173 // mode = 0; all functions are calculated
174 // mode = 1; overlap function is read from file (for Pb-Pb only)
175 // mode = 2; interaction almond functions are read from file
176 // USE THIS FOR PATH LENGTH CALC.!
186 fgWSb = new TF1("WSb", WSb, 0, fgBMax, 4);
187 fgWSb->SetParameter(0, fWSr0);
188 fgWSb->SetParameter(1, fWSd);
189 fgWSb->SetParameter(2, fWSw);
190 fgWSb->SetParameter(3, fWSn);
192 fgWSbz = new TF2("WSbz", WSbz, 0, fgBMax, 0, fgBMax, 4);
193 fgWSbz->SetParameter(0, fWSr0);
194 fgWSbz->SetParameter(1, fWSd);
195 fgWSbz->SetParameter(2, fWSw);
196 fgWSbz->SetParameter(3, fWSn);
198 fgWSz = new TF1("WSz", WSz, 0, fgBMax, 5);
199 fgWSz->SetParameter(0, fWSr0);
200 fgWSz->SetParameter(1, fWSd);
201 fgWSz->SetParameter(2, fWSw);
202 fgWSz->SetParameter(3, fWSn);
207 fgWSta = new TF1("WSta", WSta, 0., fgBMax, 0);
212 fgWStarfi = new TF2("WStarfi", WStarfi, 0., fgBMax, 0., TMath::Pi(), 1);
213 fgWStarfi->SetParameter(0, 0.);
214 fgWStarfi->SetNpx(200);
215 fgWStarfi->SetNpy(20);
218 // Participants Kernel
220 fgWKParticipants = new TF2("WKParticipants", WKParticipants, 0., fgBMax, 0., TMath::Pi(), 3);
221 fgWKParticipants->SetParameter(0, 0.);
222 fgWKParticipants->SetParameter(1, fSigmaNN);
223 fgWKParticipants->SetParameter(2, fA);
224 fgWKParticipants->SetNpx(200);
225 fgWKParticipants->SetNpy(20);
228 // Overlap and Participants
231 fgWStaa = new TF1("WStaa", WStaa, 0., fgBMax, 1);
232 fgWStaa->SetNpx(100);
233 fgWStaa->SetParameter(0,fA);
234 fgWStaa->SetNpx(100);
235 fgWParticipants = new TF1("WParticipants", WParticipants, 0., fgBMax, 2);
236 fgWParticipants->SetParameter(0, fSigmaNN);
237 fgWParticipants->SetParameter(1, fA);
238 fgWParticipants->SetNpx(100);
240 Info("Init","Reading overlap function from file %s",fName.Data());
241 TFile* f = new TFile(fName.Data());
243 Fatal("Init", "Could not open file %s",fName.Data());
245 fgWStaa = (TF1*) f->Get("WStaa");
246 fgWParticipants = (TF1*) f->Get("WParticipants");
253 fgWEnergyDensity = new TF1("WEnergyDensity", WEnergyDensity, 0., 2. * fWSr0, 1);
254 fgWEnergyDensity->SetParameter(0, fWSr0 + 1.);
257 // Geometrical Cross-Section
259 fgWSgeo = new TF1("WSgeo", WSgeo, 0., fgBMax, 1);
260 fgWSgeo->SetParameter(0,fSigmaNN); //mbarn
261 fgWSgeo->SetNpx(100);
264 // Hard cross section (binary collisions)
266 fgWSbinary = new TF1("WSbinary", WSbinary, 0., fgBMax, 1);
267 fgWSbinary->SetParameter(0, fSigmaHard); //mbarn
268 fgWSbinary->SetNpx(100);
271 // Hard collisions per event
273 fgWSN = new TF1("WSN", WSN, 0., fgBMax, 1);
277 // Almond shaped interaction region
279 fgWAlmond = new TF2("WAlmond", WAlmond, -fgBMax, fgBMax, -fgBMax, fgBMax, 1);
280 fgWAlmond->SetParameter(0, 0.);
281 fgWAlmond->SetNpx(200);
282 fgWAlmond->SetNpy(200);
285 Info("Init","Reading interaction almonds from file: %s",fName.Data());
286 Char_t almondName[100];
287 TFile* ff = new TFile(fName.Data());
288 for(Int_t k=0; k<40; k++) {
289 sprintf(almondName,"WAlmondFixedB%d",k);
290 fgWAlmondCurrent = (TF2*)ff->Get(almondName);
291 fgWAlmondFixedB[k] = fgWAlmondCurrent;
296 fgWIntRadius = new TF1("WIntRadius", WIntRadius, 0., fgBMax, 1);
297 fgWIntRadius->SetParameter(0, 0.);
300 // Path Length as a function of Phi
302 fgWPathLength0 = new TF1("WPathLength0", WPathLength0, -TMath::Pi(), TMath::Pi(), 2);
303 fgWPathLength0->SetParameter(0, 0.);
304 fgWPathLength0->SetParameter(1, 0.); //Pathlength definition
306 fgWPathLength = new TF1("WPathLength", WPathLength, -TMath::Pi(), TMath::Pi(), 3);
307 fgWPathLength->SetParameter(0, 0.); //Impact Parameter
308 fgWPathLength->SetParameter(1, 1000.); //Number of interactions used for average
309 fgWPathLength->SetParameter(2, 0); //Pathlength definition
312 void AliFastGlauber::Reset() const
315 // Reset dynamic allocated formulas
316 // in case init is called twice
318 if(fgWSb) delete fgWSb;
319 if(fgWSbz) delete fgWSbz;
320 if(fgWSz) delete fgWSz;
321 if(fgWSta) delete fgWSta;
322 if(fgWStarfi) delete fgWStarfi;
323 if(fgWAlmond) delete fgWAlmond;
324 if(fgWStaa) delete fgWStaa;
325 if(fgWSgeo) delete fgWSgeo;
326 if(fgWSbinary) delete fgWSbinary;
327 if(fgWSN) delete fgWSN;
328 if(fgWPathLength0) delete fgWPathLength0;
329 if(fgWPathLength) delete fgWPathLength;
330 if(fgWEnergyDensity) delete fgWEnergyDensity;
331 if(fgWIntRadius) delete fgWIntRadius;
332 if(fgWKParticipants) delete fgWKParticipants;
333 if(fgWParticipants) delete fgWParticipants;
336 void AliFastGlauber::DrawWSb() const
339 // Draw Wood-Saxon Nuclear Density Function
341 TCanvas *c1 = new TCanvas("c1","Wood Saxon",400,10,600,700);
343 Double_t max=fgWSb->GetMaximum(0,fgBMax)*1.01;
344 TH2F *h2f=new TH2F("h2fwsb","Wood Saxon: #rho(r) = n (1-#omega(r/r_{0})^2)/(1+exp((r-r_{0})/d)) [fm^{-3}]",2,0,fgBMax,2,0,max);
346 h2f->GetXaxis()->SetTitle("r [fm]");
347 h2f->GetYaxis()->SetNoExponent(kTRUE);
348 h2f->GetYaxis()->SetTitle("#rho [fm^{-3}]");
351 TLegend *l1a = new TLegend(0.45,0.6,.90,0.8);
352 l1a->SetFillStyle(0);
353 l1a->SetBorderSize(0);
355 sprintf(label,"r_{0} = %.2f fm",fWSr0);
356 l1a->AddEntry(fgWSb,label,"");
357 sprintf(label,"d = %.2f fm",fWSd);
358 l1a->AddEntry(fgWSb,label,"");
359 sprintf(label,"n = %.2e fm^{-3}",fWSn);
360 l1a->AddEntry(fgWSb,label,"");
361 sprintf(label,"#omega = %.2f",fWSw);
362 l1a->AddEntry(fgWSb,label,"");
367 void AliFastGlauber::DrawOverlap() const
370 // Draw Overlap Function
372 TCanvas *c2 = new TCanvas("c2","Overlap",400,10,600,700);
374 Double_t max=fgWStaa->GetMaximum(0,fgBMax)*1.01;
375 TH2F *h2f=new TH2F("h2ftaa","Overlap function: T_{AB} [mbarn^{-1}]",2,0,fgBMax,2,0, max);
377 h2f->GetXaxis()->SetTitle("b [fm]");
378 h2f->GetYaxis()->SetTitle("T_{AB} [mbarn^{-1}]");
380 fgWStaa->Draw("same");
383 void AliFastGlauber::DrawParticipants() const
386 // Draw Number of Participants Npart
388 TCanvas *c3 = new TCanvas("c3","Participants",400,10,600,700);
390 Double_t max=fgWParticipants->GetMaximum(0,fgBMax)*1.01;
391 TH2F *h2f=new TH2F("h2fpart","Number of Participants",2,0,fgBMax,2,0,max);
393 h2f->GetXaxis()->SetTitle("b [fm]");
394 h2f->GetYaxis()->SetTitle("N_{part}");
396 fgWParticipants->Draw("same");
397 TLegend *l1a = new TLegend(0.50,0.75,.90,0.9);
398 l1a->SetFillStyle(0);
399 l1a->SetBorderSize(0);
401 sprintf(label,"#sigma^{inel.}_{NN} = %.1f mbarn",fSigmaNN);
402 l1a->AddEntry(fgWParticipants,label,"");
407 void AliFastGlauber::DrawThickness() const
410 // Draw Thickness Function
412 TCanvas *c4 = new TCanvas("c4","Thickness",400,10,600,700);
414 Double_t max=fgWSta->GetMaximum(0,fgBMax)*1.01;
415 TH2F *h2f=new TH2F("h2fta","Thickness function: T_{A} [fm^{-2}]",2,0,fgBMax,2,0,max);
417 h2f->GetXaxis()->SetTitle("b [fm]");
418 h2f->GetYaxis()->SetTitle("T_{A} [fm^{-2}]");
420 fgWSta->Draw("same");
423 void AliFastGlauber::DrawGeo() const
426 // Draw Geometrical Cross-Section
428 TCanvas *c5 = new TCanvas("c5","Geometrical Cross-Section",400,10,600,700);
430 Double_t max=fgWSgeo->GetMaximum(0,fgBMax)*1.01;
431 TH2F *h2f=new TH2F("h2fgeo","Differential Geometrical Cross-Section: d#sigma^{geo}_{AB}/db [fm]",2,0,fgBMax,2,0,max);
433 h2f->GetXaxis()->SetTitle("b [fm]");
434 h2f->GetYaxis()->SetTitle("d#sigma^{geo}_{AB}/db [fm]");
436 fgWSgeo->Draw("same");
437 TLegend *l1a = new TLegend(0.10,0.8,.40,0.9);
438 l1a->SetFillStyle(0);
439 l1a->SetBorderSize(0);
441 sprintf(label,"#sigma_{NN}^{inel.} = %.1f mbarn",fSigmaNN);
442 l1a->AddEntry(fgWSgeo,label,"");
447 void AliFastGlauber::DrawBinary() const
450 // Draw Binary Cross-Section
452 TCanvas *c6 = new TCanvas("c6","Binary Cross-Section",400,10,600,700);
454 Double_t max=fgWSbinary->GetMaximum(0,fgBMax)*1.01;
455 TH2F *h2f=new TH2F("h2fbinary","Differential Binary Cross-Section: #sigma^{hard}_{NN} dT_{AB}/db [fm]",2,0,fgBMax,2,0,max);
457 h2f->GetXaxis()->SetTitle("b [fm]");
458 h2f->GetYaxis()->SetTitle("d#sigma^{hard}_{AB}/db [fm]");
460 fgWSbinary->Draw("same");
461 TLegend *l1a = new TLegend(0.50,0.8,.90,0.9);
462 l1a->SetFillStyle(0);
463 l1a->SetBorderSize(0);
465 sprintf(label,"#sigma_{NN}^{hard} = %.1f mbarn",fSigmaHard);
466 l1a->AddEntry(fgWSb,label,"");
471 void AliFastGlauber::DrawN() const
474 // Draw Binaries per event (Ncoll)
476 TCanvas *c7 = new TCanvas("c7","Binaries per event",400,10,600,700);
478 Double_t max=fgWSN->GetMaximum(0,fgBMax)*1.01;
479 TH2F *h2f=new TH2F("h2fhardcols","Number of hard collisions: T_{AB} #sigma^{hard}_{NN}/#sigma_{AB}^{geo}",2,0,fgBMax,2,0,max);
481 h2f->GetXaxis()->SetTitle("b [fm]");
482 h2f->GetYaxis()->SetTitle("N_{coll}");
485 TLegend *l1a = new TLegend(0.50,0.75,.90,0.9);
486 l1a->SetFillStyle(0);
487 l1a->SetBorderSize(0);
489 sprintf(label,"#sigma^{hard}_{NN} = %.1f mbarn",fSigmaHard);
490 l1a->AddEntry(fgWSN,label,"");
491 sprintf(label,"#sigma^{inel.}_{NN} = %.1f mbarn",fSigmaNN);
492 l1a->AddEntry(fgWSN,label,"");
497 void AliFastGlauber::DrawKernel(Double_t b) const
502 TCanvas *c8 = new TCanvas("c8","Kernel",400,10,600,700);
504 fgWStarfi->SetParameter(0, b);
505 TH2F *h2f=new TH2F("h2fkernel","Kernel of Overlap function: d^{2}T_{AB}/dr/d#phi [fm^{-3}]",2,0,fgBMax,2,0,TMath::Pi());
507 h2f->GetXaxis()->SetTitle("r [fm]");
508 h2f->GetYaxis()->SetTitle("#phi [rad]");
510 fgWStarfi->Draw("same");
511 TLegend *l1a = new TLegend(0.65,0.8,.90,0.9);
512 l1a->SetFillStyle(0);
513 l1a->SetBorderSize(0);
515 sprintf(label,"b = %.1f fm",b);
516 l1a->AddEntry(fgWStarfi,label,"");
521 void AliFastGlauber::DrawAlmond(Double_t b) const
524 // Draw Interaction Almond
526 TCanvas *c9 = new TCanvas("c9","Almond",400,10,600,700);
528 fgWAlmond->SetParameter(0, b);
529 TH2F *h2f=new TH2F("h2falmond","Interaction Almond [fm^{-4}]",2,0,fgBMax,2,0,fgBMax);
531 h2f->GetXaxis()->SetTitle("x [fm]");
532 h2f->GetYaxis()->SetTitle("y [fm]");
534 fgWAlmond->Draw("same");
535 TLegend *l1a = new TLegend(0.65,0.8,.90,0.9);
536 l1a->SetFillStyle(0);
537 l1a->SetBorderSize(0);
539 sprintf(label,"b = %.1f fm",b);
540 l1a->AddEntry(fgWAlmond,label,"");
545 void AliFastGlauber::DrawEnergyDensity() const
548 // Draw energy density
550 TCanvas *c10 = new TCanvas("c10","Energy Density",400, 10, 600, 700);
552 fgWEnergyDensity->SetMinimum(0.);
553 Double_t max=fgWEnergyDensity->GetMaximum(0,fgWEnergyDensity->GetParameter(0))*1.01;
554 TH2F *h2f=new TH2F("h2fenergydens","Energy density",2,0,fgBMax,2,0,max);
556 h2f->GetXaxis()->SetTitle("b [fm]");
557 h2f->GetYaxis()->SetTitle("fm^{-4}");
559 fgWEnergyDensity->Draw("same");
563 void AliFastGlauber::DrawPathLength0(Double_t b, Int_t iopt) const
568 TCanvas *c11 = new TCanvas("c11","Path Length",400,10,600,700);
570 fgWPathLength0->SetParameter(0, b);
571 fgWPathLength0->SetParameter(1, Double_t(iopt));
572 fgWPathLength0->SetMinimum(0.);
573 fgWPathLength0->SetMaximum(10.);
574 TH2F *h2f=new TH2F("h2fpathlength0","Path length",2,-TMath::Pi(), TMath::Pi(),2,0,10.);
576 h2f->GetXaxis()->SetTitle("#phi [rad]");
577 h2f->GetYaxis()->SetTitle("l [fm]");
579 fgWPathLength0->Draw("same");
582 void AliFastGlauber::DrawPathLength(Double_t b , Int_t ni, Int_t iopt) const
587 TCanvas *c12 = new TCanvas("c12","Path Length",400,10,600,700);
589 fgWAlmond->SetParameter(0, b);
590 fgWPathLength->SetParameter(0, b);
591 fgWPathLength->SetParameter(1, Double_t (ni));
592 fgWPathLength->SetParameter(2, Double_t (iopt));
593 fgWPathLength->SetMinimum(0.);
594 fgWPathLength->SetMaximum(10.);
595 TH2F *h2f=new TH2F("h2fpathlength","Path length",2,-TMath::Pi(), TMath::Pi(),2,0,10.);
597 h2f->GetXaxis()->SetTitle("#phi [rad]");
598 h2f->GetYaxis()->SetTitle("l [fm]");
600 fgWPathLength->Draw("same");
603 void AliFastGlauber::DrawIntRadius(Double_t b) const
606 // Draw Interaction Radius
608 TCanvas *c13 = new TCanvas("c13","Interaction Radius",400,10,600,700);
610 fgWIntRadius->SetParameter(0, b);
611 fgWIntRadius->SetMinimum(0);
612 Double_t max=fgWIntRadius->GetMaximum(0,fgBMax)*1.01;
613 TH2F *h2f=new TH2F("h2fintradius","Interaction Density",2,0.,fgBMax,2,0,max);
615 h2f->GetXaxis()->SetTitle("r [fm]");
616 h2f->GetYaxis()->SetTitle("[fm^{-3}]");
618 fgWIntRadius->Draw("same");
621 Double_t AliFastGlauber::WSb(Double_t* x, Double_t* par)
624 // Woods-Saxon Parameterisation
625 // as a function of radius (xx)
627 const Double_t kxx = x[0]; //fm
628 const Double_t kr0 = par[0]; //fm
629 const Double_t kd = par[1]; //fm
630 const Double_t kw = par[2]; //no units
631 const Double_t kn = par[3]; //fm^-3 (used to normalize integral to one)
632 Double_t y = kn * (1.+kw*(kxx/kr0)*(kxx/kr0))/(1.+TMath::Exp((kxx-kr0)/kd));
636 Double_t AliFastGlauber::WSbz(Double_t* x, Double_t* par)
639 // Wood Saxon Parameterisation
640 // as a function of z and b
642 const Double_t kbb = x[0]; //fm
643 const Double_t kzz = x[1]; //fm
644 const Double_t kr0 = par[0]; //fm
645 const Double_t kd = par[1]; //fm
646 const Double_t kw = par[2]; //no units
647 const Double_t kn = par[3]; //fm^-3 (used to normalize integral to one)
648 const Double_t kxx = TMath::Sqrt(kbb*kbb+kzz*kzz);
649 Double_t y = kn * (1.+kw*(kxx/kr0)*(kxx/kr0))/(1.+TMath::Exp((kxx-kr0)/kd));
653 Double_t AliFastGlauber::WSz(Double_t* x, Double_t* par)
656 // Wood Saxon Parameterisation
657 // as a function of z for fixed b
659 const Double_t kzz = x[0]; //fm
660 const Double_t kr0 = par[0]; //fm
661 const Double_t kd = par[1]; //fm
662 const Double_t kw = par[2]; //no units
663 const Double_t kn = par[3]; //fm^-3 (used to normalize integral to one)
664 const Double_t kbb = par[4]; //fm
665 const Double_t kxx = TMath::Sqrt(kbb*kbb+kzz*kzz);
666 Double_t y = kn * (1.+kw*(kxx/kr0)*(kxx/kr0))/(1.+TMath::Exp((kxx-kr0)/kd));
670 Double_t AliFastGlauber::WSta(Double_t* x, Double_t* /*par*/)
673 // Thickness function T_A
674 // as a function of b
676 const Double_t kb = x[0];
677 fgWSz->SetParameter(4, kb);
678 Double_t y = 2. * fgWSz->Integral(0., fgBMax);
682 Double_t AliFastGlauber::WStarfi(Double_t* x, Double_t* par)
685 // Kernel for overlap function: T_A(s)*T_A(s-b)
686 // as a function of r and phi
687 const Double_t kr1 = x[0];
688 const Double_t kphi = x[1];
689 const Double_t kb = par[0];
690 const Double_t kr2 = TMath::Sqrt(kr1*kr1 + kb*kb - 2.*kr1*kb*TMath::Cos(kphi));
691 Double_t y = kr1 * fgWSta->Eval(kr1) * fgWSta->Eval(kr2);
695 Double_t AliFastGlauber::WStaa(Double_t* x, Double_t* par)
699 // T_{AB}=Int d2s T_A(s)*T_B(s-b)
700 // as a function of b
701 // (normalized to fA*fB)
703 const Double_t kb = x[0];
704 const Double_t ka = par[0];
705 fgWStarfi->SetParameter(0, kb);
707 // root integration seems to fail
717 Double_t y = 2. * 208. * 208. * fgWStarfi->IntegralMultiple(2, al, bl, 0.001, err);
718 printf("WStaa: %.5e %.5e %.5e\n", b, y, err);
727 for (Int_t i = 0; i < fgkMCInts; i++)
730 const Double_t kphi = TMath::Pi() * gRandom->Rndm();
731 const Double_t kb1 = fgBMax * gRandom->Rndm();
732 y += fgWStarfi->Eval(kb1, kphi);
734 y *= 2. * TMath::Pi() * fgBMax / fgkMCInts; //fm^-2
735 y *= ka * ka * 0.1; //mbarn^-1
739 Double_t AliFastGlauber::WKParticipants(Double_t* x, Double_t* par)
742 // Kernel for number of participants
743 // as a function of r and phi
745 const Double_t kr1 = x[0];
746 const Double_t kphi = x[1];
747 const Double_t kb = par[0]; //fm
748 const Double_t ksig = par[1]; //mbarn
749 const Double_t ka = par[2]; //mass number
750 const Double_t kr2 = TMath::Sqrt(kr1*kr1 +kb*kb - 2.*kr1*kb*TMath::Cos(kphi));
751 const Double_t kxsi = fgWSta->Eval(kr2) * ksig * 0.1; //no units
753 Double_t y=(1-TMath::Power((1-xsi),aa))
756 Double_t sum = ka * kxsi;
758 for (Int_t i = 1; i <= ka; i++)
761 sum *= (-kxsi) * a / Float_t(i+1);
764 y *= kr1 * fgWSta->Eval(kr1);
768 Double_t AliFastGlauber::WParticipants(Double_t* x, Double_t* par)
771 // Number of Participants as
774 const Double_t kb = x[0];
775 const Double_t ksig = par[0]; //mbarn
776 const Double_t ka = par[1]; //mass number
777 fgWKParticipants->SetParameter(0, kb);
778 fgWKParticipants->SetParameter(1, ksig);
779 fgWKParticipants->SetParameter(2, ka);
785 for (Int_t i = 0; i < fgkMCInts; i++)
787 const Double_t kphi = TMath::Pi() * gRandom->Rndm();
788 const Double_t kb1 = fgBMax * gRandom->Rndm();
789 y += fgWKParticipants->Eval(kb1, kphi);
791 y *= 2. * ka * 2. * TMath::Pi() * fgBMax / fgkMCInts;
795 Double_t AliFastGlauber::WSgeo(Double_t* x, Double_t* par)
798 // Geometrical Cross-Section
799 // as a function of b
801 const Double_t kb = x[0]; //fm
802 const Double_t ksigNN = par[0]; //mbarn
803 const Double_t ktaa = fgWStaa->Eval(kb); //mbarn^-1
804 Double_t y = 2. * TMath::Pi() * kb * (1. - TMath::Exp(- ksigNN * ktaa));
808 Double_t AliFastGlauber::WSbinary(Double_t* x, Double_t* par)
811 // Number of binary hard collisions
812 // as a function of b
814 const Double_t kb = x[0]; //fm
815 const Double_t ksig = par[0]; //mbarn
816 const Double_t ktaa = fgWStaa->Eval(kb); //mbarn^-1
817 Double_t y = 2. * TMath::Pi() * kb * ksig * ktaa;
821 Double_t AliFastGlauber::WSN(Double_t* x, Double_t* /*par*/)
824 // Number of hard processes per event
825 // as a function of b
826 const Double_t kb = x[0];
827 Double_t y = fgWSbinary->Eval(kb)/fgWSgeo->Eval(kb);
831 Double_t AliFastGlauber::WEnergyDensity(Double_t* x, Double_t* par)
834 // Initial energy density
835 // as a function of the impact parameter
837 const Double_t kb = x[0];
838 const Double_t krA = par[0];
840 // Attention: area of transverse reaction zone in hard-sphere approximation !
841 const Double_t krA2=krA*krA;
842 const Double_t kb2=kb*kb;
843 const Double_t ksaa = (TMath::Pi() - 2. * TMath::ASin(kb/ 2./ krA)) * krA2
844 - kb * TMath::Sqrt(krA2 - kb2/ 4.); //fm^2
845 const Double_t ktaa = fgWStaa->Eval(kb); //mbarn^-1
846 Double_t y=ktaa/ksaa*10;
850 Double_t AliFastGlauber::WAlmond(Double_t* x, Double_t* par)
853 // Almond shaped interaction region
854 // as a function of cartesian x,y.
856 const Double_t kb = par[0];
857 const Double_t kxx = x[0] + kb/2.;
858 const Double_t kyy = x[1];
859 const Double_t kr1 = TMath::Sqrt(kxx*kxx + kyy*kyy);
860 const Double_t kphi = TMath::ATan2(kyy,kxx);
861 const Double_t kr2 = TMath::Sqrt(kr1*kr1 + kb*kb - 2.*kr1*kb*TMath::Cos(kphi));
863 // Interaction probability calculated as product of thicknesses
865 Double_t y = fgWSta->Eval(kr1) * fgWSta->Eval(kr2);
869 Double_t AliFastGlauber::WIntRadius(Double_t* x, Double_t* par)
872 // Average interaction density over radius
873 // at which interaction takes place
874 // as a function of radius
876 const Double_t kr = x[0];
877 const Double_t kb = par[0];
878 fgWAlmond->SetParameter(0, kb);
879 // Average over phi in small steps
880 const Double_t kdphi = 2. * TMath::Pi() / 100.;
883 for (Int_t i = 0; i < 100; i++) {
884 const Double_t kxx = kr * TMath::Cos(phi);
885 const Double_t kyy = kr * TMath::Sin(phi);
886 y += fgWAlmond->Eval(kxx,kyy);
889 // Result multiplied by Jacobian (2 pi r)
890 y *= 2. * TMath::Pi() * kr / 100.;
894 Double_t AliFastGlauber::WPathLength0(Double_t* x, Double_t* par)
897 // Path Length as a function of phi
898 // for interaction point fixed at (0,0)
899 // as a function of phi-direction
901 // Phi direction in Almond
902 const Double_t kphi0 = x[0];
903 const Double_t kb = par[0];
904 // Path Length definition
905 const Int_t kiopt = Int_t(par[1]);
907 // Step along radial direction phi
908 const Int_t kNp = 100; // Steps in r
909 const Double_t kDr = fgBMax/kNp;
913 for (Int_t i = 0; i < kNp; i++) {
915 // Transform into target frame
917 const Double_t kxx = r * TMath::Cos(kphi0) + kb / 2.;
918 const Double_t kyy = r * TMath::Sin(kphi0);
919 const Double_t kphi = TMath::ATan2(kyy, kxx);
920 const Double_t kr1 = TMath::Sqrt(kxx*kxx + kyy*kyy);
921 // Radius in projectile frame
922 const Double_t kr2 = TMath::Sqrt(kr1*kr1 + kb*kb - 2.*kr1*kb*TMath::Cos(kphi));
923 const Double_t ky = fgWSta->Eval(kr1) * fgWSta->Eval(kr2);
931 if (!kiopt) { // My length definition (is exact for hard disk)
932 if(w) y= 2. * rw / w;
934 const Double_t knorm=fgWSta->Eval(1e-4);
935 if(knorm) y = TMath::Sqrt(2. * rw * kDr / knorm / knorm);
940 Double_t AliFastGlauber::WPathLength(Double_t* x, Double_t* par)
943 // Path Length as a function of phi
944 // Interaction point from random distribution
945 // as a function of the phi-direction
946 const Double_t kphi0 = x[0];
947 const Double_t kb = par[0];
948 fgWAlmond->SetParameter(0, kb);
949 const Int_t kNpi = Int_t (par[1]); //Number of interactions
950 const Int_t kiopt = Int_t(par[2]); //Path Length definition
955 const Int_t kNp = 100;
956 const Double_t kDr = fgBMax/Double_t(kNp);
957 Double_t l = 0.; // Path length
958 for (Int_t in = 0; in < kNpi; in ++) {
963 fgWAlmond->GetRandom2(x0, y0);
965 const Double_t kr0 = TMath::Sqrt(x0*x0 + y0*y0);
966 const Int_t knps = Int_t ((fgBMax - kr0)/kDr) - 1;
970 for (Int_t i = 0; (i < knps ); i++) {
971 // Transform into target frame
972 const Double_t kxx = x0 + r * TMath::Cos(kphi0) + kb / 2.;
973 const Double_t kyy = y0 + r * TMath::Sin(kphi0);
974 const Double_t kphi = TMath::ATan2(kyy, kxx);
975 const Double_t kr1 = TMath::Sqrt(kxx*kxx + kyy*kyy);
976 // Radius in projectile frame
977 const Double_t kr2 = TMath::Sqrt(kr1*kr1 + kb*kb - 2.*kr1*kb*TMath::Cos(kphi));
978 const Double_t ky = fgWSta->Eval(kr1) * fgWSta->Eval(kr2);
984 // Average over interactions
986 if(w) l += (2. * rw / w);
988 const Double_t knorm=fgWSta->Eval(1e-4);
989 if(knorm) l+= 2. * rw * kDr / knorm / knorm;
996 ret=TMath::Sqrt( l / kNpi);
1000 Double_t AliFastGlauber::CrossSection(Double_t b1, Double_t b2) const
1003 // Return the geometrical cross-section integrated from b1 to b2
1005 return fgWSgeo->Integral(b1, b2)*10.; //mbarn
1008 Double_t AliFastGlauber::HardCrossSection(Double_t b1, Double_t b2) const
1011 // Return the hard cross-section integrated from b1 to b2
1013 return fgWSbinary->Integral(b1, b2)*10.; //mbarn
1016 Double_t AliFastGlauber::FractionOfHardCrossSection(Double_t b1, Double_t b2) const
1019 // Return fraction of hard cross-section integrated from b1 to b2
1021 return fgWSbinary->Integral(b1, b2)/fgWSbinary->Integral(0., 100.);
1024 Double_t AliFastGlauber::NHard(Double_t b1, Double_t b2) const
1027 // Number of binary hard collisions
1028 // as a function of b (nucl/ex/0302016 eq. 19)
1030 const Double_t kshard=HardCrossSection(b1,b2);
1031 const Double_t ksgeo=CrossSection(b1,b2);
1033 return kshard/ksgeo;
1037 Double_t AliFastGlauber::Binaries(Double_t b) const
1040 // Return number of binary hard collisions normalized to 1 at b=0
1043 return fgWSN->Eval(b)/fgWSN->Eval(1e-4);
1046 Double_t AliFastGlauber::MeanOverlap(Double_t b1, Double_t b2)
1049 // Calculate the mean overlap for impact parameter range b1 .. b2
1055 while (b < b2-0.005) {
1056 Double_t nc = GetNumberOfCollisions(b);
1057 sum += 10. * fgWStaa->Eval(b) * fgWSgeo->Eval(b) * 0.01 / (1. - TMath::Exp(-nc));
1058 sumc += 10. * fgWSgeo->Eval(b) * 0.01;
1061 return (sum / CrossSection(b1, b2));
1065 Double_t AliFastGlauber::MeanNumberOfCollisionsPerEvent(Double_t b1, Double_t b2)
1068 // Calculate the mean number of collisions per event for impact parameter range b1 .. b2
1074 while (b < b2-0.005) {
1075 Double_t nc = GetNumberOfCollisions(b);
1076 sum += nc / (1. - TMath::Exp(-nc)) * 10. * fgWSgeo->Eval(b) * 0.01;
1077 sumc += 10. * fgWSgeo->Eval(b) * 0.01;
1080 return (sum / CrossSection(b1, b2));
1084 Double_t AliFastGlauber::GetNumberOfBinaries(Double_t b) const
1087 // Return number of binary hard collisions at b
1090 return fgWSN->Eval(b);
1093 Double_t AliFastGlauber::Participants(Double_t b) const
1096 // Return the number of participants normalized to 1 at b=0
1099 return (fgWParticipants->Eval(b)/fgWParticipants->Eval(1e-4));
1102 Double_t AliFastGlauber::GetNumberOfParticipants(Double_t b) const
1105 // Return the number of participants for impact parameter b
1108 return (fgWParticipants->Eval(b));
1111 Double_t AliFastGlauber::GetNumberOfCollisions(Double_t b) const
1114 // Return the number of collisions for impact parameter b
1117 return (fgWStaa->Eval(b)*fSigmaNN);
1120 Double_t AliFastGlauber::GetNumberOfCollisionsPerEvent(Double_t b) const
1123 // Return the number of collisions per event (at least one collision)
1124 // for impact parameter b
1126 Double_t n = GetNumberOfCollisions(b);
1128 return (n / (1. - TMath::Exp(- n)));
1134 void AliFastGlauber::SimulateTrigger(Int_t n)
1137 // Simulates Trigger
1139 TH1F* mbtH = new TH1F("mbtH", "MB Trigger b-Distribution", 100, 0., 20.);
1140 TH1F* hdtH = new TH1F("hdtH", "Hard Trigger b-Distribution", 100, 0., 20.);
1141 TH1F* mbmH = new TH1F("mbmH", "MB Trigger Multiplicity Distribution", 100, 0., 8000.);
1142 TH1F* hdmH = new TH1F("hdmH", "Hard Trigger Multiplicity Distribution", 100, 0., 8000.);
1144 mbtH->SetXTitle("b [fm]");
1145 hdtH->SetXTitle("b [fm]");
1146 mbmH->SetXTitle("Multiplicity");
1147 hdmH->SetXTitle("Multiplicity");
1149 TCanvas *c0 = new TCanvas("c0","Trigger Simulation",400,10,600,700);
1151 TCanvas *c1 = new TCanvas("c1","Trigger Simulation",400,10,600,700);
1157 for (Int_t iev = 0; iev < n; iev++)
1160 GetRandom(b, p, mult);
1163 mbmH->Fill(mult, 1.);
1164 hdmH->Fill(mult, p);
1180 void AliFastGlauber::GetRandom(Float_t& b, Float_t& p, Float_t& mult)
1183 // Gives back a random impact parameter, hard trigger probability and multiplicity
1185 b = fgWSgeo->GetRandom();
1186 const Float_t kmu = fgWSN->Eval(b);
1187 p = 1.-TMath::Exp(-kmu);
1188 mult = 6000./fgWSN->Eval(1.) * kmu;
1191 void AliFastGlauber::GetRandom(Int_t& bin, Bool_t& hard)
1194 // Gives back a random impact parameter bin, and hard trigger decission
1196 const Float_t kb = fgWSgeo->GetRandom();
1197 const Float_t kmu = fgWSN->Eval(kb) * fSigmaHard;
1198 const Float_t kp = 1.-TMath::Exp(-kmu);
1201 } else if (kb < 8.6) {
1203 } else if (kb < 11.2) {
1205 } else if (kb < 13.2) {
1207 } else if (kb < 15.0) {
1213 const Float_t kr = gRandom->Rndm();
1214 if (kr < kp) hard = kTRUE;
1217 Double_t AliFastGlauber::GetRandomImpactParameter(Double_t bmin, Double_t bmax)
1220 // Gives back a random impact parameter in the range bmin .. bmax
1223 while(b < bmin || b > bmax)
1224 b = fgWSgeo->GetRandom();
1228 void AliFastGlauber::StoreFunctions() const
1231 // Store in file functions
1233 TFile* ff = new TFile(fName.Data(),"recreate");
1234 fgWStaa->Write("WStaa");
1235 fgWParticipants->Write("WParticipants");
1240 //=================== Added by A. Dainese 11/02/04 ===========================
1242 void AliFastGlauber::StoreAlmonds() const
1246 // 40 almonds for b = (0.25+k*0.5) fm (k=0->39)
1248 Char_t almondName[100];
1249 TFile* ff = new TFile(fName.Data(),"update");
1250 for(Int_t k=0; k<40; k++) {
1251 sprintf(almondName,"WAlmondFixedB%d",k);
1252 Double_t b = 0.25+k*0.5;
1253 Info("StoreAlmonds"," b = %f\n",b);
1254 fgWAlmond->SetParameter(0,b);
1255 fgWAlmond->Write(almondName);
1261 void AliFastGlauber::SetCentralityClass(Double_t xsecFrLow,Double_t xsecFrUp)
1264 // Set limits of centrality class as fractions
1265 // of the geomtrical cross section
1267 if(xsecFrLow>1. || xsecFrUp>1. || xsecFrLow>xsecFrUp) {
1268 Error("SetCentralityClass", "Please set 0 <= xsecFrLow <= xsecFrUp <= 1\n");
1272 Double_t bLow=0.,bUp=0.;
1274 const Double_t knorm=fgWSgeo->Integral(0.,100.);
1275 while(xsecFr<xsecFrLow) {
1276 xsecFr = fgWSgeo->Integral(0.,bLow)/knorm;
1280 while(xsecFr<xsecFrUp) {
1281 xsecFr = fgWSgeo->Integral(0.,bUp)/knorm;
1285 Info("SetCentralityClass", "Centrality class: %4.2f-%4.2f; %4.1f < b < %4.1f fm",
1286 xsecFrLow,xsecFrUp,bLow,bUp);
1287 fgWSbinary->SetRange(bLow,bUp);
1293 void AliFastGlauber::GetRandomBHard(Double_t& b)
1296 // Get random impact parameter according to distribution of
1297 // hard (binary) cross-section, in the range defined by the centrality class
1299 b = fgWSbinary->GetRandom();
1300 Int_t bin = 2*(Int_t)b;
1301 if( (b-(Int_t)b) > 0.5) bin++;
1302 fgWAlmondCurrent = fgWAlmondFixedB[bin];
1306 void AliFastGlauber::GetRandomXY(Double_t& x,Double_t& y)
1309 // Get random position of parton production point according to
1310 // product of thickness functions
1312 fgWAlmondCurrent->GetRandom2(x,y);
1316 void AliFastGlauber::GetRandomPhi(Double_t& phi)
1319 // Get random parton azimuthal propagation direction
1321 phi = 2.*TMath::Pi()*gRandom->Rndm();
1325 Double_t AliFastGlauber::CalculateLength(Double_t b,Double_t x0,Double_t y0,Double_t phi0)
1328 // Calculate path length for a parton with production point (x0,y0)
1329 // and propagation direction (ux=cos(phi0),uy=sin(phi0))
1330 // in a collision with impact parameter b
1333 // number of steps in l
1334 const Int_t kNp = 100;
1335 const Double_t kDl = fgBMax/Double_t(kNp);
1341 // ell = 2 * \int_0^\infty dl*l*(T_A*T_B)(x0+l*ux,y0+l*uy) /
1342 // \int_0^\infty dl*(T_A*T_B)(x0+l*ux,y0+l*uy)
1346 const Double_t kr0 = TMath::Sqrt(x0*x0 + y0*y0);
1347 const Int_t knps = Int_t ((fgBMax - kr0)/kDl) - 1;
1349 Double_t integral1 = 0.;
1350 Double_t integral2 = 0.;
1352 for (Int_t i = 0; i < knps; i++) {
1354 // Transform into target frame
1355 const Double_t kxx = x0 + l * TMath::Cos(phi0) + b / 2.;
1356 const Double_t kyy = y0 + l * TMath::Sin(phi0);
1357 const Double_t kphi = TMath::ATan2(kyy, kxx);
1358 const Double_t kr1 = TMath::Sqrt(kxx*kxx + kyy*kyy);
1359 // Radius in projectile frame
1360 const Double_t kr2 = TMath::Sqrt(kr1*kr1 + b*b - 2.*kr1*b*TMath::Cos(kphi));
1361 const Double_t kprodTATB = fgWSta->Eval(kr1) * fgWSta->Eval(kr2);
1363 integral1 += kprodTATB * l * kDl;
1364 integral2 += kprodTATB * kDl;
1370 ell = (2. * integral1 / integral2);
1372 } else if(fEllDef==2) {
1376 // ell = \int_0^\infty dl*
1377 // \Theta((T_A*T_B)(x0+l*ux,y0+l*uy)-0.5*(T_A*T_B)(0,0))
1381 const Double_t kr0 = TMath::Sqrt(x0*x0 + y0*y0);
1382 const Int_t knps = Int_t ((fgBMax - kr0)/kDl) - 1;
1383 const Double_t kprodTATBHalfMax = 0.5*fgWAlmondCurrent->Eval(0.,0.);
1386 Double_t integral = 0.;
1387 for (Int_t i = 0; i < knps; i++) {
1388 // Transform into target frame
1389 const Double_t kxx = x0 + l * TMath::Cos(phi0) + b / 2.;
1390 const Double_t kyy = y0 + l * TMath::Sin(phi0);
1391 const Double_t kphi = TMath::ATan2(kyy, kxx);
1392 const Double_t kr1 = TMath::Sqrt(kxx*kxx + kyy*kyy);
1393 // Radius in projectile frame
1394 const Double_t kr2 = TMath::Sqrt(kr1*kr1 + b*b - 2.*kr1*b*TMath::Cos(kphi));
1395 const Double_t kprodTATB = fgWSta->Eval(kr1) * fgWSta->Eval(kr2);
1396 if(kprodTATB>kprodTATBHalfMax) integral += kDl;
1399 Double_t ell = integral;
1402 Error("CalculateLength","Wrong length definition setting: %d !\n",fEllDef);
1407 void AliFastGlauber::GetLengthAndPhi(Double_t& ell,Double_t& phi,Double_t b)
1410 // Return length from random b, x0, y0, phi0
1413 Double_t x0,y0,phi0;
1414 if(b<0.) GetRandomBHard(b);
1418 ell = CalculateLength(b,x0,y0,phi0);
1422 void AliFastGlauber::GetLength(Double_t& ell,Double_t b)
1425 // Return length from random b, x0, y0, phi0
1428 GetLengthAndPhi(ell,phi,b);
1432 void AliFastGlauber::GetLengthsBackToBackAndPhi(Double_t& ell1,Double_t& ell2,Double_t &phi,Double_t b)
1435 // Return 2 lengths back to back from random b, x0, y0, phi0
1438 Double_t x0,y0,phi0;
1439 if(b<0.) GetRandomBHard(b);
1442 const Double_t kphi0plusPi = phi0+TMath::Pi();
1444 ell1 = CalculateLength(b,x0,y0,phi0);
1445 ell2 = CalculateLength(b,x0,y0,kphi0plusPi);
1449 void AliFastGlauber::GetLengthsBackToBack(Double_t& ell1,Double_t& ell2,
1453 // Return 2 lengths back to back from random b, x0, y0, phi0
1456 GetLengthsBackToBackAndPhi(ell1,ell2,phi,b);
1460 void AliFastGlauber::GetLengthsForPythia(Int_t n,Double_t* phi,Double_t* ell, Double_t b)
1463 // Returns lenghts for n partons with azimuthal angles phi[n]
1464 // from random b, x0, y0
1467 if(b < 0.) GetRandomBHard(b);
1469 for(Int_t i = 0; i< n; i++) ell[i] = CalculateLength(b,x0,y0,phi[i]);
1473 void AliFastGlauber::PlotBDistr(Int_t n)
1476 // Plot distribution of n impact parameters
1479 TH1F *hB = new TH1F("hB","dN/db",100,0,fgBMax);
1480 hB->SetXTitle("b [fm]");
1481 hB->SetYTitle("dN/db [a.u.]");
1482 hB->SetFillColor(3);
1483 for(Int_t i=0; i<n; i++) {
1487 TCanvas *cB = new TCanvas("cB","Impact parameter distribution",0,0,500,500);
1493 void AliFastGlauber::PlotLengthDistr(Int_t n,Bool_t save,const char *fname)
1496 // Plot length distribution
1499 TH1F *hEll = new TH1F("hEll","Length distribution",64,-0.5,15);
1500 hEll->SetXTitle("Transverse path length, L [fm]");
1501 hEll->SetYTitle("Probability");
1502 hEll->SetFillColor(2);
1503 for(Int_t i=0; i<n; i++) {
1507 hEll->Scale(1/(Double_t)n);
1508 TCanvas *cL = new TCanvas("cL","Length distribution",0,0,500,500);
1513 TFile *f = new TFile(fname,"recreate");
1520 void AliFastGlauber::PlotLengthB2BDistr(Int_t n,Bool_t save,const char *fname)
1523 // Plot lengths back-to-back distributions
1526 TH2F *hElls = new TH2F("hElls","Lengths back-to-back",100,0,15,100,0,15);
1527 hElls->SetXTitle("Transverse path length, L [fm]");
1528 hElls->SetYTitle("Transverse path length, L [fm]");
1529 for(Int_t i=0; i<n; i++) {
1530 GetLengthsBackToBack(ell1,ell2);
1531 hElls->Fill(ell1,ell2);
1533 hElls->Scale(1/(Double_t)n);
1534 TCanvas *cLs = new TCanvas("cLs","Length back-to-back distribution",0,0,500,500);
1535 gStyle->SetPalette(1,0);
1537 hElls->Draw("col,Z");
1539 TFile *f = new TFile(fname,"recreate");
1546 void AliFastGlauber::PlotAlmonds() const
1549 // Plot almonds for some impact parameters
1551 TCanvas *c = new TCanvas("c","Almonds",0,0,500,500);
1552 gStyle->SetPalette(1,0);
1555 fgWAlmondFixedB[0]->Draw("cont1");
1557 fgWAlmondFixedB[10]->Draw("cont1");
1559 fgWAlmondFixedB[20]->Draw("cont1");
1561 fgWAlmondFixedB[30]->Draw("cont1");
1565 //=================== Added by A. Dainese 05/03/04 ===========================
1567 void AliFastGlauber::CalculateI0I1(Double_t& integral0,Double_t& integral1,
1568 Double_t b,Double_t x0,Double_t y0,
1569 Double_t phi0,Double_t ellCut) const
1572 // Calculate integrals:
1573 // integral0 = \int_0^ellCut dl*(T_A*T_B)(x0+l*ux,y0+l*uy)
1574 // integral1 = \int_0^ellCut dl*l*(T_A*T_B)(x0+l*ux,y0+l*uy)
1576 // for a parton with production point (x0,y0)
1577 // and propagation direction (ux=cos(phi0),uy=sin(phi0))
1578 // in a collision with impact parameter b
1581 // number of steps in l
1582 const Int_t kNp = 100;
1583 const Double_t kDl = fgBMax/Double_t(kNp);
1586 const Double_t kr0 = TMath::Sqrt(x0 * x0 + y0 * y0);
1587 const Int_t knps = Int_t ((fgBMax - kr0)/kDl) - 1;
1594 while((i < knps) && (l < ellCut)) {
1595 // Transform into target frame
1596 const Double_t kxx = x0 + l * TMath::Cos(phi0) + b / 2.;
1597 const Double_t kyy = y0 + l * TMath::Sin(phi0);
1598 const Double_t kphi = TMath::ATan2(kyy, kxx);
1599 const Double_t kr1 = TMath::Sqrt(kxx*kxx + kyy*kyy);
1600 // Radius in projectile frame
1601 const Double_t kr2 = TMath::Sqrt(kr1*kr1 + b*b - 2.*kr1*b*TMath::Cos(kphi));
1602 const Double_t kprodTATB = fgWSta->Eval(kr1) * fgWSta->Eval(kr2);
1603 integral0 += kprodTATB * kDl;
1604 integral1 += kprodTATB * l * kDl;
1611 void AliFastGlauber::GetI0I1AndPhi(Double_t& integral0,Double_t& integral1,
1613 Double_t ellCut,Double_t b)
1616 // Return I0 and I1 from random b, x0, y0, phi0
1619 Double_t x0,y0,phi0;
1620 if(b<0.) GetRandomBHard(b);
1624 CalculateI0I1(integral0,integral1,b,x0,y0,phi0,ellCut);
1628 void AliFastGlauber::GetI0I1(Double_t& integral0,Double_t& integral1,
1629 Double_t ellCut,Double_t b)
1632 // Return I0 and I1 from random b, x0, y0, phi0
1635 GetI0I1AndPhi(integral0,integral1,phi,ellCut,b);
1639 void AliFastGlauber::GetI0I1BackToBackAndPhi(Double_t& integral01,Double_t& integral11,
1640 Double_t& integral02,Double_t& integral12,
1642 Double_t ellCut,Double_t b)
1645 // Return 2 pairs of I0 and I1 back to back from random b, x0, y0, phi0
1648 Double_t x0,y0,phi0;
1649 if(b<0.) GetRandomBHard(b);
1653 const Double_t kphi0plusPi = phi0+TMath::Pi();
1654 CalculateI0I1(integral01,integral11,b,x0,y0,phi0,ellCut);
1655 CalculateI0I1(integral02,integral12,b,x0,y0,kphi0plusPi,ellCut);
1659 void AliFastGlauber::GetI0I1BackToBackAndPhiAndXY(Double_t& integral01,Double_t& integral11,
1660 Double_t& integral02,Double_t& integral12,
1661 Double_t& phi,Double_t &x,Double_t &y,
1662 Double_t ellCut,Double_t b)
1665 // Return 2 pairs of I0 and I1 back to back from random b, x0, y0, phi0
1668 Double_t x0,y0,phi0;
1669 if(b<0.) GetRandomBHard(b);
1672 phi = phi0; x=x0; y=y0;
1673 const Double_t kphi0plusPi = phi0+TMath::Pi();
1674 CalculateI0I1(integral01,integral11,b,x0,y0,phi0,ellCut);
1675 CalculateI0I1(integral02,integral12,b,x0,y0,kphi0plusPi,ellCut);
1679 void AliFastGlauber::GetI0I1BackToBack(Double_t& integral01,Double_t& integral11,
1680 Double_t& integral02,Double_t& integral12,
1681 Double_t ellCut,Double_t b)
1684 // Return 2 pairs of I0 and I1 back to back from random b, x0, y0, phi0
1687 GetI0I1BackToBackAndPhi(integral01,integral11,integral02,integral12,
1692 void AliFastGlauber::GetI0I1ForPythia(Int_t n,Double_t* phi,
1693 Double_t* integral0,Double_t* integral1,
1694 Double_t ellCut,Double_t b)
1697 // Returns I0 and I1 pairs for n partons with azimuthal angles phi[n]
1698 // from random b, x0, y0
1701 if(b<0.) GetRandomBHard(b);
1703 for(Int_t i=0; i<n; i++)
1704 CalculateI0I1(integral0[i],integral1[i],b,x0,y0,phi[i],ellCut);
1708 void AliFastGlauber::GetI0I1ForPythiaAndXY(Int_t n,Double_t* phi,
1709 Double_t* integral0,Double_t* integral1,
1710 Double_t &x,Double_t& y,
1711 Double_t ellCut,Double_t b)
1714 // Returns I0 and I1 pairs for n partons with azimuthal angles phi[n]
1715 // from random b, x0, y0 and return x0,y0
1718 if(b<0.) GetRandomBHard(b);
1720 for(Int_t i=0; i<n; i++)
1721 CalculateI0I1(integral0[i],integral1[i],b,x0,y0,phi[i],ellCut);
1727 void AliFastGlauber::PlotI0I1Distr(Int_t n,Double_t ellCut,
1728 Bool_t save,const char *fname)
1731 // Plot I0-I1 distribution
1734 TH2F *hI0I1s = new TH2F("hI0I1s","I_{0} versus I_{1}",1000,0,0.001,1000,0,0.01);
1735 hI0I1s->SetXTitle("I_{0} [fm^{-3}]");
1736 hI0I1s->SetYTitle("I_{1} [fm^{-2}]");
1738 TH1F *hI0 = new TH1F("hI0","I_{0} = #hat{q}L / k",
1740 hI0->SetXTitle("I_{0} [fm^{-3}]");
1741 hI0->SetYTitle("Probability");
1742 hI0->SetFillColor(3);
1743 TH1F *hI1 = new TH1F("hI1","I_{1} = #omega_{c} / k",
1745 hI1->SetXTitle("I_{1} [fm^{-2}]");
1746 hI1->SetYTitle("Probability");
1747 hI1->SetFillColor(4);
1748 TH1F *h2 = new TH1F("h2","2 I_{1}^{2}/I_{0} = R / k",
1750 h2->SetXTitle("2 I_{1}^{2}/I_{0} [fm^{-1}]");
1751 h2->SetYTitle("Probability");
1752 h2->SetFillColor(2);
1753 TH1F *h3 = new TH1F("h3","2 I_{1}/I_{0} = L",
1755 h3->SetXTitle("2 I_{1}/I_{0} [fm]");
1756 h3->SetYTitle("Probability");
1757 h3->SetFillColor(5);
1758 TH1F *h4 = new TH1F("h4","I_{0}^{2}/(2 I_{1}) = #hat{q} / k",
1760 h4->SetXTitle("I_{0}^{2}/(2 I_{1}) [fm^{-4}]");
1761 h4->SetYTitle("Probability");
1762 h4->SetFillColor(7);
1764 for(Int_t i=0; i<n; i++) {
1765 GetI0I1(i0,i1,ellCut);
1766 hI0I1s->Fill(i0,i1);
1769 h2->Fill(2.*i1*i1/i0);
1771 h4->Fill(i0*i0/2./i1);
1773 hI0->Scale(1/(Double_t)n);
1774 hI1->Scale(1/(Double_t)n);
1775 h2->Scale(1/(Double_t)n);
1776 h3->Scale(1/(Double_t)n);
1777 h4->Scale(1/(Double_t)n);
1778 hI0I1s->Scale(1/(Double_t)n);
1780 TCanvas *cI0I1 = new TCanvas("cI0I1","I0 and I1",0,0,900,700);
1793 gStyle->SetPalette(1,0);
1794 hI0I1s->Draw("col,Z");
1797 TFile *f = new TFile(fname,"recreate");
1809 void AliFastGlauber::PlotI0I1B2BDistr(Int_t n,Double_t ellCut,
1810 Bool_t save,const char *fname)
1813 // Plot I0-I1 back-to-back distributions
1815 Double_t i01,i11,i02,i12;
1816 TH2F *hI0s = new TH2F("hI0s","I_{0}'s back-to-back",100,0,100,100,0,100);
1817 hI0s->SetXTitle("I_{0} [fm^{-3}]");
1818 hI0s->SetYTitle("I_{0} [fm^{-3}]");
1819 TH2F *hI1s = new TH2F("hI1s","I_{1}'s back-to-back",100,0,100,100,0,100);
1820 hI1s->SetXTitle("I_{1} [fm^{-2}]");
1821 hI1s->SetYTitle("I_{1} [fm^{-2}]");
1823 for(Int_t i=0; i<n; i++) {
1824 GetI0I1BackToBack(i01,i11,i02,i12,ellCut);
1825 hI0s->Fill(i01,i02);
1826 hI1s->Fill(i11,i12);
1828 hI0s->Scale(1/(Double_t)n);
1829 hI1s->Scale(1/(Double_t)n);
1831 TCanvas *cI0I1s = new TCanvas("cI0I1s","I0 and I1 back-to-back distributions",0,0,800,400);
1832 gStyle->SetPalette(1,0);
1833 cI0I1s->Divide(2,1);
1835 hI0s->Draw("col,Z");
1837 hI1s->Draw("col,Z");
1840 TFile *f = new TFile(fname,"recreate");
1848 AliFastGlauber& AliFastGlauber::operator=(const AliFastGlauber& rhs)
1850 // Assignment operator
1855 void AliFastGlauber::Copy(TObject&) const
1860 Fatal("Copy","Not implemented!\n");