]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - FASTSIM/AliFastGlauber.cxx
Ions have to be added to root particle database.
[u/mrichter/AliRoot.git] / FASTSIM / AliFastGlauber.cxx
index 6088042850e755f567e5aa6995fce79f1505ded7..7e74953784b58fd4855b6aa47916c5b9ebfca467 100644 (file)
@@ -82,7 +82,7 @@ TF1*    AliFastGlauber::fgWIntRadius     = NULL;
 TF2*    AliFastGlauber::fgWKParticipants = NULL; 
 TF1*    AliFastGlauber::fgWParticipants  = NULL; 
 TF2*    AliFastGlauber::fgWAlmondCurrent = NULL;    
-TF2     AliFastGlauber::fgWAlmondFixedB[40]; 
+TF2*    AliFastGlauber::fgWAlmondFixedB[40]; 
 const Int_t AliFastGlauber::fgkMCInts = 100000;
 Int_t AliFastGlauber::fgCounter = 0;       
 
@@ -103,6 +103,7 @@ AliFastGlauber::AliFastGlauber() : fName()
 AliFastGlauber::~AliFastGlauber()
 {
   fgCounter--;
+  for(Int_t k=0; k<40; k++) delete fgWAlmondFixedB[k];
   //if(fgCounter==0) Reset();
 }
 
@@ -148,7 +149,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);
@@ -247,7 +248,7 @@ void AliFastGlauber::Init(Int_t mode)
     for(Int_t k=0; k<40; k++) {
       sprintf(almondName,"WAlmondFixedB%d",k);
       fgWAlmondCurrent = (TF2*)ff->Get(almondName);
-      new(&fgWAlmondFixedB[k]) TF2(*fgWAlmondCurrent);
+      fgWAlmondFixedB[k] = fgWAlmondCurrent;
     }
     delete ff;
   }
@@ -972,11 +973,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 +1036,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)
 {
   //
@@ -1190,7 +1218,7 @@ void AliFastGlauber::GetRandomBHard(Double_t& b)
   b = fgWSbinary->GetRandom();
   Int_t bin = 2*(Int_t)b;
   if( (b-(Int_t)b) > 0.5) bin++;
-  fgWAlmondCurrent = &fgWAlmondFixedB[bin]; 
+  fgWAlmondCurrent = fgWAlmondFixedB[bin]; 
   return;
 }
 
@@ -1381,7 +1409,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 +1436,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
@@ -1443,13 +1471,13 @@ void AliFastGlauber::PlotAlmonds() const
   gStyle->SetPalette(1,0);
   c->Divide(2,2);
   c->cd(1);
-  fgWAlmondFixedB[0].Draw("cont1");
+  fgWAlmondFixedB[0]->Draw("cont1");
   c->cd(2);
-  fgWAlmondFixedB[10].Draw("cont1");
+  fgWAlmondFixedB[10]->Draw("cont1");
   c->cd(3);
-  fgWAlmondFixedB[20].Draw("cont1");
+  fgWAlmondFixedB[20]->Draw("cont1");
   c->cd(4);
-  fgWAlmondFixedB[30].Draw("cont1");
+  fgWAlmondFixedB[30]->Draw("cont1");
   return;
 }
 
@@ -1616,7 +1644,7 @@ void AliFastGlauber::GetI0I1ForPythiaAndXY(Int_t n,Double_t* phi,
 }
 
 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
@@ -1698,7 +1726,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