]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TDPMjet/AliGenDPMjet.cxx
Use random number service class AliDpmJetRndm.
[u/mrichter/AliRoot.git] / TDPMjet / AliGenDPMjet.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
17 // Generator using DPMJET as an external generator
18 // The main DPMJET options are accessable for the user through this interface.
19 // Uses the TDPMjet implementation of TGenerator.
20
21 #include <TDPMjet.h>
22 #include <TRandom.h>
23 #include <TArrayI.h>
24 #include <TParticle.h>
25 #include <TGraph.h>
26 #include <TDatabasePDG.h>
27 #include <TParticlePDG.h>
28 #include <TParticleClassPDG.h>
29 #include <TPDGCode.h>
30 #include <TLorentzVector.h>
31
32 #include "AliGenDPMjet.h"
33 #include "AliGenDPMjetEventHeader.h"
34 #include "AliPythia.h"
35 #include "AliRun.h"
36 #include "AliDpmJetRndm.h"
37
38  ClassImp(AliGenDPMjet)
39
40
41 //______________________________________________________________________________
42 AliGenDPMjet::AliGenDPMjet()
43                  :AliGenMC()
44 {
45 // Constructor
46     fParticles = 0;
47     fDPMjet = 0;
48     AliDpmJetRndm::SetDpmJetRandom(GetRandom());
49 }
50
51
52 //______________________________________________________________________________
53 AliGenDPMjet::AliGenDPMjet(Int_t npart)
54     :AliGenMC(npart)
55 {
56 // Default PbPb collisions at 5. 5 TeV
57 //
58     fName = "DPMJET";
59     fTitle= "Particle Generator using DPMJET";
60
61     SetBeamEnergy();
62     SetEnergyCMS();
63     SetTarget();
64     SetProjectile();
65     SetCentral();
66     SetImpactParameterRange();
67     SetBoostLHC();
68     
69     fKeep       =  0;
70     fDecaysOff  =  1;
71     fEvaluate   =  0;
72     fSelectAll  =  0;
73     fFlavor     =  0;
74     fSpectators =  1;
75     fVertex.Set(3);
76         
77     fParticles = new TClonesArray("TParticle",10000);    
78
79     // Set random number generator   
80     fDPMjet = 0;
81     // Instance AliPythia
82     AliPythia::Instance(); 
83     AliDpmJetRndm::SetDpmJetRandom(GetRandom());
84 }
85
86
87 //______________________________________________________________________________
88 AliGenDPMjet::~AliGenDPMjet()
89 {
90 // Destructor
91     delete fParticles;
92 }
93
94
95 //______________________________________________________________________________
96 void AliGenDPMjet::Init()
97 {
98 // Initialization
99     
100     SetMC(new TDPMjet(fAProjectile, fZProjectile, fATarget, fZTarget, 
101                       fBeamEn,fEnergyCMS));
102
103     fDPMjet=(TDPMjet*) fgMCEvGen;
104     //
105     // **** Flag to force central production
106     // fICentr=1. central production forced 
107     // fICentr<0 && fICentr>-100 -> bmin = fMinImpactParam, bmax = fMaxImpactParam        
108     // fICentr<-99 -> fraction of x-sec. = XSFRAC                 
109     // fICentr=-1. -> evaporation/fzc suppressed                  
110     // fICentr<-1. -> evaporation/fzc suppressed                  
111     fDPMjet->SetfFCentr(fICentr);  
112     fDPMjet->SetbRange(fMinImpactParam,fMaxImpactParam); 
113     
114 //
115 //  Initialize DPMjet  
116 //    
117     fDPMjet->Initialize();
118     //if (fEvaluate) EvaluateCrossSections();
119
120     //  Initialize AliIonPDGCodes to define ion PDG codes
121     /*AliIonPDGCodes *PDGcodes = new AliIonPDGCodes();
122     PDGcodes->AddParticlesToPdgDataBase();
123     PDGcodes->MapPDGGEant3Codes();*/
124 }
125
126
127 //______________________________________________________________________________
128 void AliGenDPMjet::Generate()
129 {
130 // Generate one event
131
132   Float_t polar[3]    =   {0,0,0};
133   Float_t origin[3]   =   {0,0,0};
134   Float_t origin0[3]  =   {0,0,0};
135   Float_t p[3], random[6];
136   Float_t tof;
137
138 //  converts from mm/c to s
139   const Float_t kconv = 0.001/2.999792458e8;
140   Int_t nt  = 0;
141   Int_t jev = 0;
142   Int_t j, kf, ks, imo;
143   kf = 0;
144         
145   fTrials = 0;
146   for (j = 0;j < 3; j++) origin0[j] = fOrigin[j];
147   if(fVertexSmear == kPerEvent) {
148       Float_t dv[3];
149       dv[2] = 1.e10;
150       while(TMath::Abs(dv[2])>fCutVertexZ*fOsigma[2]) {
151           Rndm(random,6);
152           for (j=0; j < 3; j++) {
153               dv[j] = fOsigma[j]*TMath::Cos(2*random[2*j]*TMath::Pi())*
154                   TMath::Sqrt(-2*TMath::Log(random[2*j+1]));
155           }
156       }
157       for (j=0; j < 3; j++) origin0[j] += dv[j];
158   }
159   
160   while(1)
161   {
162 //    Generate one event
163 // --------------------------------------------------------------------------
164       fSpecn = 0;  
165       fSpecp = 0;
166 // --------------------------------------------------------------------------
167       fDPMjet->GenerateEvent();
168       fTrials++;
169
170       fDPMjet->ImportParticles(fParticles,"All");      
171       if (fLHC) Boost();
172
173       // Temporaneo
174       fGenImpPar = fDPMjet->GetBImpac();
175       
176       Int_t np = fParticles->GetEntriesFast();
177       printf("\n **************************************************%d\n",np);
178       Int_t nc=0;
179       if (np==0) continue;
180       Int_t i;
181       Int_t* newPos     = new Int_t[np];
182       Int_t* pSelected  = new Int_t[np];
183
184       for (i = 0; i<np; i++) {
185           newPos[i]    = i;
186           pSelected[i] = 0;
187       }
188       
189 //      Get event vertex
190       
191       //TParticle *iparticle = (TParticle*) fParticles->At(0);
192       fVertex[0] = origin0[0];
193       fVertex[1] = origin0[1];  
194       fVertex[2] = origin0[2];
195       //for(int jj=0; jj<3; jj++) printf("      fEventVertex[%d] = %f\n",jj,fEventVertex[jj]);
196       
197 //      First select parent particles
198
199       for (i = 0; i<np; i++) {
200           TParticle *iparticle = (TParticle *) fParticles->At(i);
201           //if(!iparticle) iparticle->Print("");
202
203 // Is this a parent particle ?
204
205           if (Stable(iparticle)) continue;
206
207           Bool_t  selected             =  kTRUE;
208           Bool_t  hasSelectedDaughters =  kFALSE;
209           
210           kf = iparticle->GetPdgCode();
211           if (kf == 92) continue;
212           ks = iparticle->GetStatusCode();
213 // No initial state partons
214           if (ks==21) continue;
215             
216           if (!fSelectAll) selected = KinematicSelection(iparticle, 0) && 
217                                SelectFlavor(kf);
218           hasSelectedDaughters = DaughtersSelection(iparticle);
219
220 // Put particle on the stack if it is either selected or 
221 // it is the mother of at least one seleted particle
222
223           if (selected || hasSelectedDaughters) {
224               nc++;
225               pSelected[i] = 1;
226           } // selected
227       } // particle loop parents
228
229 // Now select the final state particles
230
231
232       for (i=0; i<np; i++) {
233           TParticle *iparticle = (TParticle *) fParticles->At(i);
234           //if(!iparticle) printf("     i = %d -> iparticle==0\n",i);
235
236 // Is this a final state particle ?
237
238           if (!Stable(iparticle)) continue;
239       
240           Bool_t  selected =  kTRUE;
241           kf = iparticle->GetPdgCode();
242           ks = iparticle->GetStatusCode();
243
244 // --------------------------------------------------------------------------
245 // Count spectator neutrons and protons (ks == 13, 14)
246           if(ks == 13 || ks == 14){
247               if(kf == kNeutron) fSpecn += 1;
248               if(kf == kProton)  fSpecp += 1;
249           }
250 // --------------------------------------------------------------------------
251
252           if (!fSelectAll) {
253               selected = KinematicSelection(iparticle,0)&&SelectFlavor(kf);
254               if (!fSpectators && selected) selected = (ks == 13 || ks == 14);
255           }
256
257 // Put particle on the stack if selected
258
259           if (selected) {
260               nc++;
261               pSelected[i] = 1;
262           } // selected
263       } // particle loop final state
264
265 // Write particles to stack
266
267       for (i = 0; i<np; i++) {
268           TParticle *  iparticle = (TParticle *) fParticles->At(i);
269           Bool_t  hasMother   = (iparticle->GetFirstMother()>=0);
270           Bool_t  hasDaughter = (iparticle->GetFirstDaughter()>=0);
271
272           if (pSelected[i]) {
273               
274               kf   = iparticle->GetPdgCode();         
275               if(kf==99999) continue;
276               ks   = iparticle->GetStatusCode();              
277               
278               p[0] = iparticle->Px();
279               p[1] = iparticle->Py();
280               p[2] = iparticle->Pz();
281               origin[0] = origin0[0]+iparticle->Vx()/10;
282               origin[1] = origin0[1]+iparticle->Vy()/10;
283               origin[2] = origin0[2]+iparticle->Vz()/10;
284               tof = kconv*iparticle->T();
285               
286               imo = -1;
287               TParticle* mother = 0;
288               if (hasMother) {
289                   imo = iparticle->GetFirstMother();
290                   mother = (TParticle *) fParticles->At(imo);
291                   imo = (mother->GetPdgCode() != 92) ? imo = newPos[imo] : -1;
292               } // if has mother   
293
294               // --- Nucleons and photons from nucleus evaporation
295               //if(ks==-1 && (kf==2212 || kf==2112 || kf==22)) imo=-1;
296               
297               // --- Offset added to nuclei, alpha particles, deuteron 
298               // --- and triton PDG codes (according to TGeant3 PDG codes)
299               //if(ks==1000 || ks==1001) kf += 10000000;
300               if(kf>10000) kf += 10000000;
301               
302               if(kf>10000 || ks==1000 || ks==1001) printf(" kf = %d, ks = %d imo = %d \n", kf,ks,imo);
303              
304
305               Bool_t tFlag = (fTrackIt && (ks == 1));
306               
307 /*            // --- Particle NOT to be tracked:
308               // --- (1) chains from valence quark system (idhkk=99999)
309               // --- (2) quark, diquarks and gluons (idhkk=1->8,...,idhkk=21)
310               if(kf==99999 || kf==1 || kf==2 || kf==3 || kf==4 || kf==5 ||
311                  kf==6 || kf==7 || kf==8 || kf==21 || 
312                  kf==1103 || kf==2101 || kf==2103 || kf==2203 || kf==3101 || 
313                  kf==3103 || kf==3201 || kf==3203 || kf==3303 || kf==4101 || 
314                  kf==4103 || kf==4201 || kf==4203 || kf==4301 || kf==4303 || 
315                  kf==4403 || kf==5101 || kf==5103 || kf==5201 || kf==5203 || 
316                  kf==5301 || kf==5303 || kf==5401 || kf==5403 || kf==5503) 
317                  tFlag=kFALSE;
318 */            
319               SetTrack(tFlag,imo,kf,p,origin,polar,tof,kPNoProcess,nt, 1., ks);
320               KeepTrack(nt);
321               newPos[i] = nt;
322           } // if selected
323       } // particle loop
324       delete[] newPos;
325       delete[] pSelected;
326       
327       printf("\n I've put %i particles on the stack \n",nc);
328       if (nc>0) {
329           jev += nc;
330           if (jev >= fNpart || fNpart == -1) {
331               printf("\n Trials: %i %i %i\n",fTrials, fNpart, jev);
332               break;
333           }
334           printf("\n Not breaking ### Trials: %i %i %i\n",fTrials, fNpart, jev);
335       }
336   } // event loop
337   MakeHeader();
338   SetHighWaterMark(nt);
339 }
340
341
342 //______________________________________________________________________________
343 void AliGenDPMjet::KeepFullEvent()
344 {
345     fKeep=1;
346 }
347
348
349 //______________________________________________________________________________
350 /*void AliGenDPMjet::EvaluateCrossSections()
351 {
352 //     Glauber Calculation of geometrical x-section
353 //
354     Float_t xTot       = 0.;          // barn
355     Float_t xTotHard   = 0.;          // barn 
356     Float_t xPart      = 0.;          // barn
357     Float_t xPartHard  = 0.;          // barn 
358     Float_t sigmaHard  = 0.1;         // mbarn
359     Float_t bMin       = 0.;
360     Float_t bMax       = fDPMjet->GetProjRadius()+fDPMjet->GetTargRadius();
361     const Float_t kdib = 0.2;
362     Int_t   kMax       = Int_t((bMax-bMin)/kdib)+1;
363
364
365     printf("\n Projectile Radius (fm): %f \n",fDPMjet->GetProjRadius());
366     printf("\n Target     Radius (fm): %f \n",fDPMjet->GetTargRadius());    
367     Int_t i;
368     Float_t oldvalue= 0.;
369
370     Float_t* b   = new Float_t[kMax];
371     Float_t* si1 = new Float_t[kMax];    
372     Float_t* si2 = new Float_t[kMax];    
373     
374     for (i = 0; i < kMax; i++)
375     {
376         Float_t xb  = bMin+i*kdib;
377         Float_t ov;
378         ov=fDPMjet->Profile(xb);
379         // ATT!->Manca la x-sec anel. nucleone-nucleone
380         Float_t gb  =  2.*0.01*fDPMjet->TMath::Pi()*kdib*xb*(1.-TMath::Exp(-fDPMjet->GetXSFrac()*ov));
381         Float_t gbh =  2.*0.01*fDPMjet->TMath::Pi()*kdib*xb*sigmaHard*ov;
382         xTot+=gb;
383         xTotHard += gbh;
384         if (xb > fMinImpactParam && xb < fMaxImpactParam)
385         {
386             xPart += gb;
387             xPartHard += gbh;
388         }
389         
390         if(oldvalue) if ((xTot-oldvalue)/oldvalue<0.0001) break;
391         oldvalue = xTot;
392         printf("\n Total cross section (barn): %d %f %f \n",i, xb, xTot);
393         printf("\n Hard  cross section (barn): %d %f %f \n\n",i, xb, xTotHard);
394         if (i>0) {
395             si1[i] = gb/kdib;
396             si2[i] = gbh/gb;
397             b[i]  = xb;
398         }
399     }
400
401     printf("\n Total cross section (barn): %f \n",xTot);
402     printf("\n Hard  cross section (barn): %f \n \n",xTotHard);
403     printf("\n Partial       cross section (barn): %f %f \n",xPart, xPart/xTot*100.);
404     printf("\n Partial  hard cross section (barn): %f %f \n",xPartHard, xPartHard/xTotHard*100.);
405
406 //  Store result as a graph
407     b[0] = 0;
408     si1[0] = 0;
409     si2[0]=si2[1];
410     
411     fDsigmaDb  = new TGraph(i, b, si1);
412     fDnDb      = new TGraph(i, b, si2);
413 }*/
414
415
416 //______________________________________________________________________________
417 Bool_t AliGenDPMjet::DaughtersSelection(TParticle* iparticle)
418 {
419 //
420 // Looks recursively if one of the daughters has been selected
421 //
422 //    printf("\n Consider daughters %d:",iparticle->GetPdgCode());
423     Int_t imin = -1;
424     Int_t imax = -1;
425     Int_t i;
426     Bool_t hasDaughters = (iparticle->GetFirstDaughter() >=0);
427     Bool_t selected = kFALSE;
428     if (hasDaughters) {
429         imin = iparticle->GetFirstDaughter();
430         imax = iparticle->GetLastDaughter();       
431         for (i = imin; i <= imax; i++){
432             TParticle *  jparticle = (TParticle *) fParticles->At(i);   
433             Int_t ip = jparticle->GetPdgCode();
434             if (KinematicSelection(jparticle,0)&&SelectFlavor(ip)) {
435                 selected=kTRUE; break;
436             }
437             if (DaughtersSelection(jparticle)) {selected=kTRUE; break; }
438         }
439     } else {
440         return kFALSE;
441     }
442     return selected;
443 }
444
445
446
447 //______________________________________________________________________________
448 Bool_t AliGenDPMjet::SelectFlavor(Int_t pid)
449 {
450 // Select flavor of particle
451 // 0: all
452 // 4: charm and beauty
453 // 5: beauty
454     Bool_t res = 0;
455     
456     if (fFlavor == 0) {
457         res = kTRUE;
458     } else {
459         Int_t ifl = TMath::Abs(pid/100);
460         if (ifl > 10) ifl/=10;
461         res = (fFlavor == ifl);
462     }
463 //
464 //  This part if gamma writing is inhibited
465     if (fNoGammas) 
466         res = res && (pid != kGamma && pid != kPi0);
467 //
468     return res;
469 }
470
471 //______________________________________________________________________________
472 Bool_t AliGenDPMjet::Stable(TParticle*  particle)
473 {
474 // Return true for a stable particle
475 //
476     
477 //    if (particle->GetFirstDaughter() < 0 ) return kTRUE;
478     if (particle->GetStatusCode() == 1) return kTRUE;
479     else return kFALSE;
480
481 }
482
483
484
485 //______________________________________________________________________________
486 /*void AliGenDPMjet::Boost(TClonesArray* particles)
487 {
488 //
489 // Boost cms into LHC lab frame
490 //
491     Double_t dy    = - 0.5 * TMath::Log(Double_t(fZProjectile) * Double_t(fATarget) / 
492                                       (Double_t(fZTarget)    * Double_t(fAProjectile)));
493     Double_t beta  = TMath::TanH(dy);
494     Double_t gamma = 1./TMath::Sqrt(1.-beta*beta);
495     Double_t gb    = gamma * beta;
496
497     printf("\n Boosting particles to lab frame %f %f %f", dy, beta, gamma);
498     
499     Int_t i;
500     Int_t np = particles->GetEntriesFast();
501     for (i = 0; i < np; i++) 
502     {
503         TParticle* iparticle = (TParticle*) fParticles->At(i);
504
505         Double_t e   = iparticle->Energy();
506         Double_t px  = iparticle->Px();
507         Double_t py  = iparticle->Py();
508         Double_t pz  = iparticle->Pz();
509
510         Double_t eb  = gamma * e -      gb * pz;
511         Double_t pzb =   -gb * e +   gamma * pz;
512
513         iparticle->SetMomentum(px, py, pzb, eb);
514     }
515 }*/
516
517
518
519 //______________________________________________________________________________
520 void AliGenDPMjet::MakeHeader()
521 {
522 // Builds the event header, to be called after each event
523     AliGenEventHeader* header = new AliGenDPMjetEventHeader("DPMJET");
524     ((AliGenDPMjetEventHeader*) header)->SetNProduced(fDPMjet->GetNumStablePc());
525     ((AliGenDPMjetEventHeader*) header)->SetImpactParameter(fDPMjet->GetBImpac());
526     ((AliGenDPMjetEventHeader*) header)->SetTotalEnergy(fDPMjet->GetTotEnergy());
527 //    ((AliGenDPMjetEventHeader*) header)->SetWounded(fDPMjet->GetProjWounded(),
528 //                                                  fDPMjet->GetTargWounded());
529     ((AliGenDPMjetEventHeader*) header)->SetParticipants(fDPMjet->GetfIp(), 
530                                                          fDPMjet->GetfIt());
531     /*((AliGenDPMjetEventHeader*) header)->SetCollisions(fDPMjet->GetN0(),
532                                                        fDPMjet->GetN01(),
533                                                        fDPMjet->GetN10(),
534                                                        fDPMjet->GetN11());*/
535     ((AliGenDPMjetEventHeader*) header)->SetSpectators(fSpecn, fSpecp);
536
537 // Bookkeeping for kinematic bias
538     ((AliGenDPMjetEventHeader*) header)->SetTrials(fTrials);
539 // Event Vertex
540     header->SetPrimaryVertex(fVertex);
541     gAlice->SetGenEventHeader(header);    
542 }
543
544
545
546 //______________________________________________________________________________
547 AliGenDPMjet& AliGenDPMjet::operator=(const  AliGenDPMjet& rhs)
548 {
549 // Assignment operator
550     return *this;
551 }
552
553
554
555 //______________________________________________________________________________