]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - EVGEN/AliGenSlowNucleons.cxx
Initialization of pointers in default constructor.
[u/mrichter/AliRoot.git] / EVGEN / AliGenSlowNucleons.cxx
index 3599590cfb5c0a3c86629cfbdaa4f883238b3849..b408a220538b0b27760c00d18bc03643fcfe0898 100644 (file)
  * provided "as is" without express or implied warranty.                  *
  **************************************************************************/
 
-/*
-$Log$
-Revision 1.3  2003/06/26 08:45:22  morsch
-- Vertex using Vertex() method.
-- Use member data in GetNumberOfSlowNucleons call.
+/* $Id$ */
 
-Revision 1.2  2003/03/25 09:55:24  morsch
-Numbers of slow nucleons either from model or user set.
-
-Revision 1.1  2003/03/24 16:49:23  morsch
-Slow nucleon generator and model.
-
-*/
-
-/*
-  Generator for slow nucluons in pA interactions. 
-  Source is modelled by a relativistic Maxwell distributions.
-  Original code by  Ferenc Sikler  <sikler@rmki.kfki.hu>
- */
+//
+//  Generator for slow nucleons in pA interactions. 
+//  Source is modelled by a relativistic Maxwell distributions.
+//  This class cooparates with AliCollisionGeometry if used inside AliGenCocktail.
+//  In this case the number of slow nucleons is determined from the number of wounded nuclei
+//  using a realisation of AliSlowNucleonModel.
+//  Original code by  Ferenc Sikler  <sikler@rmki.kfki.hu>
+// 
 
 #include <TDatabasePDG.h>
 #include <TPDGCode.h>
 #include <TH2F.h>
+#include <TH1F.h>
+#include <TF1.h>
 #include <TCanvas.h>
 
 #include "AliCollisionGeometry.h"
@@ -49,6 +42,8 @@ Slow nucleon generator and model.
 // Default constructor
     fSlowNucleonModel = 0;
     fCollisionGeometry = 0;
+    fDebugHist1 = 0;
+    fDebugHist2 = 0;
 }
 
 AliGenSlowNucleons::AliGenSlowNucleons(Int_t npart)
@@ -68,6 +63,14 @@ AliGenSlowNucleons::AliGenSlowNucleons(Int_t npart)
     fDebug = 0;
 }
 
+//____________________________________________________________
+AliGenSlowNucleons::AliGenSlowNucleons(const AliGenSlowNucleons & sn):
+    AliGenerator(sn)
+{
+// Copy constructor
+    sn.Copy(*this);
+}
+
 //____________________________________________________________
 AliGenSlowNucleons::~AliGenSlowNucleons()
 {
@@ -75,6 +78,10 @@ AliGenSlowNucleons::~AliGenSlowNucleons()
     delete  fSlowNucleonModel;
 }
 
+void AliGenSlowNucleons::SetProtonDirection(Float_t dir) {
+// Set direction of the proton to change between pA (1) and Ap (-1)
+  fProtonDirection = dir / TMath::Abs(dir);
+}
 
 void AliGenSlowNucleons::Init()
 {
@@ -87,11 +94,22 @@ void AliGenSlowNucleons::Init()
     if (fDebug) {
        fDebugHist1 = new TH2F("DebugHist1", "nu vs N_slow", 100, 0., 100., 20, 0., 20.);
        fDebugHist2 = new TH2F("DebugHist2", "b  vs N_slow", 100, 0., 100., 15, 0., 15.);
+       fCosThetaGrayHist = new TH1F("fCosThetaGrayHist", "Gray particles angle", 100, -1., 1.);
+    }
+    //
+    // non-uniform cos(theta) distribution
+    //
+    if(fThetaDistribution != 0) {
+       fCosTheta = new TF1("fCosTheta",
+                           "(2./3.14159265358979312)/(exp(2./3.14159265358979312)-exp(-2./3.14159265358979312))*exp(2.*x/3.14159265358979312)",
+                           -1., 1.);
     }
 }
 
 void AliGenSlowNucleons::FinishRun()
 {
+// End of run action
+// Show histogram for debugging if requested.
     if (fDebug) {
        TCanvas *c = new TCanvas("c","Canvas 1",400,10,600,700);
        c->Divide(2,1);
@@ -99,6 +117,8 @@ void AliGenSlowNucleons::FinishRun()
        fDebugHist1->Draw();
        c->cd(2);
        fDebugHist2->Draw();
+       c->cd(3);
+       fCosThetaGrayHist->Draw();
     }
 }
 
@@ -129,7 +149,7 @@ void AliGenSlowNucleons::Generate()
     }     
 
    //
-    Float_t p[3];
+    Float_t p[3], theta=0;
     Float_t origin[3] = {0., 0., 0.};
     Float_t polar [3] = {0., 0., 0.};    
     Int_t nt, i, j;
@@ -145,8 +165,9 @@ void AliGenSlowNucleons::Generate()
     fCharge = 1;
     kf = kProton;    
     for(i = 0; i < fNgp; i++) {
-       GenerateSlow(fCharge, fTemperatureG, fBetaSourceG, p);
-       SetTrack(fTrackIt, -1, kf, p, origin, polar,
+       GenerateSlow(fCharge, fTemperatureG, fBetaSourceG, p, theta);
+       if (fDebug) fCosThetaGrayHist->Fill(TMath::Cos(theta));
+       PushTrack(fTrackIt, -1, kf, p, origin, polar,
                 0., kPNoProcess, nt, 1.);
        KeepTrack(nt);
     }
@@ -156,8 +177,9 @@ void AliGenSlowNucleons::Generate()
     fCharge = 0;
     kf = kNeutron;    
     for(i = 0; i < fNgn; i++) {
-       GenerateSlow(fCharge, fTemperatureG, fBetaSourceG, p);
-       SetTrack(fTrackIt, -1, kf, p, origin, polar,
+       GenerateSlow(fCharge, fTemperatureG, fBetaSourceG, p, theta);
+       if (fDebug) fCosThetaGrayHist->Fill(TMath::Cos(theta));
+       PushTrack(fTrackIt, -1, kf, p, origin, polar,
                 0., kPNoProcess, nt, 1.);
        KeepTrack(nt);
     }
@@ -167,8 +189,8 @@ void AliGenSlowNucleons::Generate()
     fCharge = 1;
     kf = kProton;    
     for(i = 0; i < fNbp; i++) {
-       GenerateSlow(fCharge, fTemperatureB, fBetaSourceB, p);
-       SetTrack(fTrackIt, -1, kf, p, origin, polar,
+       GenerateSlow(fCharge, fTemperatureB, fBetaSourceB, p, theta);
+       PushTrack(fTrackIt, -1, kf, p, origin, polar,
                 0., kPNoProcess, nt, 1.);
        KeepTrack(nt);
     }
@@ -178,8 +200,8 @@ void AliGenSlowNucleons::Generate()
     fCharge = 0;
     kf = kNeutron;    
     for(i = 0; i < fNbn; i++) {
-       GenerateSlow(fCharge, fTemperatureB, fBetaSourceB, p);
-       SetTrack(fTrackIt, -1, kf, p, origin, polar,
+       GenerateSlow(fCharge, fTemperatureB, fBetaSourceB, p, theta);
+       PushTrack(fTrackIt, -1, kf, p, origin, polar,
                 0., kPNoProcess, nt, 1.);
        KeepTrack(nt);
     }
@@ -188,16 +210,17 @@ void AliGenSlowNucleons::Generate()
 
 
 
-void AliGenSlowNucleons::GenerateSlow(Int_t charge, Double_t T, Double_t beta, Float_t* q)
+void AliGenSlowNucleons::GenerateSlow(Int_t charge, Double_t T, 
+       Double_t beta, Float_t* q, Float_t &theta)
+
+{
 /* 
    Emit a slow nucleon with "temperature" T [GeV], 
    from a source moving with velocity beta         
    Three-momentum [GeV/c] is given back in q[3]    
 */
 
-{
- Double_t m, pmax, p, f, theta, phi;
+ Double_t m, pmax, p, f, phi;
  TDatabasePDG * pdg = TDatabasePDG::Instance();
  const Double_t kMassProton  = pdg->GetParticle(kProton) ->Mass();
  const Double_t kMassNeutron = pdg->GetParticle(kNeutron)->Mass();
@@ -221,8 +244,14 @@ void AliGenSlowNucleons::GenerateSlow(Int_t charge, Double_t T, Double_t beta, F
  }
  while(f < Rndm());
 
- /* Spherical symmetric emission */
- theta = TMath::ACos(2. * Rndm() - 1.);
+ /* Spherical symmetric emission for black particles (beta=0)*/
+ if(beta==0 || fThetaDistribution==0) theta = TMath::ACos(2. * Rndm() - 1.);
+ /* cos theta distributed according to experimental results for gray particles (beta=0.05)*/
+ else if(fThetaDistribution!=0){
+   Double_t costheta = fCosTheta->GetRandom();
+   theta = TMath::ACos(costheta);
+ }
+ //
  phi   = 2. * TMath::Pi() * Rndm();
 
  /* Determine momentum components in system of the moving source */
@@ -236,6 +265,7 @@ void AliGenSlowNucleons::GenerateSlow(Int_t charge, Double_t T, Double_t beta, F
 
  /* Transform to laboratory system */
  Lorentz(m, fBeta, q);
+ q[2] *= fProtonDirection; 
 }
 
 Double_t AliGenSlowNucleons::Maxwell(Double_t m, Double_t p, Double_t T)
@@ -256,6 +286,22 @@ void AliGenSlowNucleons::Lorentz(Double_t m, Double_t beta, Float_t* q)
     q[2] = gamma * (q[2] + beta*energy);
 }
 
+         
+AliGenSlowNucleons& AliGenSlowNucleons::operator=(const  AliGenSlowNucleons& rhs)
+{
+// Assignment operator
+    rhs.Copy(*this);
+    return *this;
+}
+
+void AliGenSlowNucleons::Copy(TObject&) const
+{
+    //
+    // Copy 
+    //
+    Fatal("Copy","Not implemented!\n");
+}
+