]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TUHKMgen/UHKM/InitialStateHydjet.cxx
Classes for track-based multiplicity analysis in PbPb. Initial commit.
[u/mrichter/AliRoot.git] / TUHKMgen / UHKM / InitialStateHydjet.cxx
CommitLineData
3fa37a65 1//expanding localy equilibated fireball with volume hadron radiation
2//thermal part: Blast wave model, Bjorken-like parametrization
3//hyght-pt: PYTHIA + jet quenching model PYQUEN
03896fc4 4//
5// HYDJET++
6// version 1.0:
7// InitialStateHydjet is the modified InitialStateBjorken
8// The high-pt part related with PYTHIA-PYQUEN is included
9// InitialStateBjorken (FASTMC) was used.
10//
11//
12//
13// InitialStateBjorken
14// version 2.0:
15// Ludmila Malinina malinina@lav01.sinp.msu.ru, SINP MSU/Moscow and JINR/Dubna
16// Ionut Arsene i.c.arsene@fys.uio.no, Oslo University
17// June 2007
18//
19// version 1.0:
20// Nikolai Amelin, Ludmila Malinina, Timur Pocheptsov (C) JINR/Dubna
21// amelin@sunhe.jinr.ru, malinina@sunhe.jinr.ru, pocheptsov@sunhe.jinr.ru
22// November. 2, 2005
3fa37a65 23
03896fc4 24
b1c2e580 25
b1c2e580 26
03896fc4 27#include <iostream>
28#include <fstream>
b1c2e580 29
b1c2e580 30#include <TLorentzVector.h>
31#include <TVector3.h>
b1c2e580 32#include <TMath.h>
33
b1c2e580 34#include "InitialStateHydjet.h"
b1c2e580 35#include "RandArrayFunction.h"
b1c2e580 36#include "GrandCanonical.h"
b1c2e580 37#include "StrangePotential.h"
b1c2e580 38#include "Particle.h"
b1c2e580 39#include "ParticlePDG.h"
b1c2e580 40#include "HYJET_COMMONS.h"
03896fc4 41
b1c2e580 42extern "C" void hyevnt_();
43extern "C" void myini_();
b1c2e580 44extern HYIPARCommon HYIPAR;
45extern HYFPARCommon HYFPAR;
46extern HYJPARCommon HYJPAR;
47extern HYPARTCommon HYPART;
48extern SERVICECommon SERVICE;
49
50using std::cout;
51using std::endl;
52
3fa37a65 53class ParticleAllocator;
54class TRandom3;
55
56// declaration of the static member fLastIndex
03896fc4 57Int_t Particle::fgLastIndex;
3fa37a65 58
03896fc4 59//_________________________________________________________________________________
b1c2e580 60void InitialStateHydjet::Initialize(List_t &source, ParticleAllocator & allocator) {
7b7936e9 61 // Generate initial particles from the soft and hard components
3fa37a65 62
63 // Initialize the static "last index variable"
64 Particle::InitIndexing();
65
b1c2e580 66 //----- high-pt part------------------------------
b1c2e580 67 TLorentzVector partJMom, partJPos, zeroVec;
b1c2e580 68
7b7936e9 69 // run a HYDJET event
70 hyevnt_();
7b7936e9 71
b1c2e580 72 if(fParams.fNhsel != 0) {
3fa37a65 73 //get number of particles in jets
b1c2e580 74 Int_t numbJetPart = HYPART.njp;
3fa37a65 75
76 for(Int_t i = 0; i <numbJetPart; i++) {
77 Int_t pdg = Int_t(HYPART.ppart[i][1]);
b1c2e580 78 Double_t px = HYPART.ppart[i][2];
79 Double_t py = HYPART.ppart[i][3];
80 Double_t pz = HYPART.ppart[i][4];
81 Double_t e = HYPART.ppart[i][5];
82 Double_t vx = HYPART.ppart[i][6];
83 Double_t vy = HYPART.ppart[i][7];
84 Double_t vz = HYPART.ppart[i][8];
85 Double_t vt = HYPART.ppart[i][9];
b1c2e580 86 ParticlePDG *partDef = fDatabase->GetPDGParticle(pdg);
87 Int_t type =1; //from jet
3fa37a65 88 if(partDef) {
89 partJMom.SetXYZT(px, py, pz, e);
90 partJPos.SetXYZT(vx, vy, vz, vt);
91 Particle *particle=new Particle(partDef, partJPos, partJMom, 0, 0, type, -1, zeroVec, zeroVec);
92 particle->SetIndex();
93 allocator.AddParticle(*particle, source);
94 delete particle;
95 }
b1c2e580 96 }
7b7936e9 97 } //nhsel !=0 not only hydro!
b1c2e580 98
b1c2e580 99
100 //----------HYDRO part------------------------------------------------
101 if(fParams.fNhsel < 3) {
102 const Double_t weightMax = 2*TMath::CosH(fParams.fUmax);
103 const Int_t nBins = 100;
104 Double_t probList[nBins];
105 RandArrayFunction arrayFunctDistE(nBins);
106 RandArrayFunction arrayFunctDistR(nBins);
107
108 TLorentzVector partPos, partMom, n1, p0;
109 TVector3 vec3;
b1c2e580 110 //set maximal hadron energy
111 const Double_t eMax = 5.;
112 //-------------------------------------
113 // get impact parameter
b1c2e580 114
115 //effective volume for central
116 double dYl= 2 * fParams.fYlmax; //uniform distr. [-Ylmax; Ylmax]
117 if(fParams.fEtaType >0) dYl = TMath::Sqrt(2 * TMath::Pi()) * fParams.fYlmax ; //Gaussian distr.
3eeebf49 118 Double_t volEffcent = 2 * TMath::Pi() * fParams.fTau * dYl *
119 (fParams.fR * fParams.fR)/TMath::Power((fParams.fUmax),2)*
120 ((fParams.fUmax)*TMath::SinH((fParams.fUmax))-TMath::CosH((fParams.fUmax))+ 1);
b1c2e580 121
122 //effective volume for non-central Simpson2
7b7936e9 123 Double_t volEffnoncent = fParams.fTau * dYl * SimpsonIntegrator2(0., 2.*TMath::Pi());
124 fVolEff = volEffcent * HYFPAR.npart/HYFPAR.npart0;
b1c2e580 125
7b7936e9 126 Double_t coeffRB = TMath::Sqrt(volEffcent * HYFPAR.npart/HYFPAR.npart0/volEffnoncent);
127 Double_t coeffR1 = HYFPAR.npart/HYFPAR.npart0;
128 coeffR1 = TMath::Power(coeffR1, 0.333333);
b1c2e580 129
7b7936e9 130 double veff=fVolEff;
b1c2e580 131
b1c2e580 132 //------------------------------------
133 //cycle on particles types
134 for(Int_t i = 0; i < fParams.fNPartTypes; ++i) {
7b7936e9 135 Double_t mparam = fParams.fPartMult[2 * i] * veff;
136 Int_t multiplicity = gRandom->Poisson(mparam);
b1c2e580 137 const Int_t encoding = fParams.fPartEnc[i];
138
139 if(multiplicity > 0) {
140 ParticlePDG *partDef = fDatabase->GetPDGParticle(encoding);
141 if(!partDef) {
142 Error("InitialStateHydjet::Initialize", "No particle with encoding %d", encoding);
143 continue;
144 }
145 //no charm now !
03896fc4 146 if(TMath::Abs(partDef->GetCharmQNumber())>0 || TMath::Abs(partDef->GetCharmAQNumber())>0){
b1c2e580 147 continue;
148 }
149
150 //compute chemical potential for single f.o. mu==mu_ch
151 //compute chemical potential for thermal f.o.
152 Double_t mu = fParams.fPartMu[2 * i];
153
154 //choose Bose-Einstein or Fermi-Dirac statistics
155 const Double_t d = !(Int_t(2*partDef->GetSpin()) & 1) ? -1 : 1;
156 const Double_t mass = partDef->GetMass();
157
158 //prepare histogram to sample hadron energy:
159 Double_t h = (eMax - mass) / nBins;
160 Double_t x = mass + 0.5 * h;
7b7936e9 161 Int_t ii;
162 for(ii = 0; ii < nBins; ++ii) {
163 if(x>=mu && fParams.fThFO>0)probList[ii] = x * TMath::Sqrt(x * x - mass * mass) /
b1c2e580 164 (TMath::Exp((x - mu) / (fParams.fThFO)) + d);
7b7936e9 165 if(x>=mu && fParams.fThFO<=0)probList[ii] = x * TMath::Sqrt(x * x - mass * mass) /
b1c2e580 166 (TMath::Exp((x - mu) / (fParams.fT)) + d);
7b7936e9 167 if(x<mu)probList[ii] = 0.;
b1c2e580 168 x += h;
169 }
170 arrayFunctDistE.PrepareTable(probList);
171
172 //prepare histogram to sample hadron transverse radius:
173 h = (fParams.fR) / nBins;
174 x = 0.5 * h;
175 Double_t param = (fParams.fUmax) / (fParams.fR);
7b7936e9 176 for(ii = 0; ii < nBins; ++ii) {
177 probList[ii] = x * TMath::CosH(param*x);
b1c2e580 178 x += h;
179 }
180 arrayFunctDistR.PrepareTable(probList);
181
182 //loop over hadrons, assign hadron coordinates and momenta
183 Double_t weight = 0., yy = 0., px0 = 0., py0 = 0., pz0 = 0.;
184 Double_t e = 0., x0 = 0., y0 = 0., z0 = 0., t0 = 0., etaF = 0.;
7b7936e9 185 Double_t r, rB, phiF;
b1c2e580 186
187 for(Int_t j = 0; j < multiplicity; ++j) {
188 do {
189 fParams.fEtaType <=0 ? etaF = fParams.fYlmax * (2. * gRandom->Rndm() - 1.)
190 : etaF = (fParams.fYlmax) * (gRandom->Gaus());
191 n1.SetXYZT(0.,0.,TMath::SinH(etaF),TMath::CosH(etaF));
192 if(TMath::Abs(etaF)>5.)continue;
193
7b7936e9 194 rB = fParams.fR * coeffRB * coeffR1;
b1c2e580 195
196 Double_t rho = TMath::Sqrt(gRandom->Rndm());
197 Double_t phi = TMath::TwoPi() * gRandom->Rndm();
7b7936e9 198 Double_t rx = TMath::Sqrt(1-fParams.fEpsilon)*rB;
199 Double_t ry = TMath::Sqrt(1+fParams.fEpsilon)*rB;
b1c2e580 200
7b7936e9 201 x0 = rx * rho * TMath::Cos(phi);
202 y0 = ry * rho * TMath::Sin(phi);
b1c2e580 203 r = TMath::Sqrt(x0*x0+y0*y0);
204 phiF = TMath::Abs(TMath::ATan(y0/x0));
205
206 if(x0<0&&y0>0)phiF = TMath::Pi()-phiF;
207 if(x0<0&&y0<0)phiF = TMath::Pi()+phiF;
208 if(x0>0&&y0<0)phiF = 2.*TMath::Pi()-phiF;
209
210 //proper time with emission duration
7b7936e9 211 Double_t tau = coeffR1 * fParams.fTau + sqrt(2.) * fParams.fSigmaTau * coeffR1 * (gRandom->Gaus());
b1c2e580 212 z0 = tau * TMath::SinH(etaF);
213 t0 = tau * TMath::CosH(etaF);
214
7b7936e9 215 Double_t rhou = fParams.fUmax * r / rB;
b1c2e580 216
b1c2e580 217
218 Double_t uxf = TMath::SinH(rhou)*TMath::Sqrt(1+fParams.fDelta)*TMath::Cos(phiF);
219 Double_t uyf = TMath::SinH(rhou)*TMath::Sqrt(1-fParams.fDelta)*TMath::Sin(phiF);
220 Double_t utf = TMath::CosH(etaF) * TMath::CosH(rhou) *
221 TMath::Sqrt(1+fParams.fDelta*TMath::Cos(2*phiF)*TMath::TanH(rhou)*TMath::TanH(rhou));
222 Double_t uzf = TMath::SinH(etaF) * TMath::CosH(rhou) *
223 TMath::Sqrt(1+fParams.fDelta*TMath::Cos(2*phiF)*TMath::TanH(rhou)*TMath::TanH(rhou));
224
225 vec3.SetXYZ(uxf / utf, uyf / utf, uzf / utf);
226 n1.Boost(-vec3);
227
228 yy = weightMax * gRandom->Rndm();
229
230 Double_t php0 = TMath::TwoPi() * gRandom->Rndm();
231 Double_t ctp0 = 2. * gRandom->Rndm() - 1.;
bfd20868 232 Double_t stp0 = TMath::Sqrt((1.-ctp0)*(1.+ctp0));
b1c2e580 233 e = mass + (eMax - mass) * arrayFunctDistE();
bfd20868 234 Double_t pp0 = TMath::Sqrt((e-mass)*(e+mass));
b1c2e580 235 px0 = pp0 * stp0 * TMath::Sin(php0);
236 py0 = pp0 * stp0 * TMath::Cos(php0);
237 pz0 = pp0 * ctp0;
238 p0.SetXYZT(px0, py0, pz0, e);
239
240 //weight for rdr
241 weight = (n1 * p0) /e; // weight for rdr gammar: weight = (n1 * p0) / n1[3] / e;
242 } while(yy >= weight);
243
b1c2e580 244 partMom.SetXYZT(px0, py0, pz0, e);
245 partPos.SetXYZT(x0, y0, z0, t0);
246 partMom.Boost(vec3);
247 Int_t type =0; //hydro
7b7936e9 248
3fa37a65 249 Particle *particle=new Particle(partDef, partPos, partMom, 0., 0, type, -1, zeroVec, zeroVec);
250 particle->SetIndex();
251 allocator.AddParticle(*particle, source);
252 delete particle;
b1c2e580 253 } //nhsel==4 , no hydro part
254 }
255 }
256 }
257
b1c2e580 258}
259
03896fc4 260//_________________________________________________________________________________
b1c2e580 261Bool_t InitialStateHydjet::ReadParams() {
7b7936e9 262 // Read parameters from an input file in ascii
263
b1c2e580 264 Float_t par[200] = {0.};
265 Int_t i = 0;
266 std::string s(40,' ');
267 std::ifstream input("RunInputHydjet");
268 if (!input) {
269 Error("Ukm::ReadParams", "Cannot open RunInputHydjet");
270 return kFALSE;
271 }
272
273 while (std::getline(input, s)) {
274 input>>par[i];
275 if (i < 140)
276 std::cout<<s<<" = "<<par[i]<<std::endl;
277 ++i;
278 std::getline(input,s);
279 }
280
281 std::cout<<"\nFor output use the files RunOutput.root \n\n"<< std::endl;
282
283 fParams.fNevnt = Int_t(par[0]); //number of events
284 fParams.fSqrtS = par[1]; //cms energy per nucleon
285 fParams.fAw = par[2]; // atomic number of colliding nuclei
286 fParams.fIfb = Int_t(par[3]); // flag of type of centrality generation (=0 is fixed by fBfix, not 0
287 //impact parameter is generated in each event between fBfmin
288 //and fBmax according with Glauber model (f-la 30)
289 fParams.fBmin = par[4]; //minimum impact parameter in units of nuclear radius RA
290 fParams.fBmax = par[5]; //maximum impact parameter in units of nuclear radius RA
291 fParams.fBfix = par[6]; //fix impact parameter in units of nuclear radius RA
292
293 fParams.fSeed = Int_t(par[7]); //parameter to set the random nuber seed (=0 the current time is used
294 //to set the random generator seed, !=0 the value fSeed is
295 //used to set the random generator seed and then the state of random
296 //number generator in PYTHIA MRPY(1)=fSeed
297
298 fParams.fT = par[8]; //chemical freeze-out temperature in GeV
299 fParams.fMuB = par[9]; //baryon potential
300 fParams.fMuS = par[10]; //strangeness potential
301 fParams.fMuI3 = par[11]; //isospin potential
302 fParams.fThFO = par[12]; //thermal freeze-out temperature T^th in GeV
303 fParams.fMu_th_pip = par[13]; // effective chemical potential of positivly charged pions at thermal in GeV
304
305
306 fParams.fTau = par[14]; //proper time value
307 fParams.fSigmaTau = par[15]; //its standart deviation (emission duration)
308 fParams.fR = par[16]; //maximal transverse radius
309 fParams.fYlmax = par[17]; //maximal longitudinal rapidity
310 fParams.fUmax = par[18]; //maximal transverse velocity multiplaed on \gamma_r
311 fParams.fDelta = par[19]; //momentum asymmetry parameter
312 fParams.fEpsilon = par[20]; //coordinate asymmetry parameter
313
314 fParams.fDecay = Int_t(par[21]); // flag to switch on/off hadron decays<0: decays off,>=0: decays on, (default: 0)
315 fParams.fWeakDecay = Int_t(par[22]); //flag to switch on/off weak hadron decays <0: decays off, >0: decays on, (default: 0)
316
317 fParams.fEtaType = Int_t(par[23]); // flag to choose rapidity distribution, if fEtaType<=0,
318 //then uniform rapidity distribution in [-fYlmax,fYlmax] if fEtaType>0,
319 //then Gaussian with dispertion = fYlmax
320
321 fParams.fTMuType = Int_t(par[24]); // flag to use calculated chemical freeze-out temperature,
322 //baryon potential and strangeness potential as a function of fSqrtS
323
324 fParams.fCorrS = par[25]; // flag and value to include strangeness supression factor
325 fParams.fNhsel = Int_t(par[26]); //flag to switch on/off jet and hydro-state production (0: jet
326 // production off and hydro on, 1: jet production on and jet quenching
327 // off and hydro on, 2: jet production on and jet quenching on and
328 // hydro on, 3: jet production on and jet quenching off and hydro
329 // off, 4: jet production on and jet quenching on and hydro off
330
331 fParams.fIshad= Int_t(par[27]); //flag to switch on/off impact parameter dependent nuclear
332 // shadowing for gluons and light sea quarks (u,d,s) (0: shadowing off,
333 // 1: shadowing on for fAw=207, 197, 110, 40, default: 1
334
335 fParams.fPtmin = par[28]; //minimal transverse momentum transfer p_T of hard
336 // parton-parton scatterings in GeV (the PYTHIA parameter ckin(3)=fPtmin)
337
338 // PYQUEN energy loss model parameters:
339
340 fParams.fT0 = par[29]; // initial temperature (in GeV) of QGP for
341 //central Pb+Pb collisions at mid-rapidity (initial temperature for other
342 //centralities and atomic numbers will be calculated automatically) (allowed range is 0.2<fT0<2)
343
344 fParams.fTau0= par[30]; //proper QGP formation time in fm/c (0.01<fTau0<10)
345 fParams.fNf= Int_t(par[31]); //number of active quark flavours N_f in QGP fNf=0, 1,2 or 3
346 fParams.fIenglu= Int_t(par[32]); // flag to fix type of in-medium partonic energy loss
347 //(0: radiative and collisional loss, 1: radiative loss only, 2:
348 //collisional loss only) (default: 0);
349 fParams.fIanglu= Int_t(par[33]); //flag to fix type of angular distribution of in-medium emitted
350 // gluons (0: small-angular, 1: wide-angular, 2:collinear) (default: 0).
351
352
353 //PYTHIA parameters:
354 Int_t jj;
355 for (Int_t j = 0; j <25; ++j) {
356 jj= 35+j;
357 SERVICE.parPYTH[j]=par[jj];
358 }
359
360 // Set Random Number seed
361
362 gRandom->SetSeed(fParams.fSeed); //Set 0 to use the current time
363 //to send seed in PYTHIA
364 SERVICE.iseed_fromC=gRandom->GetSeed();
365 std::cout<<"Seed for random number generation= "<<gRandom->GetSeed()<<std::endl;
366
367 fParams.fNPartTypes = 0; //counter of hadron species
368 return kTRUE;
369}
370
03896fc4 371//_________________________________________________________________________________
b1c2e580 372Bool_t InitialStateHydjet::MultIni() {
7b7936e9 373 // Calculate average multiplicities, chemical potentials (if necessary),
374 // initialize pyquen
375
b1c2e580 376 //check and redefine input parameters
377 if(fParams.fTMuType>0 && fParams.fSqrtS > 2.24) {
378 if(fParams.fSqrtS < 2.24){
379 Error("InitialStateHydjet::MultIni", "SqrtS<2.24 not allowed with fParams.fTMuType>0");
380 return 0;
381 }
382
383 //sqrt(s) = 2.24 ==> T_kin = 0.8 GeV
384 //see J. Cleymans, H. Oeschler, K. Redlich,S. Wheaton, Phys Rev. C73 034905 (2006)
385 fParams.fMuB = 1.308/(1. + fParams.fSqrtS*0.273);
386 fParams.fT = 0.166 - 0.139*fParams.fMuB*fParams.fMuB - 0.053*fParams.fMuB*fParams.fMuB*
387 fParams.fMuB*fParams.fMuB;
388 fParams.fMuI3 = 0.;
389 fParams.fMuS = 0.;
390 //create strange potential object and set strangeness density 0
03896fc4 391 StrangePotential* psp = new StrangePotential(0., fDatabase);
b1c2e580 392 psp->SetBaryonPotential(fParams.fMuB);
393 psp->SetTemperature(fParams.fT);
394 //compute strangeness potential
395 if(fParams.fMuB > 0.01)
396 fParams.fMuS = psp->CalculateStrangePotential();
b1c2e580 397
398 //if user choose fYlmax larger then allowed by kinematics at the specified beam energy sqrt(s)
399 if(fParams.fYlmax > TMath::Log(fParams.fSqrtS/0.94)){
400 Error("InitialStateHydjet::MultIni", "fParams.fYlmax > TMath::Log(fParams.fSqrtS/0.94)!!! ");
401 return 0;
402 }
403
404
405 if(fParams.fCorrS <= 0.) {
406 //see F. Becattini, J. Mannien, M. Gazdzicki, Phys Rev. C73 044905 (2006)
407 fParams.fCorrS = 1. - 0.386* TMath::Exp(-1.23*fParams.fT/fParams.fMuB);
b1c2e580 408
409 }
410 std::cout<<"The phenomenological f-la J. Cleymans et al. PRC73 034905 (2006) for Tch mu_B was used." << std::endl;
411 std::cout<<"The simulation will be done with the calculated parameters:" << std::endl;
412 std::cout<<"Baryon chemical potential = "<<fParams.fMuB<< " [GeV]" << std::endl;
413 std::cout<<"Strangeness chemical potential = "<<fParams.fMuS<< " [GeV]" << std::endl;
414 std::cout<<"Isospin chemical potential = "<<fParams.fMuI3<< " [GeV]" << std::endl;
415 std::cout<<"Strangeness suppression parameter = "<<fParams.fCorrS << std::endl;
416 std::cout<<"Eta_max = "<<fParams.fYlmax<< std::endl;
417 std::cout << std::endl;
418
419 }
420
421
b1c2e580 422 //initialisation of high-pt part
423
424 HYJPAR.nhsel = fParams.fNhsel;
425 HYJPAR.ptmin = fParams.fPtmin;
426 HYJPAR.ishad = fParams.fIshad;
427 HYIPAR.bminh = fParams.fBmin;
428 HYIPAR.bmaxh = fParams.fBmax;
429 HYIPAR.AW = fParams.fAw;
430
431 HYPYIN.ifb = fParams.fIfb;
432 HYPYIN.bfix = fParams.fBfix;
433 HYPYIN.ene = fParams.fSqrtS;
434
435 PYQPAR.T0 = fParams.fT0;
436 PYQPAR.tau0 = fParams.fTau0;
437 PYQPAR.nf = fParams.fNf;
438 PYQPAR.ienglu = fParams.fIenglu;
439 PYQPAR.ianglu = fParams.fIanglu;
440
441
7b7936e9 442 myini_(); //
b1c2e580 443
444
445 // calculation of multiplicities of different particle species
446 // according to the grand canonical approach
447 GrandCanonical gc(15, fParams.fT, fParams.fMuB, fParams.fMuS, fParams.fMuI3);
3fa37a65 448 GrandCanonical gcCh(15, fParams.fT, fParams.fMuB, fParams.fMuS, fParams.fMuI3);
449 GrandCanonical gcPiTh(15, fParams.fThFO, 0., 0., fParams.fMu_th_pip);
450 GrandCanonical gcTh0(15, fParams.fThFO, 0., 0., 0.);
b1c2e580 451
452 //effective volume for central
453 double dYl= 2 * fParams.fYlmax; //uniform distr. [-Ylmax; Ylmax]
454 if (fParams.fEtaType >0) dYl = TMath::Sqrt(2 * TMath::Pi()) * fParams.fYlmax ; //Gaussian distr.
455 fVolEff = 2 * TMath::Pi() * fParams.fTau * dYl * (fParams.fR * fParams.fR)/TMath::Power((fParams.fUmax),2) *
456 ((fParams.fUmax)*TMath::SinH((fParams.fUmax))-TMath::CosH((fParams.fUmax))+ 1);
b1c2e580 457
3fa37a65 458 Double_t particleDensityPiCh=0;
459 Double_t particleDensityPiTh=0;
b1c2e580 460
461 if(fParams.fThFO != fParams.fT && fParams.fThFO > 0){
3fa37a65 462 particleDensityPiCh = gcCh.ParticleNumberDensity(fDatabase->GetPDGParticle(211));
463 particleDensityPiTh = gcPiTh.ParticleNumberDensity(fDatabase->GetPDGParticle(211));
b1c2e580 464 }
465
466 for(Int_t particleIndex = 0; particleIndex < fDatabase->GetNParticles(); particleIndex++) {
467 ParticlePDG *currParticle = fDatabase->GetPDGParticleByIndex(particleIndex);
468 Int_t encoding = currParticle->GetPDG();
469 //strangeness supression
470 Double_t gammaS = 1;
7b7936e9 471 Int_t s = Int_t(currParticle->GetStrangeness());
472 if(encoding == 333)
473 s = 2;
474 if(fParams.fCorrS < 1. && s != 0)
475 gammaS = TMath::Power(fParams.fCorrS,-TMath::Abs(s));
b1c2e580 476 //average densities
7b7936e9 477 Double_t particleDensity = gc.ParticleNumberDensity(currParticle)/gammaS;
b1c2e580 478
479 //compute chemical potential for single f.o. mu==mu_ch
480 Double_t mu = fParams.fMuB * Int_t(currParticle->GetBaryonNumber()) +
481 fParams.fMuS * Int_t(currParticle->GetStrangeness()) +
482 fParams.fMuI3 * Int_t(currParticle->GetElectricCharge());
483
484 //thermal f.o.
485 if(fParams.fThFO != fParams.fT && fParams.fThFO > 0){
3fa37a65 486 Double_t particleDensityCh = gcCh.ParticleNumberDensity(currParticle);
487 Double_t particleDensityTh0 = gcTh0.ParticleNumberDensity(currParticle);
488 Double_t numbDensBolt = particleDensityPiTh*particleDensityCh/particleDensityPiCh;
489 mu = fParams.fThFO*TMath::Log(numbDensBolt/particleDensityTh0);
b1c2e580 490 if(abs(encoding)==211 || encoding==111)mu= fParams.fMu_th_pip;
3fa37a65 491 particleDensity = numbDensBolt;
b1c2e580 492 }
493
7b7936e9 494 // set particle number density to zero for some species
495 // photons
496 if(abs(encoding)==22)
497 particleDensity=0;
498 // K0L and K0S
499 if(abs(encoding)==130 || abs(encoding)==310) {
500 particleDensity=0;
501 }
b1c2e580 502
503 if(particleDensity > 0.) {
504 // outMult<<encoding<< " " <<particleDensity<< " "<<mu<<std::endl;
505 fParams.fPartEnc[fParams.fNPartTypes] = encoding;
506 fParams.fPartMult[2 * fParams.fNPartTypes] = particleDensity;
507 fParams.fPartMu[2 * fParams.fNPartTypes] = mu;
508 ++fParams.fNPartTypes;
509 if(fParams.fNPartTypes > 1000)
510 Error("in Bool_t MultIni:", "fNPartTypes is too large %d", fParams.fNPartTypes);
511 }
512 }
513 return kTRUE;
514}
515
03896fc4 516//_________________________________________________________________________________
b1c2e580 517Double_t InitialStateHydjet::SimpsonIntegrator2(Double_t a, Double_t b) {
7b7936e9 518 // Simpson integration
b1c2e580 519 Int_t nsubIntervals=10000;
520 Double_t h = (b - a)/nsubIntervals; //0-pi, phi
521 Double_t s=0;
b1c2e580 522 Double_t x = 0; //phi
523 for(Int_t j = 1; j < nsubIntervals; j++) {
524 x += h; // phi
525 Double_t e = fParams.fEpsilon;
3fa37a65 526 Double_t rSB = fParams.fR; //test: podstavit' *coefff_RB
527 Double_t rB = rSB *(TMath::Sqrt(1-e*e)/TMath::Sqrt(1+e*TMath::Cos(2*x))); //f-la7 rB
7b7936e9 528 Double_t sr = SimpsonIntegrator(0,rB,x);
b1c2e580 529 s += sr;
530 }
531 return s*h;
532
533}
534
03896fc4 535//_________________________________________________________________________________
b1c2e580 536Double_t InitialStateHydjet::SimpsonIntegrator(Double_t a, Double_t b, Double_t phi) {
7b7936e9 537 // Simpson integration
b1c2e580 538 Int_t nsubIntervals=100;
539 Double_t h = (b - a)/nsubIntervals;
03896fc4 540 Double_t s = F2(phi,a + 0.5*h);
541 Double_t t = 0.5*(F2(phi,a) + F2(phi,b));
b1c2e580 542 Double_t x = a;
543 Double_t y = a + 0.5*h;
544 for(Int_t i = 1; i < nsubIntervals; i++) {
545 x += h;
546 y += h;
03896fc4 547 s += F2(phi,y);
548 t += F2(phi,x);
b1c2e580 549 }
550 t += 2.0*s;
551 return t*h/3.0;
552}
553
554
555//f2=f(phi,r)
03896fc4 556//_________________________________________________________________________________
557Double_t InitialStateHydjet::F2(Double_t x, Double_t y) {
7b7936e9 558 // formula
3fa37a65 559 Double_t rSB = fParams.fR; //test: podstavit' *coefff_RB
560 Double_t rhou = fParams.fUmax * y / rSB;
b1c2e580 561 Double_t ff = y*TMath::CosH(rhou)*
562 TMath::Sqrt(1+fParams.fDelta*TMath::Cos(2*x)*TMath::TanH(rhou)*TMath::TanH(rhou));
563//n_mu u^mu f-la 20
564 return ff;
565}
566
03896fc4 567//_________________________________________________________________________________
b1c2e580 568Double_t InitialStateHydjet::MidpointIntegrator2(Double_t a, Double_t b) {
7b7936e9 569 // Perform integration through the mid-point method
b1c2e580 570 Int_t nsubIntervals=2000;
571 Int_t nsubIntervals2=1;
572 Double_t h = (b - a)/nsubIntervals; //0-pi , phi
7b7936e9 573 Double_t h2 = (fParams.fR)/nsubIntervals; //0-R maximal rB ?
b1c2e580 574 Double_t x = a + 0.5*h;
575 Double_t y = 0;
03896fc4 576 Double_t t = F2(x,y);
b1c2e580 577 Double_t e = fParams.fEpsilon;
b1c2e580 578 for(Int_t j = 1; j < nsubIntervals; j++) {
579 x += h; // integr phi
3fa37a65 580 Double_t rSB = fParams.fR; //test: podstavit' *coefff_RB
581 Double_t rB = rSB *(TMath::Sqrt(1-e*e)/TMath::Sqrt(1+e*TMath::Cos(2*x))); //f-la7 rB
7b7936e9 582 nsubIntervals2 = Int_t(rB / h2)+1;
b1c2e580 583 // integr R
584 y=0;
585 for(Int_t i = 1; i < nsubIntervals2; i++)
03896fc4 586 t += F2(x,(y += h2));
b1c2e580 587 }
588 return t*h*h2;
589}
590