1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 // $Id: AliCollider.cxx,v 1.12 2004/05/04 15:33:04 nick Exp $
18 ///////////////////////////////////////////////////////////////////////////
20 // Pythia based universal physics event generator.
21 // This class is derived from TPythia6 and has some extensions to
22 // support also generation of nucleus-nucleus interactions and to allow
23 // investigation of the effect of detector resolving power.
24 // Furthermore, the produced event information is provided in a format
25 // using the AliEvent structure.
26 // For the produced AliTrack objects, the particle ID code is set to the
27 // Pythia KF value, which is compatible with the PDG identifier.
28 // This will allow a direct analysis of the produced data using the
29 // Ralice physics analysis tools.
31 // For further details concerning the produced output structure,
32 // see the docs of the memberfunctions SetVertexMode and SetResolution.
34 // Example job of minimum biased Pb+Pb interactions :
35 // --------------------------------------------------
37 // gSystem->Load("libEG");
38 // gSystem->Load("libEGPythia6");
39 // gSystem->Load("ralice");
41 // AliCollider* gen=new AliCollider();
43 // gen->SetOutputFile("test.root");
44 // gen->SetVertexMode(3);
45 // gen->SetResolution(1e-6); // 1 micron vertex resolution
47 // gen->SetRunNumber(1);
54 // gen->Init("fixt",zp,ap,zt,at,158);
56 // gen->SetTitle("SPS Pb-Pb collision at 158A GeV/c beam energy");
61 // Float_t* rans=new Float_t[nevents];
62 // rndm.Uniform(rans,nevents,2,ap+at);
64 // for (Int_t i=0; i<nevents; i++)
67 // gen->MakeEvent(npart);
69 // AliEvent* evt=gen->GetEvent();
78 // Example job of a cosmic nu+p atmospheric interaction.
79 // -----------------------------------------------------
81 // gSystem->Load("libEG");
82 // gSystem->Load("libEGPythia6");
83 // gSystem->Load("ralice");
85 // AliCollider* gen=new AliCollider();
87 // gen->SetOutputFile("test.root");
89 // gen->SetRunNumber(1);
91 // gen->Init("fixt","nu_mu","p",1e11);
93 // gen->SetTitle("Atmospheric nu_mu-p interaction at 1e20 eV");
97 // for (Int_t i=0; i<nevents; i++)
99 // gen->MakeEvent(0,1);
101 // AliEvent* evt=gen->GetEvent();
110 //--- Author: Nick van Eijndhoven 22-nov-2002 Utrecht University
111 //- Modified: NvE $Date: 2004/05/04 15:33:04 $ Utrecht University
112 ///////////////////////////////////////////////////////////////////////////
115 #include "AliCollider.h"
116 #include "Riostream.h"
118 ClassImp(AliCollider) // Class implementation to enable ROOT I/O
120 AliCollider::AliCollider() : TPythia6()
122 // Default constructor.
123 // All variables initialised to default values.
125 // Some Pythia default MC parameters are automatically modified to provide
126 // more suitable running conditions for soft processes in view of
127 // nucleus-nucleus interactions and astrophysical processes.
128 // The user may initialise the generator with all the default Pythia
129 // parameters and obtain full user control to modify the settings by means
130 // of the SetUserControl memberfunction.
132 // Refer to the SetElastic memberfunction for the inclusion of elastic
133 // and diffractive processes.
134 // By default these processes are not included.
136 fVertexmode=0; // No vertex structure creation
137 fResolution=1e-7; // Standard resolution is 0.1 micron
141 fUserctrl=0; // Automatic optimisation of some MC parameters
142 fElastic=0; // No elastic and diffractive processes
171 ///////////////////////////////////////////////////////////////////////////
172 AliCollider::~AliCollider()
174 // Default destructor
196 ///////////////////////////////////////////////////////////////////////////
197 void AliCollider::SetOutputFile(TString s)
199 // Create the output file containing all the data in ROOT output format.
205 fOutFile=new TFile(s.Data(),"RECREATE","AliCollider data");
212 fOutTree=new TTree("T","AliCollider event data");
216 fOutTree->Branch("Events","AliEvent",&fEvent,bsize,split);
218 ///////////////////////////////////////////////////////////////////////////
219 void AliCollider::SetVertexMode(Int_t mode)
221 // Set the mode of the vertex structure creation.
223 // By default all generated tracks will only appear in the AliEvent
224 // structure without any primary (and secondary) vertex structure.
225 // The user can build the vertex structure if he/she wants by means
226 // of the beginpoint location of each AliTrack.
228 // However, one can also let AliCollider automatically create
229 // the primary (and secondary) vertex structure(s).
230 // In this case the primary vertex is given Id=1 and all sec. vertices
231 // are given Id's 2,3,4,....
232 // All vertices are created as standalone entities in the AliEvent structure
233 // without any linking between the various vertices.
234 // For this automated process, the user-selected resolution
235 // (see SetResolution) is used to decide whether or not certain vertex
236 // locations can be resolved.
237 // In case no vertex creation is selected (i.e. the default mode=0),
238 // the value of the resolution is totally irrelevant.
240 // The user can also let AliCollider automatically connect the sec. vertices
241 // to the primary vertex (i.e. mode=3). This process will also automatically
242 // generate the tracks connecting the vertices.
243 // Note that the result of the mode=3 operation may be very sensitive to
244 // the resolution parameter. Therefore, no attempt is made to distinguish
245 // between secondary, tertiary etc... vertices. All sec. vertices are
246 // linked to the primary one.
248 // Irrespective of the selected mode, all generated tracks can be obtained
249 // directly from the AliEvent structure.
250 // In case (sec.) vertex creation is selected, all generated vertices can
251 // also be obtained directly from the AliEvent structure.
252 // These (sec.) vertices contain only the corresponding pointers to the various
253 // tracks which are stored in the AliEvent structure.
255 // Overview of vertex creation modes :
256 // -----------------------------------
257 // mode = 0 ==> No vertex structure will be created
258 // 1 ==> Only primary vertex structure will be created
259 // 2 ==> Unconnected primary and secondary vertices will be created
260 // 3 ==> Primary and secondary vertices will be created where all the
261 // sec. vertices will be connected to the primary vertex.
262 // Also the vertex connecting tracks will be automatically
265 if (mode<0 || mode >3)
267 cout << " *AliCollider::SetVertexMode* Invalid argument mode : " << mode << endl;
275 ///////////////////////////////////////////////////////////////////////////
276 Int_t AliCollider::GetVertexMode() const
278 // Provide the current mode for vertex structure creation.
281 ///////////////////////////////////////////////////////////////////////////
282 void AliCollider::SetResolution(Double_t res)
284 // Set the resolution (in meter) for resolving (sec.) vertices.
285 // By default this resolution is set to 0.1 micron.
286 // Note : In case no vertex creation has been selected, the value of
287 // the resolution is totally irrelevant.
288 fResolution=fabs(res);
290 ///////////////////////////////////////////////////////////////////////////
291 Double_t AliCollider::GetResolution() const
293 // Provide the current resolution (in meter) for resolving (sec.) vertices.
296 ///////////////////////////////////////////////////////////////////////////
297 void AliCollider::SetRunNumber(Int_t run)
299 // Set the user defined run number.
300 // By default the run number is set to 0.
303 ///////////////////////////////////////////////////////////////////////////
304 Int_t AliCollider::GetRunNumber() const
306 // Provide the user defined run number.
309 ///////////////////////////////////////////////////////////////////////////
310 void AliCollider::SetPrintFreq(Int_t n)
312 // Set the print frequency for every 'n' events.
313 // By default the printfrequency is set to 1 (i.e. every event).
316 ///////////////////////////////////////////////////////////////////////////
317 Int_t AliCollider::GetPrintFreq() const
319 // Provide the user selected print frequency.
322 ///////////////////////////////////////////////////////////////////////////
323 void AliCollider::SetUserControl(Int_t flag)
325 // Set the user control flag w.r.t. disabling automatic optimisation
326 // of some Pythia default MC parameters for soft interactions in view of
327 // nucleus-nucleus collisions and astrophysical processes.
328 // Flag = 0 : Limited user control (automatic optimisation enabled)
329 // 1 : Full user control (automatic optimisation disabled)
330 // By default the user control is set to 0 (i.e. automatic optimisation).
331 // See the Init() memberfunctions for further details w.r.t. the optimisations.
334 ///////////////////////////////////////////////////////////////////////////
335 Int_t AliCollider::GetUserControl() const
337 // Provide the value of the user control flag.
340 ///////////////////////////////////////////////////////////////////////////
341 void AliCollider::SetElastic(Int_t flag)
343 // Set the flag w.r.t. inclusion of elastic and diffractive processes.
344 // By default these processes are not included.
345 // Flag = 0 : Do not include elastic and diffractive processes
346 // 1 : Elastic and diffractive processes will be included
349 ///////////////////////////////////////////////////////////////////////////
350 Int_t AliCollider::GetElastic() const
352 // Provide the value of the control flag for elastic and diffractive processes.
355 ///////////////////////////////////////////////////////////////////////////
356 void AliCollider::Init(char* frame,char* beam,char* target,Float_t win)
358 // Initialisation of the underlying Pythia generator package.
359 // The event number is reset to 0.
360 // This routine just invokes TPythia6::Initialize(...) and the arguments
361 // have the corresponding meaning.
362 // Some Pythia default MC parameters are automatically modified to provide
363 // more suitable running conditions for soft processes in view of
364 // astrophysical processes.
365 // The optimisations consist of :
366 // * Usage of real photons for photon beams or targets
367 // * Minimum CMS energy of 3 GeV for the event
368 // * Activation of the default K factor values
369 // with separate settings for ordinary and color annihilation graphs.
370 // The user may initialise the generator with all the default Pythia
371 // parameters and obtain full user control to modify the settings by means
372 // of invoking the SetUserControl memberfunction before this initialisation.
373 // Note that the inclusion of elastic and diffractive processes is controlled
374 // by invokation of the SetElastic memberfunction before this initialisation,
375 // irrespective of the UserControl selection.
377 if (!fUserctrl) // Optimisation of some MC parameters
379 SetMSTP(14,10); // Real photons for photon beams or targets
380 SetPARP(2,3.); // Minimum CMS energy for the event
381 SetMSTP(33,2); // Activate K factor. Separate for ordinary and color annih. graphs
384 if (fElastic) SetMSEL(2); // Include low-Pt, elastic and diffractive events
390 Initialize(frame,beam,target,win);
393 cout << " *AliCollider::Init* Standard Pythia initialisation." << endl;
394 cout << " Beam particle : " << beam << " Target particle : " << target
395 << " Frame = " << frame << " Energy = " << win
398 ///////////////////////////////////////////////////////////////////////////
399 void AliCollider::Init(char* frame,Int_t zp,Int_t ap,Int_t zt,Int_t at,Float_t win)
401 // Initialisation of the underlying Pythia generator package for the generation
402 // of nucleus-nucleus interactions.
403 // The event number is reset to 0.
404 // In addition to the Pythia standard arguments 'frame' and 'win', the user
405 // can specify here (Z,A) values of the projectile and target nuclei.
407 // Note : The 'win' value denotes either the cms energy per nucleon-nucleon collision
408 // (i.e. frame="cms") or the momentum per nucleon in all other cases.
410 // Some Pythia default MC parameters are automatically modified to provide
411 // more suitable running conditions for soft processes in view of
412 // nucleus-nucleus interactions and astrophysical processes.
413 // The optimisations consist of :
414 // * Minimum CMS energy of 3 GeV for the event
415 // * Activation of the default K factor values
416 // with separate settings for ordinary and color annihilation graphs.
417 // The user may initialise the generator with all the default Pythia
418 // parameters and obtain full user control to modify the settings by means
419 // of invoking the SetUserControl memberfunction before this initialisation.
420 // Note that the inclusion of elastic and diffractive processes is controlled
421 // by invokation of the SetElastic memberfunction before this initialisation,
422 // irrespective of the UserControl selection.
424 if (!fUserctrl) // Optimisation of some MC parameters
426 SetPARP(2,3.); // Minimum CMS energy for the event
427 SetMSTP(33,2); // Activate K factor. Separate for ordinary and color annih. graphs
430 if (fElastic) SetMSEL(2); // Include low-Pt, elastic and diffractive events
445 if (ap<1 || at<1 || zp>ap || zt>at)
448 cout << " *AliCollider::Init* Invalid input value(s). Zproj = " << zp
449 << " Aproj = " << ap << " Ztarg = " << zt << " Atarg = " << at << endl;
459 cout << " *AliCollider::Init* Nucleus-Nucleus generator initialisation." << endl;
460 cout << " Zproj = " << zp << " Aproj = " << ap << " Ztarg = " << zt << " Atarg = " << at
461 << " Frame = " << frame << " Energy = " << win
464 ///////////////////////////////////////////////////////////////////////////
465 void AliCollider::GetFractions(Float_t zp,Float_t ap,Float_t zt,Float_t at)
467 // Determine the fractions for the various N-N collision processes.
468 // The various processes are : p+p, n+p, p+n and n+n.
479 fFracpp=(zp/ap)*(zt/at);
480 fFracnp=(1.-zp/ap)*(zt/at);
481 fFracpn=(zp/ap)*(1.-zt/at);
482 fFracnn=(1.-zp/ap)*(1.-zt/at);
485 ///////////////////////////////////////////////////////////////////////////
486 void AliCollider::MakeEvent(Int_t npt,Int_t mlist,Int_t medit)
488 // Generate one event.
489 // In case of a nucleus-nucleus interaction, the argument 'npt' denotes
490 // the number of participant nucleons.
491 // Normally also the spectator tracks will be stored into the event structure.
492 // The spectator tracks have a negative user Id to distinguish them from the
493 // ordinary generated tracks.
494 // In case the user has selected the creation of vertex structures, the spectator
495 // tracks will be linked to the primary vertex.
496 // However, specification of npt<0 will suppress the storage of spectator tracks.
497 // In the latter case abs(npt) will be taken as the number of participants.
498 // In case of a standard Pythia run for 'elementary' particle interactions,
499 // the value of npt is totally irrelevant.
501 // The argument 'mlist' denotes the list mode used for Pylist().
502 // Note : mlist<0 suppresses the invokation of Pylist().
503 // By default, no listing is produced (i.e. mlist=-1).
505 // The argument 'medit' denotes the edit mode used for Pyedit().
506 // Note : medit<0 suppresses the invokation of Pyedit().
507 // By default, only 'stable' final particles are kept (i.e. medit=1).
509 // In the case of a standard Pythia run concerning 'elementary' particle
510 // interactions, the projectile and target particle ID's for the created
511 // event structure are set to the corresponding Pythia KF codes.
512 // All the A and Z values are in that case set to zero.
513 // In case of a nucleus-nucleus interaction, the proper A and Z values for
514 // the projectile and target particles are set in the event structure.
515 // However, in this case both particle ID's are set to zero.
517 // Note : Only in case an event passed the selection criteria as specified
518 // via SelectEvent(), the event will appear on the output file.
529 // Counters for the various (proj,targ) combinations : p+p, n+p, p+n and n+n
530 Int_t ncols[4]={0,0,0,0};
540 if (npt<1 || npt>(fAproj+fAtarg))
542 cout << " *AliCollider::MakeEvent* Invalid input value. npt = " << npt
543 << " Aproj = " << fAproj << " Atarg = " << fAtarg << endl;
547 // Determine the number of nucleon-nucleon collisions
549 if (npt%2 && fRan.Uniform()>0.5) ncol+=1;
551 // Determine the number of the various types of N+N interactions
556 Int_t maxa=2; // Indicator whether proj (1) or target (2) has maximal A left
558 Float_t* rans=new Float_t[ncol];
559 fRan.Uniform(rans,ncol);
561 for (Int_t i=0; i<ncol; i++)
563 GetFractions(zp,ap,zt,at);
565 if (rndm<=fFracpp) // p+p interaction
579 if (rndm>fFracpp && rndm<=(fFracpp+fFracnp)) // n+p interaction
592 if (rndm>(fFracpp+fFracnp) && rndm<=(fFracpp+fFracnp+fFracpn)) // p+n interaction
605 if (rndm>(fFracpp+fFracnp+fFracpn)) // n+n interaction
621 if (!(fEventnum%fPrintfreq))
623 cout << " *AliCollider::MakeEvent* Run : " << fRunnum << " Event : " << fEventnum
627 cout << " npart = " << npt << " ncol = " << ncol
628 << " ncolpp = " << ncols[0] << " ncolnp = " << ncols[1]
629 << " ncolpn = " << ncols[2] << " ncolnn = " << ncols[3] << endl;
635 fEvent=new AliEvent();
637 fEvent->SetName(GetName());
638 fEvent->SetTitle(GetTitle());
642 fEvent->SetRunNumber(fRunnum);
643 fEvent->SetEventNumber(fEventnum);
650 Ali3Vector pproj,ptarg;
654 // Make sure the primary vertex gets correct location and Id=1
658 r.SetPosition(v,"car");
662 r.SetPositionErrors(v,"car");
665 vert.SetTrackCopy(0);
666 vert.SetVertexCopy(0);
668 fEvent->AddVertex(vert,0);
672 Float_t charge=0,mass=0;
677 // Singular settings for a normal Pythia elementary particle interation
684 // Generate all the various collisions
685 fSelect=0; // Flag to indicate whether the total event is selected or not
686 Int_t select=0; // Flag to indicate whether the sub-event is selected or not
687 Int_t first=1; // Flag to indicate the first collision process
691 for (Int_t itype=0; itype<ntypes; itype++)
695 if (itype==0 && ncols[itype]) Initialize(fFrame,"p","p",fWin);
696 if (itype==1 && ncols[itype]) Initialize(fFrame,"n","p",fWin);
697 if (itype==2 && ncols[itype]) Initialize(fFrame,"p","n",fWin);
698 if (itype==3 && ncols[itype]) Initialize(fFrame,"n","n",fWin);
700 for (Int_t jcol=0; jcol<ncols[itype]; jcol++)
705 if (select) fSelect=1;
707 if (first) // Store generator parameter information in the event structure
709 // Enter generator parameters as a device in the event
711 params.SetNameTitle("AliCollider","AliCollider generator parameters");
712 params.SetSlotName("Medit",1);
713 params.SetSlotName("Vertexmode",2);
714 params.SetSlotName("Resolution",3);
715 params.SetSlotName("Userctrl",4);
716 params.SetSlotName("Elastic",5);
718 params.SetSignal(medit,1);
719 params.SetSignal(fVertexmode,2);
720 params.SetSignal(fResolution,3);
721 params.SetSignal(fUserctrl,4);
722 params.SetSignal(fElastic,5);
724 // Store projectile and target information in the event structure
730 pproj.SetVector(v,"car");
731 pnucl=pproj.GetNorm();
732 fEvent->SetProjectile(fAproj,fZproj,pnucl);
736 ptarg.SetVector(v,"car");
737 pnucl=ptarg.GetNorm();
738 fEvent->SetTarget(fAtarg,fZtarg,pnucl);
740 params.AddNamedSlot("specmode");
741 params.AddNamedSlot("Specpmin");
742 params.AddNamedSlot("npart");
743 params.AddNamedSlot("ncolpp");
744 params.AddNamedSlot("ncolnp");
745 params.AddNamedSlot("ncolpn");
746 params.AddNamedSlot("ncolnn");
748 params.SetSignal(specmode,"specmode");
749 params.SetSignal(fSpecpmin,"Specpmin");
750 params.SetSignal(npt,"npart");
751 params.SetSignal(ncols[0],"ncolpp");
752 params.SetSignal(ncols[1],"ncolnp");
753 params.SetSignal(ncols[2],"ncolpn");
754 params.SetSignal(ncols[3],"ncolnn");
761 pnucl=sqrt(v[0]*v[0]+v[1]*v[1]+v[2]*v[2]);
763 fEvent->SetProjectile(0,0,pnucl,kf);
767 pnucl=sqrt(v[0]*v[0]+v[1]*v[1]+v[2]*v[2]);
769 fEvent->SetTarget(0,0,pnucl,kf);
772 fEvent->AddDevice(params);
777 if (medit >= 0) Pyedit(medit); // Define which particles are to be kept
779 if (mlist>=0 && select)
786 for (Int_t jpart=1; jpart<=npart; jpart++)
789 charge=Pychge(kf)/3.;
793 // 3-momentum in GeV/c
797 p.SetVector(v,"car");
799 // Production location in meter.
800 v[0]=GetV(jpart,1)/1000.;
801 v[1]=GetV(jpart,2)/1000.;
802 v[2]=GetV(jpart,3)/1000.;
803 r.SetPosition(v,"car");
809 t.SetParticleCode(kf);
810 t.SetName(name.Data());
818 // Build the vertex structures if requested
821 // Check if track belongs within the resolution to an existing vertex
823 for (Int_t jv=1; jv<=fEvent->GetNvertices(); jv++)
825 AliVertex* vx=fEvent->GetVertex(jv);
828 rx=vx->GetPosition();
829 dist=rx.GetDistance(r);
830 if (dist < fResolution)
832 AliTrack* tx=fEvent->GetIdTrack(ntk);
840 if (add) break; // No need to look further for vertex candidates
843 // If track was not close enough to an existing vertex
844 // a new secondary vertex is created
845 if (!add && fVertexmode>1)
847 AliTrack* tx=fEvent->GetIdTrack(ntk);
853 r.SetPositionErrors(v,"car");
855 vert.SetTrackCopy(0);
856 vert.SetVertexCopy(0);
857 vert.SetId((fEvent->GetNvertices())+1);
860 fEvent->AddVertex(vert,0);
864 } // End of loop over the produced particles for each collision
865 } // End of loop over number of collisions for each type
866 } // End of loop over collision types
868 // Link sec. vertices to primary if requested
869 // Note that also the connecting tracks are automatically created
872 AliVertex* vp=fEvent->GetIdVertex(1); // Primary vertex
875 for (Int_t i=2; i<=fEvent->GetNvertices(); i++)
877 AliVertex* vx=fEvent->GetVertex(i);
880 if (vx->GetId() != 1) vp->AddVertex(vx);
886 // Include the spectator tracks in the event structure.
887 if (fNucl && specmode)
892 r.SetPosition(v,"car");
894 zp=fZproj-(ncols[0]+ncols[2]);
896 ap=fAproj-(ncols[0]+ncols[1]+ncols[2]+ncols[3]);
898 zt=fZtarg-(ncols[0]+ncols[1]);
900 at=fAtarg-(ncols[0]+ncols[1]+ncols[2]+ncols[3]);
905 if (pproj.GetNorm() > fSpecpmin)
907 kf=2212; // Projectile spectator protons
908 charge=Pychge(kf)/3.;
909 mass=GetPMAS(Pycomp(kf),1);
911 for (Int_t iprojp=1; iprojp<=zp; iprojp++)
916 t.SetParticleCode(kf);
917 t.SetName(name.Data());
918 t.SetTitle("Projectile spectator proton");
921 t.Set3Momentum(pproj);
927 kf=2112; // Projectile spectator neutrons
928 charge=Pychge(kf)/3.;
929 mass=GetPMAS(Pycomp(kf),1);
931 for (Int_t iprojn=1; iprojn<=(ap-zp); iprojn++)
936 t.SetParticleCode(kf);
937 t.SetName(name.Data());
938 t.SetTitle("Projectile spectator neutron");
941 t.Set3Momentum(pproj);
948 if (ptarg.GetNorm() > fSpecpmin)
950 kf=2212; // Target spectator protons
951 charge=Pychge(kf)/3.;
952 mass=GetPMAS(Pycomp(kf),1);
954 for (Int_t itargp=1; itargp<=zt; itargp++)
959 t.SetParticleCode(kf);
960 t.SetName(name.Data());
961 t.SetTitle("Target spectator proton");
964 t.Set3Momentum(ptarg);
970 kf=2112; // Target spectator neutrons
971 charge=Pychge(kf)/3.;
972 mass=GetPMAS(Pycomp(kf),1);
974 for (Int_t itargn=1; itargn<=(at-zt); itargn++)
979 t.SetParticleCode(kf);
980 t.SetName(name.Data());
981 t.SetTitle("Target spectator neutron");
984 t.Set3Momentum(ptarg);
991 // Link the spectator tracks to the primary vertex.
994 AliVertex* vp=fEvent->GetIdVertex(1);
997 for (Int_t ispec=1; ispec<=nspec; ispec++)
999 AliTrack* tx=fEvent->GetIdTrack(-ispec);
1000 if (tx) vp->AddTrack(tx);
1006 if (!(fEventnum%fPrintfreq) && (mlist || fEvent))
1010 cout << " Number of tracks in the event structure : "
1011 << fEvent->GetNtracks() << endl;
1013 cout << endl; // Create empty output line after the event
1016 if (fOutTree && fSelect) fOutTree->Fill();
1018 ///////////////////////////////////////////////////////////////////////////
1019 AliEvent* AliCollider::GetEvent(Int_t select) const
1021 // Provide pointer to the generated event structure.
1023 // select = 0 : Always return the pointer to the generated event.
1024 // 1 : Only return the pointer to the generated event in case
1025 // the event passed the selection criteria as specified via
1026 // SelectEvent(). Otherwise the value 0 will be returned.
1028 // By invoking GetEvent() the default of select=0 will be used.
1030 if (!select || fSelect)
1039 ///////////////////////////////////////////////////////////////////////////
1040 void AliCollider::EndRun()
1042 // Properly close the output file (if needed).
1047 cout << " *AliCollider::EndRun* Output file correctly closed." << endl;
1050 ///////////////////////////////////////////////////////////////////////////
1051 void AliCollider::SetStable(Int_t id,Int_t mode)
1053 // Declare whether a particle must be regarded as stable or not.
1054 // The parameter "id" indicates the Pythia KF particle code, which
1055 // basically is the PDG particle identifier code.
1056 // The parameter "mode" indicates the action to be taken.
1058 // mode = 0 : Particle will be able to decay
1059 // 1 : Particle will be regarded as stable.
1061 // In case the user does NOT explicitly invoke this function, the standard
1062 // Pythia settings for the decay tables are used.
1064 // When this function is invoked without the "mode" argument, then the
1065 // default of mode=1 will be used for the specified particle.
1069 // 1) This function should be invoked after the initialisation call
1070 // to AliCollider::Init.
1071 // 2) Due to the internals of Pythia, there is no need to specify particles
1072 // and their corresponding anti-particles separately as (un)stable.
1073 // Once a particle has been declared (un)stable, the corresponding
1074 // anti-particle will be treated in the same way.
1076 if (mode==0 || mode==1)
1078 Int_t kc=Pycomp(id);
1082 SetMDCY(kc,1,decay);
1086 cout << " *AliCollider::SetStable* Unknown particle code. id = " << id << endl;
1091 cout << " *AliCollider::SetStable* Invalid parameter. mode = " << mode << endl;
1094 ///////////////////////////////////////////////////////////////////////////
1095 void AliCollider::SelectEvent(Int_t id)
1097 // Add a particle to the event selection list.
1098 // The parameter "id" indicates the Pythia KF particle code, which
1099 // basically is the PDG particle identifier code.
1100 // In case the user has built a selection list via this procedure, only the
1101 // events in which one of the particles specified in the list was generated
1103 // The investigation of the generated particles takes place when the complete
1104 // event is in memory, including all (shortlived) mother particles and resonances.
1105 // So, the settings of the various particle decay modes have no influence on
1106 // the event selection described here.
1108 // If no list has been specified, all events will be accepted.
1110 // Note : id=0 will delete the selection list.
1112 // Be aware of the fact that severe selection criteria (i.e. selecting only
1113 // rare events) may result in long runtimes before an event sample has been
1126 Int_t kc=Pycomp(id);
1129 fSelections=new TArrayI(1);
1130 fSelections->AddAt(kc,0);
1135 Int_t size=fSelections->GetSize();
1136 for (Int_t i=0; i<size; i++)
1138 if (kc==fSelections->At(i))
1147 fSelections->Set(size+1);
1148 fSelections->AddAt(kc,size);
1153 ///////////////////////////////////////////////////////////////////////////
1154 Int_t AliCollider::GetSelectionFlag() const
1156 // Return the value of the selection flag for the total event.
1157 // When the event passed the selection criteria as specified via
1158 // SelectEvent() the value 1 is returned, otherwise the value 0 is returned.
1161 ///////////////////////////////////////////////////////////////////////////
1162 Int_t AliCollider::IsSelected()
1164 // Check whether the generated (sub)event contains one of the particles
1165 // specified in the selection list via SelectEvent().
1166 // If this is the case or when no selection list is present, the value 1
1167 // will be returned, indicating the event is selected to be kept.
1168 // Otherwise the value 0 will be returned.
1170 if (!fSelections) return 1;
1172 Int_t nsel=fSelections->GetSize();
1177 for (Int_t jpart=1; jpart<=npart; jpart++)
1181 for (Int_t i=0; i<nsel; i++)
1183 if (kc==fSelections->At(i))
1193 ///////////////////////////////////////////////////////////////////////////
1194 void AliCollider::SetSpectatorPmin(Float_t pmin)
1196 // Set minimal momentum in GeV/c for spectator tracks to be stored.
1197 // Spectator tracks with a momentum below this threshold will not be stored
1198 // in the (output) event structure.
1199 // This facility allows to minimise the output file size.
1200 // Note that when the user wants to boost the event into another reference
1201 // frame these spectator tracks might have got momenta above the threshold.
1202 // However, when the spectator tracks were not stored in the event structure
1203 // in the original frame, there is no way to retreive them anymore.
1206 ///////////////////////////////////////////////////////////////////////////
1207 Float_t AliCollider::GetSpectatorPmin() const
1209 // Provide the minimal spectator momentum in GeV/c.
1212 ///////////////////////////////////////////////////////////////////////////
1213 TString AliCollider::GetPyname(Int_t kf)
1215 // Provide the correctly truncated Pythia particle name for PGD code kf
1217 // The TPythia6::Pyname returned name is copied into a TString and truncated
1218 // at the first blank to prevent funny trailing characters due to incorrect
1219 // stripping of empty characters in TPythia6::Pyname.
1220 // The truncation at the first blank is allowed due to the Pythia convention
1221 // that particle names never contain blanks.
1226 for (Int_t i=1; i<16; i++)
1228 if (name[i]==' ') break;
1229 sname=sname+name[i];
1233 ///////////////////////////////////////////////////////////////////////////