]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Gray particle generator, first commit.
authormorsch <morsch@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 8 Oct 2002 13:53:17 +0000 (13:53 +0000)
committermorsch <morsch@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 8 Oct 2002 13:53:17 +0000 (13:53 +0000)
EVGEN/AliGenGrayParticles.cxx [new file with mode: 0644]
EVGEN/AliGenGrayParticles.h [new file with mode: 0644]
EVGEN/EVGENLinkDef.h
EVGEN/libEVGEN.pkg

diff --git a/EVGEN/AliGenGrayParticles.cxx b/EVGEN/AliGenGrayParticles.cxx
new file mode 100644 (file)
index 0000000..39d904a
--- /dev/null
@@ -0,0 +1,162 @@
+/**************************************************************************
+ * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ *                                                                        *
+ * Author: The ALICE Off-line Project.                                    *
+ * Contributors are mentioned in the code where appropriate.              *
+ *                                                                        *
+ * Permission to use, copy, modify and distribute this software and its   *
+ * documentation strictly for non-commercial purposes is hereby granted   *
+ * without fee, provided that the above copyright notice appears in all   *
+ * copies and that both the copyright notice and this permission notice   *
+ * appear in the supporting documentation. The authors make no claims     *
+ * about the suitability of this software for any purpose. It is          *
+ * provided "as is" without express or implied warranty.                  *
+ **************************************************************************/
+
+/*
+$Log$
+*/
+
+/*
+  Generator for gray nucluons in pA interactions. 
+  Source is modelled by a relativistic Maxwell distributions.
+  Original code by  Ferenc Sikler  <sikler@rmki.kfki.hu>
+ */
+#include "AliGenGrayParticles.h"
+#include "AliPDG.h"
+#include <TDatabasePDG.h>
+
+ ClassImp(AliGenGrayParticles)
+    
+ AliGenGrayParticles::AliGenGrayParticles():AliGenerator(-1)
+{
+// Default constructor
+}
+
+AliGenGrayParticles::AliGenGrayParticles(Int_t npart)
+    :AliGenerator(npart)
+{
+// Constructor
+    fName  = "GrayParticles";
+    fTitle = "Generator for gray particles in pA collisions";
+    SetPmax();
+    SetTarget();
+    SetNominalCmsEnergy();
+    SetTemperature();
+}
+
+//____________________________________________________________
+AliGenGrayParticles::~AliGenGrayParticles()
+{
+// Destructor
+}
+
+
+void AliGenGrayParticles::Init()
+{
+  //
+  // Initialization
+  //
+    Float_t kMass  = TDatabasePDG::Instance()->GetParticle(kProton)->Mass();
+    fMomentum = fCMS/2. * fZTarget / fATarget;
+    fBeta     = fMomentum / TMath::Sqrt(kMass * kMass + fMomentum * fMomentum);
+}
+
+
+void AliGenGrayParticles::Generate()
+{
+  //
+  // Generate one event
+  //
+    Float_t p[3];
+    Float_t origin[3] = {0., 0., 0.};
+    Float_t polar [3] = {0., 0., 0.};    
+    Int_t nt, i;
+    for(i = 0;i < fNpart; i++) {
+       Int_t kf = kProton;
+       GenerateSlow(1, fTemperature, 0., p);
+       
+       SetTrack(fTrackIt, -1, kf, p, origin, polar,
+                0., kPNoProcess, nt, 1.);
+
+       KeepTrack(nt);
+    }
+}
+
+
+
+
+void AliGenGrayParticles::GenerateSlow(Int_t charge, Double_t T, Double_t beta, Float_t* q)
+/* 
+   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;
+ TDatabasePDG * pdg = TDatabasePDG::Instance();
+ const Double_t kMassProton  = pdg->GetParticle(kProton) ->Mass();
+ const Double_t kMassNeutron = pdg->GetParticle(kNeutron)->Mass();
+ /* Select nucleon type */
+ if(charge == 0) m = kMassNeutron;
+ else m = kMassProton;
+
+ /* Momentum at maximum of Maxwell-distribution */
+
+ pmax = TMath::Sqrt(2*T*(T+sqrt(T*T+m*m)));
+
+ /* Try until proper momentum                                  */
+ /* for lack of primitive function of the Maxwell-distribution */
+ /* a brute force trial-accept loop, normalized at pmax        */
+
+ do
+ {
+     p = Rndm() * fPmax;
+     f = Maxwell(m, p, T) / Maxwell(m , pmax, T);
+ }
+ while(f < Rndm());
+
+ /* Spherical symmetric emission */
+ theta = TMath::ACos(2. * Rndm() - 1.);
+ phi   = 2. * TMath::Pi() * Rndm();
+
+ /* Determine momentum components in system of the moving source */
+ q[0] = p * TMath::Sin(theta) * TMath::Cos(phi);
+ q[1] = p * TMath::Sin(theta) * TMath::Sin(phi);
+ q[2] = p * TMath::Cos(theta);
+
+ /* Transform to system of the target nucleus                             */
+ /* beta is passed as negative, because the gray nucleons are slowed down */
+ Lorentz(m, -beta, q);
+
+ /* Transform to laboratory system */
+ Lorentz(m, fBeta, q);
+}
+
+Double_t AliGenGrayParticles::Maxwell(Double_t m, Double_t p, Double_t T)
+{
+/* Relativistic Maxwell-distribution */
+    Double_t ekin;
+    ekin = TMath::Sqrt(p*p+m*m)-m;
+    return (p*p * exp(-ekin/T));
+}
+
+
+void AliGenGrayParticles::Lorentz(Double_t m, Double_t beta, Float_t* q)
+{
+/* Lorentz transform in the direction of q[2] */
+    Double_t gamma  = 1/sqrt(1-beta*beta); 
+    Double_t energy = sqrt(m*m + q[0]*q[0] + q[1]*q[1] + q[2]*q[2]);
+    q[2] = gamma * (q[2] + beta*energy);
+}
+
+
+
+
+
+
+
diff --git a/EVGEN/AliGenGrayParticles.h b/EVGEN/AliGenGrayParticles.h
new file mode 100644 (file)
index 0000000..6dcec20
--- /dev/null
@@ -0,0 +1,43 @@
+#ifndef ALIGENGRAYPARTICLES_H
+#define ALIGENGRAYPARTICLES_H
+/* Copyright(c) 198-1999, ALICE Experiment at CERN, All rights reserved. *
+ * See cxx source for full Copyright notice                               */
+
+/* $Id$ */
+
+#include "AliGenerator.h"
+
+class AliGenGrayParticles : public AliGenerator
+{
+public:
+    AliGenGrayParticles();
+    AliGenGrayParticles(Int_t npart);
+    virtual ~AliGenGrayParticles();
+    virtual void Init();
+    virtual void Generate();
+    virtual void SetPmax(Float_t pmax = 10.) {fPmax = pmax;}
+    virtual void SetNominalCmsEnergy(Float_t energy = 14000.) {fCMS = energy;}
+    virtual void SetTarget(Float_t a=208, Float_t z=82) {fATarget = a; fZTarget = z;}
+    virtual void SetTemperature(Double_t t = 0.05) {fTemperature = t;}
+    
+ protected:
+    void     GenerateSlow(Int_t charge, Double_t T, Double_t beta, Float_t* q);
+    Double_t Maxwell(Double_t m, Double_t p, Double_t t);
+    void     Lorentz(Double_t m, Double_t beta, Float_t* q);
+ protected:
+    Float_t  fCMS;         // Center of mass energy
+    Float_t  fMomentum;    // Target nucleus momentum
+    Float_t  fBeta;        // Target nucleus beta
+    Float_t  fPmax;        // Maximum slow nucleon momentum
+    Float_t  fATarget;     // Target nucleus mass number
+    Float_t  fZTarget;     // Target nucleus charge number
+    Float_t  fTemperature; // Source Temperature
+  ClassDef(AliGenGrayParticles,1) // Gray Particle Generator
+};
+#endif
+
+
+
+
+
+
index ce13aef9f5d7e89bc17d2478cec3bf8770c44dcc..705846eddb223d77f86e5119e80f735bf45c5e39 100644 (file)
@@ -58,6 +58,7 @@
 #pragma link C++ class  AliGenAfterBurnerFlow+;
 #pragma link C++ class  AliPartonicEnergyLoss+;
 #pragma link C++ class  AliStructFuncType+;
+#pragma link C++ class  AliGenGrayParticles+;
 #endif
 
 
index 2f4a0f152b4420d512f8e4dd2be725e95d79a65f..53d10ee4fa3917796836e7f011860148e33ae71b 100644 (file)
@@ -20,7 +20,7 @@ SRCS          = AliGenHIJINGpara.cxx AliGenBox.cxx AliGenFixed.cxx \
                AliGenThetaSlice.cxx AliGenSTRANGElib.cxx AliGenBeamGas.cxx\
                AliGenAfterBurnerFlow.cxx \
                AliPartonicEnergyLoss.cxx AliGenHerwig.cxx \
-               AliStructFuncType.cxx
+               AliStructFuncType.cxx AliGenGrayParticles.cxx
 
 # Headerfiles for this particular package (Path respect to own directory)
 HDRS= $(SRCS:.cxx=.h)