]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVGEN/AliSlowNucleonModelExp.cxx
increase cell time diff histogram binning
[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     fLCPparam(0.585)
43 {
44   //
45   // Default constructor
46   //
47   //
48   fSlownparam[0] = 60.;
49   fSlownparam[1] = 469.2;
50   fSlownparam[2] = 8.762;
51   printf("\n\n ******** Initializing slow nucleon model with parameters:\n");
52   printf(" \t alpha_{gray} %1.2f  alpha_{black} %1.2f\n",fAlphaGray, fAlphaBlack);
53   printf(" \t SATURATION %d w. %d (gray) %d (black) \n\n",fApplySaturation,fnGraySaturation,fnBlackSaturation);
54   printf(" \t LCP parameter %f   Slown parameters = {%f, %f, %f}\n\n",fLCPparam,fSlownparam[0],fSlownparam[1],fSlownparam[2]); 
55 }
56
57
58 void AliSlowNucleonModelExp::GetNumberOfSlowNucleons(AliCollisionGeometry* geo, 
59                                                       Int_t& ngp, Int_t& ngn, Int_t & nbp, Int_t & nbn) const
60 {
61 //
62 // Return the number of black and gray nucleons
63 //
64 // Number of collisions
65
66     Float_t nu = geo->NN() + geo->NwN() + geo->NNw(); 
67
68 // Mean number of gray nucleons 
69
70     Float_t nGray         = fAlphaGray * nu;
71     Float_t nGrayNeutrons = nGray * fN / (fN + fP);
72     Float_t nGrayProtons  = nGray - nGrayNeutrons;
73
74 // Mean number of black nucleons 
75     Float_t nBlack  = 0.;
76     if(!fApplySaturation || (fApplySaturation && nGray<fnGraySaturation)) nBlack = fAlphaBlack * nu;
77     else if(fApplySaturation && nGray>=fnGraySaturation) nBlack = fnBlackSaturation;
78     Float_t nBlackNeutrons = nBlack * 0.84;
79     Float_t nBlackProtons  = nBlack - nBlackNeutrons;
80
81 // Actual number (including fluctuations) from binomial distribution
82     Double_t p;
83
84 //  gray neutrons
85     p =  nGrayNeutrons/fN;
86     ngn = gRandom->Binomial((Int_t) fN, p);
87     
88 //  gray protons
89     p =  nGrayProtons/fP;
90     ngp = gRandom->Binomial((Int_t) fP, p);
91
92 //  black neutrons
93     p =  nBlackNeutrons/fN;
94     nbn = gRandom->Binomial((Int_t) fN, p);
95     
96 //  black protons
97     p =  nBlackProtons/fP;
98     nbp = gRandom->Binomial((Int_t) fP, p);
99
100 }
101
102 void AliSlowNucleonModelExp::GetNumberOfSlowNucleons2(AliCollisionGeometry* geo, 
103                                                       Int_t& ngp, Int_t& ngn, Int_t & nbp, Int_t & nbn) const
104 {
105 //
106 // Return the number of black and gray nucleons
107 //
108 // Number of collisions
109
110    // based on E910 model ================================================================
111
112    Float_t nu = (Float_t) (geo->NN() + geo->NwN() + geo->NNw()); 
113    //
114    nu = nu+1.*gRandom->Rndm();
115    //
116    Float_t  poverpd = 0.843; 
117    Float_t  zAu2zPb = 82./79.;
118    Float_t  nGrayp = (-0.27 + 0.63 * nu - 0.0008 *nu *nu)*poverpd*zAu2zPb;
119
120 //  gray protons
121     Double_t p;
122     p =  nGrayp/fP;
123     ngp = gRandom->Binomial((Int_t) fP, p);
124     //ngp = gRandom->Gaus(nGrayp, TMath::Sqrt(fP*p*(1-p)));
125     if(nGrayp<0.) ngp=0;
126     
127     //Float_t blackovergray = 3./7.;// from spallation
128     Float_t blackovergray = 0.65; // from COSY
129     Float_t nBlackp  = blackovergray*nGrayp; 
130
131 //  black protons
132     p =  nBlackp/fP;
133     nbp = gRandom->Binomial((Int_t) fP, p);
134     //nbp = gRandom->Gaus(nBlackp, TMath::Sqrt(fP*p*(1-p)));
135     if(nBlackp<0.) nbp=0;
136     
137     if(nu<3.){
138       nGrayp = -0.836 + 0.9112 *nu - 0.05381 *nu *nu;
139       nBlackp  = blackovergray*nGrayp; 
140     }
141     
142     //printf(" \t Using LCP parameter %f   Slown parameters = {%f, %f, %f}\n\n",fLCPparam,fSlownparam[0],fSlownparam[1],fSlownparam[2]); 
143     Float_t nGrayNeutrons = 0.;
144     Float_t nBlackNeutrons = 0.;
145     Float_t cp = (nGrayp+nBlackp)/fLCPparam;
146     
147     if(cp>0.){
148       Float_t nSlow      = fSlownparam[0]+fSlownparam[1]/(-fSlownparam[2]-cp);
149       Float_t paramRetta = fSlownparam[0]+fSlownparam[1]/(-fSlownparam[2]-3);
150       if(cp<3.) nSlow = 0.+(paramRetta-0.)/(3.-0.)*(cp-0.);
151     
152       nGrayNeutrons = nSlow * 0.1; 
153       nBlackNeutrons = nSlow - nGrayNeutrons;
154     }
155     else{
156       // Sikler "pasturato" (qui non entra mai!!!!)
157       nGrayNeutrons = 0.47 * fAlphaGray *  nu; 
158       nBlackNeutrons = 0.88 * fAlphaBlack * nu;      
159       printf("nslowp=0 -> ncoll = %1.0f -> ngrayn = %1.0f  nblackn = %1.0f \n", nu, nGrayNeutrons, nBlackNeutrons);
160     }
161     
162 //  gray neutrons
163     p =  nGrayNeutrons/fN;
164 //    ngn = gRandom->Binomial((Int_t) fN, p);
165     ngn = gRandom->Gaus(nGrayNeutrons, TMath::Sqrt(fN*p*(1-p)));
166
167 //  black neutrons
168     p =  nBlackNeutrons/fN;
169 //    nbn = gRandom->Binomial((Int_t) fN, p);
170     nbn = gRandom->Gaus(nBlackNeutrons, TMath::Sqrt(fN*p*(1-p)));
171     
172     
173 }
174
175 void AliSlowNucleonModelExp::SetParameters(Float_t alpha1, Float_t alpha2)
176 {
177     // Set the model parameters
178     fAlphaGray  = alpha1;
179     fAlphaBlack = alpha2;
180 }
181