ROOT logo
/**************************************************************************
 * 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.                  *
 **************************************************************************/

/* $Id$ */

//
// Experimental data inspired Gray Particle Model for p-Pb collisions
// The number of gray nucleons  is proportional to the number of collisions.
// The number of black nucleons is proportional to the number of collisions
// Fluctuations are calculated from a binomial distribution.
// Author: A.Morsch
//

#include "AliSlowNucleonModelExp.h"
#include "AliCollisionGeometry.h"
#include <TRandom.h>
#include <TMath.h>

ClassImp(AliSlowNucleonModelExp)


AliSlowNucleonModelExp::AliSlowNucleonModelExp():
    fP(82),
    fN (126),
    fAlphaGray(2.3),
    fAlphaBlack(3.6),
    fApplySaturation(kTRUE),
    fnGraySaturation(15),
    fnBlackSaturation(28),
    fLCPparam(0.585),
    fSigmaSmear(0.25)
{
  //
  // Default constructor
  //
  //
  fSlownparam[0] = 60.;
  fSlownparam[1] = 469.2;
  fSlownparam[2] = 8.762;
  /*printf("\n\n ******** Initializing slow nucleon model with parameters:\n");
  printf(" \t alpha_{gray} %1.2f  alpha_{black} %1.2f\n",fAlphaGray, fAlphaBlack);
  printf(" \t SATURATION %d w. %d (gray) %d (black) \n\n",fApplySaturation,fnGraySaturation,fnBlackSaturation);
  printf(" \t LCP parameter %f   Slown parameters = {%f, %f,
  %f}\n\n",fLCPparam,fSlownparam[0],fSlownparam[1],fSlownparam[2]); */
}


void AliSlowNucleonModelExp::GetNumberOfSlowNucleons(AliCollisionGeometry* geo, 
						      Int_t& ngp, Int_t& ngn, Int_t & nbp, Int_t & nbn) const
{
//
// Return the number of black and gray nucleons
//
// Number of collisions

    Float_t nu = geo->NN() + geo->NwN() + geo->NNw(); 

// Mean number of gray nucleons 

    Float_t nGray         = fAlphaGray * nu;
    Float_t nGrayNeutrons = nGray * fN / (fN + fP);
    Float_t nGrayProtons  = nGray - nGrayNeutrons;

// Mean number of black nucleons 
    Float_t nBlack  = 0.;
    if(!fApplySaturation || (fApplySaturation && nGray<fnGraySaturation)) nBlack = fAlphaBlack * nu;
    else if(fApplySaturation && nGray>=fnGraySaturation) nBlack = fnBlackSaturation;
    Float_t nBlackNeutrons = nBlack * 0.84;
    Float_t nBlackProtons  = nBlack - nBlackNeutrons;

// Actual number (including fluctuations) from binomial distribution
    Double_t p;

//  gray neutrons
    p =  nGrayNeutrons/fN;
    ngn = gRandom->Binomial((Int_t) fN, p);
    
//  gray protons
    p =  nGrayProtons/fP;
    ngp = gRandom->Binomial((Int_t) fP, p);

//  black neutrons
    p =  nBlackNeutrons/fN;
    nbn = gRandom->Binomial((Int_t) fN, p);
    
//  black protons
    p =  nBlackProtons/fP;
    nbp = gRandom->Binomial((Int_t) fP, p);

}

void AliSlowNucleonModelExp::GetNumberOfSlowNucleons2(AliCollisionGeometry* geo, 
						      Int_t& ngp, Int_t& ngn, Int_t & nbp, Int_t & nbn) const
{
//
// Return the number of black and gray nucleons
//
// Number of collisions

   // based on E910 model ================================================================

   Float_t nu = (Float_t) (geo->NN() + geo->NwN() + geo->NNw()); 
   //
   //nu = nu+1.*gRandom->Rndm();
   nu = gRandom->Gaus(nu, 0.5);
   if(nu<0.) nu=0.;
   //
   Float_t  poverpd = 0.843; 
   Float_t  zAu2zPb = 82./79.;
   Float_t  nGrayp = (-0.27 + 0.63 * nu - 0.0008 *nu *nu)*poverpd*zAu2zPb;

//  gray protons
    Double_t p;
    p =  nGrayp/fP;
    ngp = gRandom->Binomial((Int_t) fP, p);
    //ngp = gRandom->Gaus(nGrayp, TMath::Sqrt(fP*p*(1-p)));
    if(nGrayp<0.) ngp=0;
    
    //Float_t blackovergray = 3./7.;// from spallation
    Float_t blackovergray = 0.65; // from COSY
    Float_t nBlackp  = blackovergray*nGrayp; 

//  black protons
    p =  nBlackp/fP;
    nbp = gRandom->Binomial((Int_t) fP, p);
    //nbp = gRandom->Gaus(nBlackp, TMath::Sqrt(fP*p*(1-p)));
    if(nBlackp<0.) nbp=0;
    
    if(nu<3.){
      nGrayp = -0.836 + 0.9112 *nu - 0.05381 *nu *nu;
      nBlackp  = blackovergray*nGrayp; 
    }
    
    //printf(" \t Using LCP parameter %f   Slown parameters = {%f, %f, %f}\n\n",fLCPparam,fSlownparam[0],fSlownparam[1],fSlownparam[2]); 
    Float_t nGrayNeutrons = 0.;
    Float_t nBlackNeutrons = 0.;
    Float_t cp = (nGrayp+nBlackp)/fLCPparam;
    
    if(cp>0.){
      Float_t nSlow      = fSlownparam[0]+fSlownparam[1]/(-fSlownparam[2]-cp);
      Float_t paramRetta = fSlownparam[0]+fSlownparam[1]/(-fSlownparam[2]-3);
      if(cp<3.) nSlow = 0.+(paramRetta-0.)/(3.-0.)*(cp-0.);
    
      nGrayNeutrons = nSlow * 0.1; 
      nBlackNeutrons = nSlow - nGrayNeutrons;
    }
    else{
      // Sikler "pasturato" (qui non entra mai!!!!)
      nGrayNeutrons = 0.47 * fAlphaGray *  nu; 
      nBlackNeutrons = 0.88 * fAlphaBlack * nu;      
      //printf("nslowp=0 -> ncoll = %1.0f -> ngrayn = %1.0f  nblackn = %1.0f \n", nu, nGrayNeutrons, nBlackNeutrons);
    }
    
//  gray neutrons
    p =  nGrayNeutrons/fN;
//    ngn = gRandom->Binomial((Int_t) fN, p);
    ngn = gRandom->Gaus(nGrayNeutrons, TMath::Sqrt(fN*p*(1-p)));

//  black neutrons
    p =  nBlackNeutrons/fN;
//    nbn = gRandom->Binomial((Int_t) fN, p);
    nbn = gRandom->Gaus(nBlackNeutrons, TMath::Sqrt(fN*p*(1-p)));
    
    
}

void AliSlowNucleonModelExp::GetNumberOfSlowNucleons2s(AliCollisionGeometry* geo, 
						      Int_t& ngp, Int_t& ngn, Int_t & nbp, Int_t & nbn) const
{
//
// Return the number of black and gray nucleons
//
// Number of collisions

   // based on E910 model ================================================================

   Float_t nu = (Float_t) (geo->NN() + geo->NwN() + geo->NNw()); 
   //
   Float_t  poverpd = 0.843; 
   Float_t  zAu2zPb = 82./79.;
   Float_t  grayp = (-0.27 + 0.63 * nu - 0.0008 *nu *nu)*poverpd*zAu2zPb;
   Float_t  nGrayp = gRandom->Gaus(grayp, fSigmaSmear);
   if(nGrayp<0.) nGrayp=0.;

//  gray protons
    Double_t p=0.;
    p = nGrayp/fP;
    ngp = gRandom->Binomial((Int_t) fP, p);
    //ngp = gRandom->Gaus(nGrayp, TMath::Sqrt(fP*p*(1-p)));
    if(nGrayp<0.) ngp=0;
    
    //Float_t blackovergray = 3./7.;// from spallation
    Float_t blackovergray = 0.65; // from COSY
    //Float_t blackp  = blackovergray*grayp; 
    //Float_t nBlackp = gRandom->Gaus(nblackp, fSigmaSmear);
    Float_t nBlackp = blackovergray*nGrayp;
    if(nBlackp<0.) nBlackp=0.;

//  black protons
    p =  nBlackp/fP;
    nbp = gRandom->Binomial((Int_t) fP, p);
    //nbp = gRandom->Gaus(nBlackp, TMath::Sqrt(fP*p*(1-p)));
    if(nBlackp<0.) nbp=0;
    
    Float_t nGrayNeutrons = 0.;
    Float_t nBlackNeutrons = 0.;
    Float_t cp = (nGrayp+nBlackp)/fLCPparam;
    
    if(cp>0.){
      Float_t nSlow = fSlownparam[0]+fSlownparam[1]/(-fSlownparam[2]-cp);
      
      nGrayNeutrons = nSlow * 0.1; 
      nBlackNeutrons = nSlow - nGrayNeutrons;
    }
    else{
      // Sikler "pasturato" (qui non entra mai!!!!)
      nGrayNeutrons = 0.47 * fAlphaGray *  nu; 
      nBlackNeutrons = 0.88 * fAlphaBlack * nu;      
      //printf("nslowp=0 -> ncoll = %1.0f -> ngrayn = %1.0f  nblackn = %1.0f \n", nu, nGrayNeutrons, nBlackNeutrons);
    }
    //
    if(nGrayNeutrons<0.) nGrayNeutrons=0.;
    if(nBlackNeutrons<0.) nBlackNeutrons=0.;
    
//  gray neutrons
    p =  nGrayNeutrons/fN;
//    ngn = gRandom->Binomial((Int_t) fN, p);
    ngn = gRandom->Gaus(nGrayNeutrons, TMath::Sqrt(fN*p*(1-p)));
    if(nGrayNeutrons<0.) ngn=0;

//  black neutrons
    p =  nBlackNeutrons/fN;
//    nbn = gRandom->Binomial((Int_t) fN, p);
    nbn = gRandom->Gaus(nBlackNeutrons, TMath::Sqrt(fN*p*(1-p)));
    if(nBlackNeutrons<0.) nbn=0;
    
}

void AliSlowNucleonModelExp::SetParameters(Float_t alpha1, Float_t alpha2)
{
    // Set the model parameters
    fAlphaGray  = alpha1;
    fAlphaBlack = alpha2;
}

 AliSlowNucleonModelExp.cxx:1
 AliSlowNucleonModelExp.cxx:2
 AliSlowNucleonModelExp.cxx:3
 AliSlowNucleonModelExp.cxx:4
 AliSlowNucleonModelExp.cxx:5
 AliSlowNucleonModelExp.cxx:6
 AliSlowNucleonModelExp.cxx:7
 AliSlowNucleonModelExp.cxx:8
 AliSlowNucleonModelExp.cxx:9
 AliSlowNucleonModelExp.cxx:10
 AliSlowNucleonModelExp.cxx:11
 AliSlowNucleonModelExp.cxx:12
 AliSlowNucleonModelExp.cxx:13
 AliSlowNucleonModelExp.cxx:14
 AliSlowNucleonModelExp.cxx:15
 AliSlowNucleonModelExp.cxx:16
 AliSlowNucleonModelExp.cxx:17
 AliSlowNucleonModelExp.cxx:18
 AliSlowNucleonModelExp.cxx:19
 AliSlowNucleonModelExp.cxx:20
 AliSlowNucleonModelExp.cxx:21
 AliSlowNucleonModelExp.cxx:22
 AliSlowNucleonModelExp.cxx:23
 AliSlowNucleonModelExp.cxx:24
 AliSlowNucleonModelExp.cxx:25
 AliSlowNucleonModelExp.cxx:26
 AliSlowNucleonModelExp.cxx:27
 AliSlowNucleonModelExp.cxx:28
 AliSlowNucleonModelExp.cxx:29
 AliSlowNucleonModelExp.cxx:30
 AliSlowNucleonModelExp.cxx:31
 AliSlowNucleonModelExp.cxx:32
 AliSlowNucleonModelExp.cxx:33
 AliSlowNucleonModelExp.cxx:34
 AliSlowNucleonModelExp.cxx:35
 AliSlowNucleonModelExp.cxx:36
 AliSlowNucleonModelExp.cxx:37
 AliSlowNucleonModelExp.cxx:38
 AliSlowNucleonModelExp.cxx:39
 AliSlowNucleonModelExp.cxx:40
 AliSlowNucleonModelExp.cxx:41
 AliSlowNucleonModelExp.cxx:42
 AliSlowNucleonModelExp.cxx:43
 AliSlowNucleonModelExp.cxx:44
 AliSlowNucleonModelExp.cxx:45
 AliSlowNucleonModelExp.cxx:46
 AliSlowNucleonModelExp.cxx:47
 AliSlowNucleonModelExp.cxx:48
 AliSlowNucleonModelExp.cxx:49
 AliSlowNucleonModelExp.cxx:50
 AliSlowNucleonModelExp.cxx:51
 AliSlowNucleonModelExp.cxx:52
 AliSlowNucleonModelExp.cxx:53
 AliSlowNucleonModelExp.cxx:54
 AliSlowNucleonModelExp.cxx:55
 AliSlowNucleonModelExp.cxx:56
 AliSlowNucleonModelExp.cxx:57
 AliSlowNucleonModelExp.cxx:58
 AliSlowNucleonModelExp.cxx:59
 AliSlowNucleonModelExp.cxx:60
 AliSlowNucleonModelExp.cxx:61
 AliSlowNucleonModelExp.cxx:62
 AliSlowNucleonModelExp.cxx:63
 AliSlowNucleonModelExp.cxx:64
 AliSlowNucleonModelExp.cxx:65
 AliSlowNucleonModelExp.cxx:66
 AliSlowNucleonModelExp.cxx:67
 AliSlowNucleonModelExp.cxx:68
 AliSlowNucleonModelExp.cxx:69
 AliSlowNucleonModelExp.cxx:70
 AliSlowNucleonModelExp.cxx:71
 AliSlowNucleonModelExp.cxx:72
 AliSlowNucleonModelExp.cxx:73
 AliSlowNucleonModelExp.cxx:74
 AliSlowNucleonModelExp.cxx:75
 AliSlowNucleonModelExp.cxx:76
 AliSlowNucleonModelExp.cxx:77
 AliSlowNucleonModelExp.cxx:78
 AliSlowNucleonModelExp.cxx:79
 AliSlowNucleonModelExp.cxx:80
 AliSlowNucleonModelExp.cxx:81
 AliSlowNucleonModelExp.cxx:82
 AliSlowNucleonModelExp.cxx:83
 AliSlowNucleonModelExp.cxx:84
 AliSlowNucleonModelExp.cxx:85
 AliSlowNucleonModelExp.cxx:86
 AliSlowNucleonModelExp.cxx:87
 AliSlowNucleonModelExp.cxx:88
 AliSlowNucleonModelExp.cxx:89
 AliSlowNucleonModelExp.cxx:90
 AliSlowNucleonModelExp.cxx:91
 AliSlowNucleonModelExp.cxx:92
 AliSlowNucleonModelExp.cxx:93
 AliSlowNucleonModelExp.cxx:94
 AliSlowNucleonModelExp.cxx:95
 AliSlowNucleonModelExp.cxx:96
 AliSlowNucleonModelExp.cxx:97
 AliSlowNucleonModelExp.cxx:98
 AliSlowNucleonModelExp.cxx:99
 AliSlowNucleonModelExp.cxx:100
 AliSlowNucleonModelExp.cxx:101
 AliSlowNucleonModelExp.cxx:102
 AliSlowNucleonModelExp.cxx:103
 AliSlowNucleonModelExp.cxx:104
 AliSlowNucleonModelExp.cxx:105
 AliSlowNucleonModelExp.cxx:106
 AliSlowNucleonModelExp.cxx:107
 AliSlowNucleonModelExp.cxx:108
 AliSlowNucleonModelExp.cxx:109
 AliSlowNucleonModelExp.cxx:110
 AliSlowNucleonModelExp.cxx:111
 AliSlowNucleonModelExp.cxx:112
 AliSlowNucleonModelExp.cxx:113
 AliSlowNucleonModelExp.cxx:114
 AliSlowNucleonModelExp.cxx:115
 AliSlowNucleonModelExp.cxx:116
 AliSlowNucleonModelExp.cxx:117
 AliSlowNucleonModelExp.cxx:118
 AliSlowNucleonModelExp.cxx:119
 AliSlowNucleonModelExp.cxx:120
 AliSlowNucleonModelExp.cxx:121
 AliSlowNucleonModelExp.cxx:122
 AliSlowNucleonModelExp.cxx:123
 AliSlowNucleonModelExp.cxx:124
 AliSlowNucleonModelExp.cxx:125
 AliSlowNucleonModelExp.cxx:126
 AliSlowNucleonModelExp.cxx:127
 AliSlowNucleonModelExp.cxx:128
 AliSlowNucleonModelExp.cxx:129
 AliSlowNucleonModelExp.cxx:130
 AliSlowNucleonModelExp.cxx:131
 AliSlowNucleonModelExp.cxx:132
 AliSlowNucleonModelExp.cxx:133
 AliSlowNucleonModelExp.cxx:134
 AliSlowNucleonModelExp.cxx:135
 AliSlowNucleonModelExp.cxx:136
 AliSlowNucleonModelExp.cxx:137
 AliSlowNucleonModelExp.cxx:138
 AliSlowNucleonModelExp.cxx:139
 AliSlowNucleonModelExp.cxx:140
 AliSlowNucleonModelExp.cxx:141
 AliSlowNucleonModelExp.cxx:142
 AliSlowNucleonModelExp.cxx:143
 AliSlowNucleonModelExp.cxx:144
 AliSlowNucleonModelExp.cxx:145
 AliSlowNucleonModelExp.cxx:146
 AliSlowNucleonModelExp.cxx:147
 AliSlowNucleonModelExp.cxx:148
 AliSlowNucleonModelExp.cxx:149
 AliSlowNucleonModelExp.cxx:150
 AliSlowNucleonModelExp.cxx:151
 AliSlowNucleonModelExp.cxx:152
 AliSlowNucleonModelExp.cxx:153
 AliSlowNucleonModelExp.cxx:154
 AliSlowNucleonModelExp.cxx:155
 AliSlowNucleonModelExp.cxx:156
 AliSlowNucleonModelExp.cxx:157
 AliSlowNucleonModelExp.cxx:158
 AliSlowNucleonModelExp.cxx:159
 AliSlowNucleonModelExp.cxx:160
 AliSlowNucleonModelExp.cxx:161
 AliSlowNucleonModelExp.cxx:162
 AliSlowNucleonModelExp.cxx:163
 AliSlowNucleonModelExp.cxx:164
 AliSlowNucleonModelExp.cxx:165
 AliSlowNucleonModelExp.cxx:166
 AliSlowNucleonModelExp.cxx:167
 AliSlowNucleonModelExp.cxx:168
 AliSlowNucleonModelExp.cxx:169
 AliSlowNucleonModelExp.cxx:170
 AliSlowNucleonModelExp.cxx:171
 AliSlowNucleonModelExp.cxx:172
 AliSlowNucleonModelExp.cxx:173
 AliSlowNucleonModelExp.cxx:174
 AliSlowNucleonModelExp.cxx:175
 AliSlowNucleonModelExp.cxx:176
 AliSlowNucleonModelExp.cxx:177
 AliSlowNucleonModelExp.cxx:178
 AliSlowNucleonModelExp.cxx:179
 AliSlowNucleonModelExp.cxx:180
 AliSlowNucleonModelExp.cxx:181
 AliSlowNucleonModelExp.cxx:182
 AliSlowNucleonModelExp.cxx:183
 AliSlowNucleonModelExp.cxx:184
 AliSlowNucleonModelExp.cxx:185
 AliSlowNucleonModelExp.cxx:186
 AliSlowNucleonModelExp.cxx:187
 AliSlowNucleonModelExp.cxx:188
 AliSlowNucleonModelExp.cxx:189
 AliSlowNucleonModelExp.cxx:190
 AliSlowNucleonModelExp.cxx:191
 AliSlowNucleonModelExp.cxx:192
 AliSlowNucleonModelExp.cxx:193
 AliSlowNucleonModelExp.cxx:194
 AliSlowNucleonModelExp.cxx:195
 AliSlowNucleonModelExp.cxx:196
 AliSlowNucleonModelExp.cxx:197
 AliSlowNucleonModelExp.cxx:198
 AliSlowNucleonModelExp.cxx:199
 AliSlowNucleonModelExp.cxx:200
 AliSlowNucleonModelExp.cxx:201
 AliSlowNucleonModelExp.cxx:202
 AliSlowNucleonModelExp.cxx:203
 AliSlowNucleonModelExp.cxx:204
 AliSlowNucleonModelExp.cxx:205
 AliSlowNucleonModelExp.cxx:206
 AliSlowNucleonModelExp.cxx:207
 AliSlowNucleonModelExp.cxx:208
 AliSlowNucleonModelExp.cxx:209
 AliSlowNucleonModelExp.cxx:210
 AliSlowNucleonModelExp.cxx:211
 AliSlowNucleonModelExp.cxx:212
 AliSlowNucleonModelExp.cxx:213
 AliSlowNucleonModelExp.cxx:214
 AliSlowNucleonModelExp.cxx:215
 AliSlowNucleonModelExp.cxx:216
 AliSlowNucleonModelExp.cxx:217
 AliSlowNucleonModelExp.cxx:218
 AliSlowNucleonModelExp.cxx:219
 AliSlowNucleonModelExp.cxx:220
 AliSlowNucleonModelExp.cxx:221
 AliSlowNucleonModelExp.cxx:222
 AliSlowNucleonModelExp.cxx:223
 AliSlowNucleonModelExp.cxx:224
 AliSlowNucleonModelExp.cxx:225
 AliSlowNucleonModelExp.cxx:226
 AliSlowNucleonModelExp.cxx:227
 AliSlowNucleonModelExp.cxx:228
 AliSlowNucleonModelExp.cxx:229
 AliSlowNucleonModelExp.cxx:230
 AliSlowNucleonModelExp.cxx:231
 AliSlowNucleonModelExp.cxx:232
 AliSlowNucleonModelExp.cxx:233
 AliSlowNucleonModelExp.cxx:234
 AliSlowNucleonModelExp.cxx:235
 AliSlowNucleonModelExp.cxx:236
 AliSlowNucleonModelExp.cxx:237
 AliSlowNucleonModelExp.cxx:238
 AliSlowNucleonModelExp.cxx:239
 AliSlowNucleonModelExp.cxx:240
 AliSlowNucleonModelExp.cxx:241
 AliSlowNucleonModelExp.cxx:242
 AliSlowNucleonModelExp.cxx:243
 AliSlowNucleonModelExp.cxx:244
 AliSlowNucleonModelExp.cxx:245
 AliSlowNucleonModelExp.cxx:246
 AliSlowNucleonModelExp.cxx:247
 AliSlowNucleonModelExp.cxx:248
 AliSlowNucleonModelExp.cxx:249
 AliSlowNucleonModelExp.cxx:250
 AliSlowNucleonModelExp.cxx:251
 AliSlowNucleonModelExp.cxx:252
 AliSlowNucleonModelExp.cxx:253
 AliSlowNucleonModelExp.cxx:254
 AliSlowNucleonModelExp.cxx:255
 AliSlowNucleonModelExp.cxx:256
 AliSlowNucleonModelExp.cxx:257