newly implemented model for slow nucleon production
[u/mrichter/AliRoot.git] / EVGEN / AliSlowNucleonModelExp.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 /* $Id$ */
17
18 //
19 // Experimental data inspired Gray Particle Model for p-Pb collisions
20 // The number of gray nucleons  is proportional to the number of collisions.
21 // The number of black nucleons is proportional to the number of collisions
22 // Fluctuations are calculated from a binomial distribution.
23 // Author: A.Morsch
24 //
25
26 #include "AliSlowNucleonModelExp.h"
27 #include "AliCollisionGeometry.h"
28 #include <TRandom.h>
29 #include <TMath.h>
30
31 ClassImp(AliSlowNucleonModelExp)
32
33
34 AliSlowNucleonModelExp::AliSlowNucleonModelExp():
35     fP(82),
36     fN (126),
37     fAlphaGray(2.3),
38     fAlphaBlack(3.6),
39     fApplySaturation(kTRUE),
40     fnGraySaturation(15),
41     fnBlackSaturation(28)
42 {
43   //
44   // Default constructor
45   //
46   //
47   printf("\n\nInitializing slow nucleon model with parameters:\n");
48   printf(" \t alpha_{gray} %1.2f  alpha_{black} %1.2f\n",fAlphaGray, fAlphaBlack);
49   //printf(" \t SATURATION %d w. %d (gray) %d (black) \n\n",fApplySaturation,fnGraySaturation,fnBlackSaturation);
50 }
51
52
53 void AliSlowNucleonModelExp::GetNumberOfSlowNucleons(AliCollisionGeometry* geo, 
54                                                       Int_t& ngp, Int_t& ngn, Int_t & nbp, Int_t & nbn) const
55 {
56 //
57 // Return the number of black and gray nucleons
58 //
59 // Number of collisions
60
61     Float_t nu = geo->NN() + geo->NwN() + geo->NNw(); 
62
63 // Mean number of gray nucleons 
64
65     Float_t nGray         = fAlphaGray * nu;
66     Float_t nGrayNeutrons = nGray * fN / (fN + fP);
67     Float_t nGrayProtons  = nGray - nGrayNeutrons;
68
69 // Mean number of black nucleons 
70     Float_t nBlack  = 0.;
71     if(!fApplySaturation || (fApplySaturation && nGray<fnGraySaturation)) nBlack = fAlphaBlack * nu;
72     else if(fApplySaturation && nGray>=fnGraySaturation) nBlack = fnBlackSaturation;
73     Float_t nBlackNeutrons = nBlack * 0.84;
74     Float_t nBlackProtons  = nBlack - nBlackNeutrons;
75
76 // Actual number (including fluctuations) from binomial distribution
77     Double_t p;
78
79 //  gray neutrons
80     p =  nGrayNeutrons/fN;
81     ngn = gRandom->Binomial((Int_t) fN, p);
82     
83 //  gray protons
84     p =  nGrayProtons/fP;
85     ngp = gRandom->Binomial((Int_t) fP, p);
86
87 //  black neutrons
88     p =  nBlackNeutrons/fN;
89     nbn = gRandom->Binomial((Int_t) fN, p);
90     
91 //  black protons
92     p =  nBlackProtons/fP;
93     nbp = gRandom->Binomial((Int_t) fP, p);
94
95 }
96
97 void AliSlowNucleonModelExp::GetNumberOfSlowNucleons2(AliCollisionGeometry* geo, 
98                                                       Int_t& ngp, Int_t& ngn, Int_t & nbp, Int_t & nbn) const
99 {
100 //
101 // Return the number of black and gray nucleons
102 //
103 // Number of collisions
104
105    // based on E910 model ================================================================
106
107    Float_t nu = (Float_t) (geo->NN() + geo->NwN() + geo->NNw()); 
108    //
109    nu = nu+1.*gRandom->Rndm();
110    //
111    Float_t  poverpd = 0.843; 
112    Float_t  zAu2zPb = 82./79.;
113    Float_t  nGrayp = (-0.27 + 0.63 * nu - 0.0008 *nu *nu)*poverpd*zAu2zPb;
114
115 //  gray protons
116     Double_t p;
117     p =  nGrayp/fP;
118     ngp = gRandom->Binomial((Int_t) fP, p);
119     //ngp = gRandom->Gaus(nGrayp, TMath::Sqrt(fP*p*(1-p)));
120     if(nGrayp<0.) ngp=0;
121     
122     //Float_t blackovergray = 3./7.;// from spallation
123     Float_t blackovergray = 0.65; // from COSY
124     Float_t nBlackp  = blackovergray*nGrayp; 
125
126 //  black protons
127     p =  nBlackp/fP;
128     nbp = gRandom->Binomial((Int_t) fP, p);
129     //nbp = gRandom->Gaus(nBlackp, TMath::Sqrt(fP*p*(1-p)));
130     if(nBlackp<0.) nbp=0;
131     
132     if(nu<3.){
133       nGrayp = -0.836 + 0.9112 *nu - 0.05381 *nu *nu;
134       nBlackp  = blackovergray*nGrayp; 
135     }
136     
137     Float_t nGrayNeutrons = 0.;
138     Float_t nBlackNeutrons = 0.;
139     Float_t cp = (nGrayp+nBlackp)/0.24;
140     
141     if(cp>0.){
142       Float_t nSlow  = 51.5+469.2/(-8.762-cp);
143       //if(cp<2.5) nSlow = 1+(9.9-1)/(2.5-0)*(cp-0);
144       if(cp<3.) nSlow = 0.+(11.6-0.)/(3.-0.)*(cp-0.);
145     
146       nGrayNeutrons = nSlow * 0.1; 
147       nBlackNeutrons = nSlow - nGrayNeutrons;
148     }
149     else{
150       // Sikler "pasturato"
151       nGrayNeutrons = 0.47 * fAlphaGray *  nu; 
152       nBlackNeutrons = 0.88 * fAlphaBlack * nu;      
153       printf("nslowp=0 -> ncoll = %1.0f -> ngrayn = %1.0f  nblackn = %1.0f \n", nu, nGrayNeutrons, nBlackNeutrons);
154     }
155     
156 //  gray neutrons
157     p =  nGrayNeutrons/fN;
158 //    ngn = gRandom->Binomial((Int_t) fN, p);
159     ngn = gRandom->Gaus(nGrayNeutrons, TMath::Sqrt(fN*p*(1-p)));
160
161 //  black neutrons
162     p =  nBlackNeutrons/fN;
163 //    nbn = gRandom->Binomial((Int_t) fN, p);
164     nbn = gRandom->Gaus(nBlackNeutrons, TMath::Sqrt(fN*p*(1-p)));
165     
166     
167 }
168
169 void AliSlowNucleonModelExp::SetParameters(Float_t alpha1, Float_t alpha2)
170 {
171     // Set the model parameters
172     fAlphaGray  = alpha1;
173     fAlphaBlack = alpha2;
174 }
175