]> git.uio.no Git - u/mrichter/AliRoot.git/blame - RALICE/AliCollider.cxx
removing temporary file
[u/mrichter/AliRoot.git] / RALICE / AliCollider.cxx
CommitLineData
fdbea0ce 1// $Id$
2
3///////////////////////////////////////////////////////////////////////////
4// Class AliCollider
5// Pythia based universal physics event generator.
6// This event class is derived from TPythia6 and has some extensions to
7// support also generation of nucleus-nucleus interactions and to allow
8// investigation of the effect of detector resolving power.
9// Furthermore, the produced event information is provided in a format
10// using the AliEvent structure.
11// For the produced AliTrack objects, the particle ID code is set to the
12// Pythia KF value, which is compatible with the PDG identifier.
13// This will allow a direct analysis of the produced data using the
14// Ralice physics analysis tools.
15//
16// For further details concerning the produced output structure,
17// see the docs of the memberfunctions SetVertexMode and SetResolution.
18//
19// Example job of minimum biased Pb+Pb interactions :
20// --------------------------------------------------
21// {
22// gSystem->Load("libEG");
23// gSystem->Load("libEGPythia6");
24// gSystem->Load("ralice");
25//
26// AliCollider* gen=new AliCollider();
27//
28// gen->SetOutputFile("test.root");
29// gen->SetVertexMode(3);
30// gen->SetResolution(1e-4); // 1 micron vertex resolution
31//
32// gen->SetRunNumber(1);
33//
34// Int_t zp=82;
35// Int_t ap=208;
36// Int_t zt=82;
37// Int_t at=208;
38//
39// gen->Init("fixt",zp,ap,zt,at,158);
40//
41// Int_t nevents=5;
42//
43// AliRandom rndm;
44// Float_t* rans=new Float_t[nevents];
45// rndm.Uniform(rans,nevents,2,ap+at);
46// Int_t npart;
47// for (Int_t i=0; i<nevents; i++)
48// {
49// npart=rans[i];
50// gen->MakeEvent(npart);
51//
52// AliEvent* evt=gen->GetEvent();
53//
54// evt->List();
55// }
56//
57// gen->EndRun();
58// }
59//
60//
61// Example job of a cosmic nu+p atmospheric interaction.
62// -----------------------------------------------------
63// {
64// gSystem->Load("libEG");
65// gSystem->Load("libEGPythia6");
66// gSystem->Load("ralice");
67//
68// AliCollider* gen=new AliCollider();
69//
70// gen->SetOutputFile("test.root");
71//
72// gen->SetRunNumber(1);
73//
74// gen->Init("fixt","nu_mu","p",1e11);
75//
76// Int_t nevents=10;
77//
78// for (Int_t i=0; i<nevents; i++)
79// {
80// gen->MakeEvent(0,1);
81//
82// AliEvent* evt=gen->GetEvent();
83//
84// evt->Info();
85// }
86//
87// gen->EndRun();
88// }
89//
90//
91//--- Author: Nick van Eijndhoven 22-nov-2002 Utrecht University
92//- Modified: NvE $Date$ Utrecht University
93///////////////////////////////////////////////////////////////////////////
94
95#include "AliCollider.h"
96
97ClassImp(AliCollider) // Class implementation to enable ROOT I/O
98
99AliCollider::AliCollider()
100{
101// Default constructor.
102// All variables initialised to default values.
103 fVertexmode=0; // No vertex structure creation
104 fResolution=1e-5; // Standard resolution is 0.1 micron
105 fRunnum=0;
106 fEventnum=0;
107 fPrintfreq=1;
108
109 fEvent=0;
110
111 fFrame="none";
112 fWin=0;
113
114 fNucl=0;
115 fZproj=0;
116 fAproj=0;
117 fZtarg=0;
118 fAtarg=0;
119 fFracpp=0;
120 fFracnp=0;
121 fFracpn=0;
122 fFracnn=0;
123
124 fOutFile=0;
125 fOutTree=0;
126}
127///////////////////////////////////////////////////////////////////////////
128AliCollider::~AliCollider()
129{
130// Default destructor
131 if (fEvent)
132 {
133 delete fEvent;
134 fEvent=0;
135 }
136 if (fOutFile)
137 {
138 delete fOutFile;
139 fOutFile=0;
140 }
141 if (fOutTree)
142 {
143 delete fOutTree;
144 fOutTree=0;
145 }
146}
147///////////////////////////////////////////////////////////////////////////
148void AliCollider::SetOutputFile(TString s)
149{
150// Create the output file containing all the data in ROOT output format.
151 if (fOutFile)
152 {
153 delete fOutFile;
154 fOutFile=0;
155 }
156 fOutFile=new TFile(s.Data(),"RECREATE","AliCollider data");
157
158 if (fOutTree)
159 {
160 delete fOutTree;
161 fOutTree=0;
162 }
163 fOutTree=new TTree("T","AliCollider event data");
164
165 Int_t bsize=32000;
166 Int_t split=0;
167 fOutTree->Branch("Events","AliEvent",&fEvent,bsize,split);
168}
169///////////////////////////////////////////////////////////////////////////
170void AliCollider::SetVertexMode(Int_t mode)
171{
172// Set the mode of the vertex structure creation.
173//
174// By default all generated tracks will only appear in the AliEvent
175// structure without any primary (and secondary) vertex structure.
176// The user can build the vertex structure if he/she wants by means
177// of the beginpoint location of each AliTrack.
178//
179// However, one can also let AliCollider automatically create
180// the primary (and secondary) vertex structure(s).
181// In this case the primary vertex is given Id=1 and all sec. vertices
182// are given Id's 2,3,4,....
183// All vertices are created as standalone entities in the AliEvent structure
184// without any linking between the various vertices.
185// For this automated process, the user-selected resolution
186// (see SetResolution) is used to decide whether or not certain vertex
187// locations can be resolved.
188// In case no vertex creation is selected (i.e. the default mode=0),
189// the value of the resolution is totally irrelevant.
190//
191// The user can also let AliCollider automatically connect the sec. vertices
192// to the primary vertex (i.e. mode=3). This process will also automatically
193// generate the tracks connecting the vertices.
194// Note that the result of the mode=3 operation may be very sensitive to
195// the resolution parameter. Therefore, no attempt is made to distinguish
196// between secondary, tertiary etc... vertices. All sec. vertices are
197// linked to the primary one.
198//
199// Irrespective of the selected mode, all generated tracks can be obtained
200// directly from the AliEvent structure.
201// In case (sec.) vertex creation is selected, all generated vertices can
202// also be obtained directly from the AliEvent structure.
203// These (sec.) vertices contain only the corresponding pointers to the various
204// tracks which are stored in the AliEvent structure.
205//
206// Overview of vertex creation modes :
207// -----------------------------------
208// mode = 0 ==> No vertex structure will be created
209// 1 ==> Only primary vertex structure will be created
210// 2 ==> Unconnected primary and secondary vertices will be created
211// 3 ==> Primary and secondary vertices will be created where all the
212// sec. vertices will be connected to the primary vertex.
213// Also the vertex connecting tracks will be automatically
214// generated.
215//
216 if (mode<0 || mode >3)
217 {
218 cout << " *AliCollider::SetVertexMode* Invalid argument mode : " << mode << endl;
219 fVertexmode=0;
220 }
221 else
222 {
223 fVertexmode=mode;
224 }
225}
226///////////////////////////////////////////////////////////////////////////
227Int_t AliCollider::GetVertexMode()
228{
229// Provide the current mode for vertex structure creation.
230 return fVertexmode;
231}
232///////////////////////////////////////////////////////////////////////////
233void AliCollider::SetResolution(Double_t res)
234{
235// Set the resolution (in cm) for resolving (sec.) vertices.
236// By default this resolution is set to 0.1 micron.
237// Note : In case no vertex creation has been selected, the value of
238// the resolution is totally irrelevant.
239 fResolution=fabs(res);
240}
241///////////////////////////////////////////////////////////////////////////
242Double_t AliCollider::GetResolution()
243{
244// Provide the current resolution (in cm) for resolving (sec.) vertices.
245 return fResolution;
246}
247///////////////////////////////////////////////////////////////////////////
248void AliCollider::SetRunNumber(Int_t run)
249{
250// Set the user defined run number.
251// By default the run number is set to 0.
252 fRunnum=run;
253}
254///////////////////////////////////////////////////////////////////////////
255Int_t AliCollider::GetRunNumber()
256{
257// Provide the user defined run number.
258 return fRunnum;
259}
260///////////////////////////////////////////////////////////////////////////
261void AliCollider::SetPrintFreq(Int_t n)
262{
263// Set the print frequency for every 'n' events.
264// By default the printfrequency is set to 1 (i.e. every event).
265 fPrintfreq=n;
266}
267///////////////////////////////////////////////////////////////////////////
268Int_t AliCollider::GetPrintFreq()
269{
270// Provide the user selected print frequency.
271 return fPrintfreq;
272}
273///////////////////////////////////////////////////////////////////////////
274void AliCollider::Init(char* frame,char* beam,char* target,Float_t win)
275{
276// Initialisation of the underlying Pythia generator package.
277// This routine just invokes TPythia6::Initialize(...) and the arguments
278// have the corresponding meaning.
279// The event number is reset to 0.
280 fEventnum=0;
281 fNucl=0;
282 fFrame=frame;
283 fWin=win;
284 Initialize(frame,beam,target,win);
285}
286///////////////////////////////////////////////////////////////////////////
287void AliCollider::Init(char* frame,Int_t zp,Int_t ap,Int_t zt,Int_t at,Float_t win)
288{
289// Initialisation of the underlying Pythia generator package for the generation
290// of nucleus-nucleus interactions.
291// In addition to the Pythia standard arguments 'frame' and 'win', the user
292// can specify here (Z,A) values of the projectile and target nuclei and the number
293// 'npart' of the participant nucleons for this collision.
294// The event number is reset to 0.
295 fEventnum=0;
296 fNucl=1;
297 fFrame=frame;
298 fWin=win;
299 fZproj=0;
300 fAproj=0;
301 fZtarg=0;
302 fAtarg=0;
303 fFracpp=0;
304 fFracnp=0;
305 fFracpn=0;
306 fFracnn=0;
307
308 if (ap<1 || at<1 || zp>ap || zt>at)
309 {
310 cout << " *AliCollider::Init* Invalid input value(s). Zproj = " << zp
311 << " Aproj = " << ap << " Ztarg = " << zt << " Atarg = " << at << endl;
312 return;
313 }
314
315 fZproj=zp;
316 fAproj=ap;
317 fZtarg=zt;
318 fAtarg=at;
319
320 cout << " *AliCollider::Init* Nucleus-Nucleus generator initialisation." << endl;
321 cout << " Zproj = " << zp << " Aproj = " << ap << " Ztarg = " << zt << " Atarg = " << at
322 << " Frame = " << frame << " Energy = " << win
323 << endl;
324}
325///////////////////////////////////////////////////////////////////////////
326void AliCollider::GetFractions(Float_t zp,Float_t ap,Float_t zt,Float_t at)
327{
328// Determine the fractions for the various N-N collision processes.
329// The various processes are : p+p, n+p, p+n and n+n.
330 if (zp<0) zp=0;
331 if (zt<0) zt=0;
332
333 fFracpp=0;
334 fFracnp=0;
335 fFracpn=0;
336 fFracnn=0;
337
338 if (ap>0 && at>0)
339 {
340 fFracpp=(zp/ap)*(zt/at);
341 fFracnp=(1.-zp/ap)*(zt/at);
342 fFracpn=(zp/ap)*(1.-zt/at);
343 fFracnn=(1.-zp/ap)*(1.-zt/at);
344 }
345}
346///////////////////////////////////////////////////////////////////////////
347void AliCollider::MakeEvent(Int_t npt,Int_t mlist,Int_t medit)
348{
349// Generate one event.
350// In case of a nucleus-nucleus interaction, the argument 'npt' denotes
351// the number of participant nucleons.
352// In case of a standard Pythia run for 'elementary' particle interactions,
353// the value of npt is totally irrelevant.
354// The argument 'medit' denotes the edit mode used for Pyedit().
355// By default, only 'stable' final particles are kept (i.e. medit=1).
356// The argument 'mlist' denotes the list mode used for Pylist().
357// By default, no listing is produced (i.e. mlist=0).
358
359 fEventnum++;
360
361 // Counters for the various (proj,targ) combinations : p+p, n+p, p+n and n+n
362 Int_t ncols[4]={0,0,0,0};
363
364 if (fNucl)
365 {
366 if (npt<1 || npt>(fAproj+fAtarg))
367 {
368 cout << " *AliCollider::MakeEvent* Invalid input value. npt = " << npt
369 << " Aproj = " << fAproj << " Atarg = " << fAtarg << endl;
370 return;
371 }
372
373 // Determine the number of nucleon-nucleon collisions
374 Int_t ncol=npt/2.;
375 if (npt%2 && fRan.Uniform()>0.5) ncol+=1;
376
377 // Determine the number of the various types of N+N interactions
378 Int_t zp=fZproj;
379 Int_t ap=fAproj;
380 Int_t zt=fZtarg;
381 Int_t at=fAtarg;
382 Int_t maxa=2; // Indicator whether proj (1) or target (2) has maximal A left
383 if (ap>at) maxa=1;
384 Float_t* rans=new Float_t[ncol];
385 fRan.Uniform(rans,ncol);
386 Float_t rndm=0;
387 for (Int_t i=0; i<ncol; i++)
388 {
389 GetFractions(zp,ap,zt,at);
390 rndm=rans[i];
391 if (rndm<=fFracpp) // p+p interaction
392 {
393 ncols[0]++;
394 if (maxa=2)
395 {
396 at--;
397 zt--;
398 }
399 else
400 {
401 ap--;
402 zp--;
403 }
404 }
405 if (rndm>fFracpp && rndm<=(fFracpp+fFracnp)) // n+p interaction
406 {
407 ncols[1]++;
408 if (maxa=2)
409 {
410 at--;
411 zt--;
412 }
413 else
414 {
415 ap--;
416 }
417 }
418 if (rndm>(fFracpp+fFracnp) && rndm<=(fFracpp+fFracnp+fFracpn)) // p+n interaction
419 {
420 ncols[2]++;
421 if (maxa=2)
422 {
423 at--;
424 }
425 else
426 {
427 ap--;
428 zp--;
429 }
430 }
431 if (rndm>(fFracpp+fFracnp+fFracpn)) // n+n interaction
432 {
433 ncols[3]++;
434 if (maxa=2)
435 {
436 at--;
437 }
438 else
439 {
440 ap--;
441 }
442 }
443 }
444 delete [] rans;
445
446 if (!(fEventnum%fPrintfreq))
447 {
448 cout << " *AliCollider::MakeEvent* Run : " << fRunnum << " Event : " << fEventnum
449 << endl;
450 cout << " npart = " << npt << " ncol = " << ncol
451 << " ncolpp = " << ncols[0] << " ncolnp = " << ncols[1]
452 << " ncolpn = " << ncols[2] << " ncolnn = " << ncols[3] << endl;
453 }
454
455 }
456
457 if (!fEvent)
458 {
459 fEvent=new AliEvent();
460 fEvent->SetOwner();
461 }
462
463 fEvent->Reset();
464 fEvent->SetRunNumber(fRunnum);
465 fEvent->SetEventNumber(fEventnum);
466
467 AliVertex vert;
468 if (fVertexmode)
469 {
470 // Make sure the primary vertex gets correct location and Id=1
471 vert.SetId(1);
472 vert.SetTrackCopy(0);
473 vert.SetVertexCopy(0);
474 fEvent->AddVertex(vert,0);
475 }
476
477 AliTrack t;
478 Ali3Vector p;
479 AliPosition r,rx;
480 Float_t v[3];
481
482 Int_t kf=0,kc=0;
483 Float_t charge=0,mass=0;
484
485 TMCParticle* part=0;
486
487 Int_t ntypes=4;
488
489 // Singular settings for a normal Pythia elementary particle interation
490 if (!fNucl)
491 {
492 ntypes=1;
493 ncols[0]=1;
494 }
495
496 // Generate all the various collisions
497 Int_t npart=0,ntk=0;
498 Double_t dist=0;
499 for (Int_t itype=0; itype<ntypes; itype++)
500 {
501 if (fNucl)
502 {
503 if (itype==0 && ncols[itype]) Initialize(fFrame,"p","p",fWin);
504 if (itype==1 && ncols[itype]) Initialize(fFrame,"n","p",fWin);
505 if (itype==2 && ncols[itype]) Initialize(fFrame,"p","n",fWin);
506 if (itype==3 && ncols[itype]) Initialize(fFrame,"n","n",fWin);
507 }
508 for (Int_t jcol=0; jcol<ncols[itype]; jcol++)
509 {
510 GenerateEvent();
511
512 Pyedit(medit); // Define which particles are to be kept
513
514 if (mlist) Pylist(mlist);
515
516 ImportParticles();
517 npart=0;
518 if (fParticles) npart=fParticles->GetEntries();
519
520 for (Int_t jpart=0; jpart<npart; jpart++)
521 {
522 part=(TMCParticle*)fParticles->At(jpart);
523 if (!part) continue;
524
525 kf=part->GetKF();
526 kc=Pycomp(kf);
527
528 charge=GetKCHG(kc,1)/3.;
529 if (kf<0) charge*=-1;
530 mass=GetPMAS(kc,1);
531
532 // 3-momentum in GeV/c
533 v[0]=part->GetPx();
534 v[1]=part->GetPy();
535 v[2]=part->GetPz();
536 p.SetVector(v,"car");
537
538 // Production location in cm.
539 v[0]=(part->GetVx())/10;
540 v[1]=(part->GetVy())/10;
541 v[2]=(part->GetVz())/10;
542 r.SetVector(v,"car");
543
544 ntk++;
545
546 t.Reset();
547 t.SetId(ntk);
548 t.SetParticleCode(kf);
549 t.SetCharge(charge);
550 t.SetMass(mass);
551 t.Set3Momentum(p);
552 t.SetBeginPoint(r);
553
554 fEvent->AddTrack(t);
555
556 // Build the vertex structures if requested
557 if (fVertexmode)
558 {
559 // Check if track belongs within the resolution to an existing vertex
560 Int_t add=0;
561 for (Int_t jv=1; jv<=fEvent->GetNvertices(); jv++)
562 {
563 AliVertex* vx=fEvent->GetVertex(jv);
564 if (vx)
565 {
566 rx=vx->GetPosition();
567 dist=rx.GetDistance(r);
568 if (dist < fResolution)
569 {
570 AliTrack* tx=fEvent->GetIdTrack(ntk);
571 if (tx)
572 {
573 vx->AddTrack(tx);
574 add=1;
575 }
576 }
577 }
578 if (add) break; // No need to look further for vertex candidates
579 }
580
581 // If track was not close enough to an existing vertex
582 // a new secondary vertex is created
583 if (!add && fVertexmode>1)
584 {
585 AliTrack* tx=fEvent->GetIdTrack(ntk);
586 if (tx)
587 {
588 vert.Reset();
589 vert.SetTrackCopy(0);
590 vert.SetVertexCopy(0);
591 vert.SetId((fEvent->GetNvertices())+1);
592 vert.SetPosition(r);
593 vert.AddTrack(tx);
594 fEvent->AddVertex(vert,0);
595 }
596 }
597 }
598 } // End of loop over the produced particles for each collision
599 } // End of loop over number of collisions for each type
600 } // End of loop over collision types
601
602 // Link sec. vertices to primary if requested
603 // Note that also the connecting tracks are automatically created
604 if (fVertexmode>2)
605 {
606 AliVertex* vp=fEvent->GetIdVertex(1); // Primary vertex
607 if (vp)
608 {
609 for (Int_t i=2; i<=fEvent->GetNvertices(); i++)
610 {
611 AliVertex* vx=fEvent->GetVertex(i);
612 if (vx)
613 {
614 if (vx->GetId() != 1) vp->AddVertex(vx);
615 }
616 }
617 }
618 }
619
620 if (mlist) cout << endl; // Create empty output line after the event
621 if (fOutTree) fOutTree->Fill();
622}
623///////////////////////////////////////////////////////////////////////////
624AliEvent* AliCollider::GetEvent()
625{
626// Provide pointer to the generated event structure.
627 return fEvent;
628}
629///////////////////////////////////////////////////////////////////////////
630void AliCollider::EndRun()
631{
632// Properly close the output file (if needed).
633 if (fOutFile)
634 {
635 fOutFile->Write();
636 fOutFile->Close();
637 cout << " *AliCollider::EndRun* Output file correctly closed." << endl;
638 }
639}
640///////////////////////////////////////////////////////////////////////////