]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - FASTSIM/AliFastGlauber.cxx
Moving the functions used to initialize TF1 and TF2 to the pivate part of the class
[u/mrichter/AliRoot.git] / FASTSIM / AliFastGlauber.cxx
index 20d9983016416240735ad01b9063526c802e3955..804f0e52da55b5b21f91f0c7d1089d503487f2a8 100644 (file)
@@ -148,7 +148,7 @@ void AliFastGlauber::Init(Int_t mode)
   fgWSb->SetParameter(2, fWSw);
   fgWSb->SetParameter(3, fWSn);
 
-  fgWSbz = new TF2("WSbz", WSbz, 0, fgBMax, 4);
+  fgWSbz = new TF2("WSbz", WSbz, 0, fgBMax, 0, fgBMax, 4);
   fgWSbz->SetParameter(0, fWSr0);
   fgWSbz->SetParameter(1, fWSd);
   fgWSbz->SetParameter(2, fWSw);
@@ -972,11 +972,24 @@ Double_t AliFastGlauber::HardCrossSection(Double_t b1, Double_t b2) const
 Double_t AliFastGlauber::FractionOfHardCrossSection(Double_t b1, Double_t b2) const
 {
   //
-  // Return raction of hard cross-section integrated from b1 to b2 
+  // Return fraction of hard cross-section integrated from b1 to b2 
   //
   return fgWSbinary->Integral(b1, b2)/fgWSbinary->Integral(0., 100.);
 }
 
+Double_t AliFastGlauber::NHard(Double_t b1, Double_t b2) const
+{
+  //
+  //  Number of binary hard collisions 
+  //  as a function of b (nucl/ex/0302016 eq. 19)
+  //
+  const Double_t kshard=HardCrossSection(b1,b2);
+  const Double_t ksgeo=CrossSection(b1,b2); 
+  if(ksgeo>0)
+    return kshard/ksgeo;
+  else return -1; 
+}
+
 Double_t AliFastGlauber::Binaries(Double_t b) const
 {
   //
@@ -1022,6 +1035,20 @@ Double_t AliFastGlauber::GetNumberOfCollisions(Double_t  b) const
   return (fgWStaa->Eval(b)*fSigmaNN);
 }
 
+Double_t AliFastGlauber::GetNumberOfCollisionsPerEvent(Double_t  b) const
+{
+  //
+  // Return the number of collisions per event (at least one collision)
+  // for impact parameter b
+  //
+    Double_t n = GetNumberOfCollisions(b);
+    if (n > 0.) {
+       return (n / (1. - TMath::Exp(- n)));
+    } else {
+       return (0.);
+    }
+}
+
 void AliFastGlauber::SimulateTrigger(Int_t n)
 {
   //
@@ -1381,7 +1408,7 @@ void AliFastGlauber::PlotBDistr(Int_t n)
   return;
 }
 
-void AliFastGlauber::PlotLengthDistr(Int_t n,Bool_t save,Char_t *fname)
+void AliFastGlauber::PlotLengthDistr(Int_t n,Bool_t save,const char *fname)
 {
   //
   // Plot length distribution
@@ -1408,7 +1435,7 @@ void AliFastGlauber::PlotLengthDistr(Int_t n,Bool_t save,Char_t *fname)
   return;
 }
 
-void AliFastGlauber::PlotLengthB2BDistr(Int_t n,Bool_t save,Char_t *fname)
+void AliFastGlauber::PlotLengthB2BDistr(Int_t n,Bool_t save,const char *fname)
 {
   //
   // Plot lengths back-to-back distributions
@@ -1547,6 +1574,26 @@ void AliFastGlauber::GetI0I1BackToBackAndPhi(Double_t& integral01,Double_t& inte
   return;
 }
 
+void AliFastGlauber::GetI0I1BackToBackAndPhiAndXY(Double_t& integral01,Double_t& integral11,
+                                                 Double_t& integral02,Double_t& integral12,
+                                                 Double_t& phi,Double_t &x,Double_t &y,
+                                                 Double_t ellCut,Double_t b)
+{
+  //
+  // Return 2 pairs of I0 and I1 back to back from random b, x0, y0, phi0 
+  // Return also phi0
+  //
+  Double_t x0,y0,phi0;
+  if(b<0.) GetRandomBHard(b);
+  GetRandomXY(x0,y0);
+  GetRandomPhi(phi0);
+  phi = phi0; x=x0; y=y0;
+  const Double_t kphi0plusPi = phi0+TMath::Pi();
+  CalculateI0I1(integral01,integral11,b,x0,y0,phi0,ellCut);
+  CalculateI0I1(integral02,integral12,b,x0,y0,kphi0plusPi,ellCut);
+  return;
+}
+
 void AliFastGlauber::GetI0I1BackToBack(Double_t& integral01,Double_t& integral11,
                                       Double_t& integral02,Double_t& integral12,
                                       Double_t ellCut,Double_t b)
@@ -1576,8 +1623,27 @@ void AliFastGlauber::GetI0I1ForPythia(Int_t n,Double_t* phi,
   return;
 }
 
+void AliFastGlauber::GetI0I1ForPythiaAndXY(Int_t n,Double_t* phi,
+                                     Double_t* integral0,Double_t* integral1,
+                                     Double_t &x,Double_t& y,
+                                     Double_t ellCut,Double_t b)
+{
+  //
+  // Returns I0 and I1 pairs for n partons with azimuthal angles phi[n] 
+  // from random b, x0, y0 and return x0,y0
+  //
+  Double_t x0,y0;
+  if(b<0.) GetRandomBHard(b);
+  GetRandomXY(x0,y0);
+  for(Int_t i=0; i<n; i++) 
+    CalculateI0I1(integral0[i],integral1[i],b,x0,y0,phi[i],ellCut);
+  x=x0;
+  y=y0;
+  return;
+}
+
 void AliFastGlauber::PlotI0I1Distr(Int_t n,Double_t ellCut,
-                                  Bool_t save,Char_t *fname)
+                                  Bool_t save,const char *fname)
 {
   //
   // Plot I0-I1 distribution
@@ -1659,7 +1725,7 @@ void AliFastGlauber::PlotI0I1Distr(Int_t n,Double_t ellCut,
 }
 
 void AliFastGlauber::PlotI0I1B2BDistr(Int_t n,Double_t ellCut,
-                                     Bool_t save,Char_t *fname)
+                                     Bool_t save,const char *fname)
 {
   //
   // Plot I0-I1 back-to-back distributions