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.11 2004/03/12 09:18:24 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-4); // 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/03/12 09:18:24 $ Utrecht University
112 ///////////////////////////////////////////////////////////////////////////
114 #include "AliCollider.h"
115 #include "Riostream.h"
117 ClassImp(AliCollider) // Class implementation to enable ROOT I/O
119 AliCollider::AliCollider() : TPythia6()
121 // Default constructor.
122 // All variables initialised to default values.
123 fVertexmode=0; // No vertex structure creation
124 fResolution=1e-5; // Standard resolution is 0.1 micron
156 ///////////////////////////////////////////////////////////////////////////
157 AliCollider::~AliCollider()
159 // Default destructor
181 ///////////////////////////////////////////////////////////////////////////
182 void AliCollider::SetOutputFile(TString s)
184 // Create the output file containing all the data in ROOT output format.
190 fOutFile=new TFile(s.Data(),"RECREATE","AliCollider data");
197 fOutTree=new TTree("T","AliCollider event data");
201 fOutTree->Branch("Events","AliEvent",&fEvent,bsize,split);
203 ///////////////////////////////////////////////////////////////////////////
204 void AliCollider::SetVertexMode(Int_t mode)
206 // Set the mode of the vertex structure creation.
208 // By default all generated tracks will only appear in the AliEvent
209 // structure without any primary (and secondary) vertex structure.
210 // The user can build the vertex structure if he/she wants by means
211 // of the beginpoint location of each AliTrack.
213 // However, one can also let AliCollider automatically create
214 // the primary (and secondary) vertex structure(s).
215 // In this case the primary vertex is given Id=1 and all sec. vertices
216 // are given Id's 2,3,4,....
217 // All vertices are created as standalone entities in the AliEvent structure
218 // without any linking between the various vertices.
219 // For this automated process, the user-selected resolution
220 // (see SetResolution) is used to decide whether or not certain vertex
221 // locations can be resolved.
222 // In case no vertex creation is selected (i.e. the default mode=0),
223 // the value of the resolution is totally irrelevant.
225 // The user can also let AliCollider automatically connect the sec. vertices
226 // to the primary vertex (i.e. mode=3). This process will also automatically
227 // generate the tracks connecting the vertices.
228 // Note that the result of the mode=3 operation may be very sensitive to
229 // the resolution parameter. Therefore, no attempt is made to distinguish
230 // between secondary, tertiary etc... vertices. All sec. vertices are
231 // linked to the primary one.
233 // Irrespective of the selected mode, all generated tracks can be obtained
234 // directly from the AliEvent structure.
235 // In case (sec.) vertex creation is selected, all generated vertices can
236 // also be obtained directly from the AliEvent structure.
237 // These (sec.) vertices contain only the corresponding pointers to the various
238 // tracks which are stored in the AliEvent structure.
240 // Overview of vertex creation modes :
241 // -----------------------------------
242 // mode = 0 ==> No vertex structure will be created
243 // 1 ==> Only primary vertex structure will be created
244 // 2 ==> Unconnected primary and secondary vertices will be created
245 // 3 ==> Primary and secondary vertices will be created where all the
246 // sec. vertices will be connected to the primary vertex.
247 // Also the vertex connecting tracks will be automatically
250 if (mode<0 || mode >3)
252 cout << " *AliCollider::SetVertexMode* Invalid argument mode : " << mode << endl;
260 ///////////////////////////////////////////////////////////////////////////
261 Int_t AliCollider::GetVertexMode() const
263 // Provide the current mode for vertex structure creation.
266 ///////////////////////////////////////////////////////////////////////////
267 void AliCollider::SetResolution(Double_t res)
269 // Set the resolution (in cm) for resolving (sec.) vertices.
270 // By default this resolution is set to 0.1 micron.
271 // Note : In case no vertex creation has been selected, the value of
272 // the resolution is totally irrelevant.
273 fResolution=fabs(res);
275 ///////////////////////////////////////////////////////////////////////////
276 Double_t AliCollider::GetResolution() const
278 // Provide the current resolution (in cm) for resolving (sec.) vertices.
281 ///////////////////////////////////////////////////////////////////////////
282 void AliCollider::SetRunNumber(Int_t run)
284 // Set the user defined run number.
285 // By default the run number is set to 0.
288 ///////////////////////////////////////////////////////////////////////////
289 Int_t AliCollider::GetRunNumber() const
291 // Provide the user defined run number.
294 ///////////////////////////////////////////////////////////////////////////
295 void AliCollider::SetPrintFreq(Int_t n)
297 // Set the print frequency for every 'n' events.
298 // By default the printfrequency is set to 1 (i.e. every event).
301 ///////////////////////////////////////////////////////////////////////////
302 Int_t AliCollider::GetPrintFreq() const
304 // Provide the user selected print frequency.
307 ///////////////////////////////////////////////////////////////////////////
308 void AliCollider::Init(char* frame,char* beam,char* target,Float_t win)
310 // Initialisation of the underlying Pythia generator package.
311 // This routine just invokes TPythia6::Initialize(...) and the arguments
312 // have the corresponding meaning.
313 // The event number is reset to 0.
318 Initialize(frame,beam,target,win);
321 cout << " *AliCollider::Init* Standard Pythia initialisation." << endl;
322 cout << " Beam particle : " << beam << " Target particle : " << target
323 << " Frame = " << frame << " Energy = " << win
326 ///////////////////////////////////////////////////////////////////////////
327 void AliCollider::Init(char* frame,Int_t zp,Int_t ap,Int_t zt,Int_t at,Float_t win)
329 // Initialisation of the underlying Pythia generator package for the generation
330 // of nucleus-nucleus interactions.
331 // In addition to the Pythia standard arguments 'frame' and 'win', the user
332 // can specify here (Z,A) values of the projectile and target nuclei.
334 // Note : The 'win' value denotes either the cms energy per nucleon-nucleon collision
335 // (i.e. frame="cms") or the momentum per nucleon in all other cases.
337 // The event number is reset to 0.
351 if (ap<1 || at<1 || zp>ap || zt>at)
354 cout << " *AliCollider::Init* Invalid input value(s). Zproj = " << zp
355 << " Aproj = " << ap << " Ztarg = " << zt << " Atarg = " << at << endl;
365 cout << " *AliCollider::Init* Nucleus-Nucleus generator initialisation." << endl;
366 cout << " Zproj = " << zp << " Aproj = " << ap << " Ztarg = " << zt << " Atarg = " << at
367 << " Frame = " << frame << " Energy = " << win
370 ///////////////////////////////////////////////////////////////////////////
371 void AliCollider::GetFractions(Float_t zp,Float_t ap,Float_t zt,Float_t at)
373 // Determine the fractions for the various N-N collision processes.
374 // The various processes are : p+p, n+p, p+n and n+n.
385 fFracpp=(zp/ap)*(zt/at);
386 fFracnp=(1.-zp/ap)*(zt/at);
387 fFracpn=(zp/ap)*(1.-zt/at);
388 fFracnn=(1.-zp/ap)*(1.-zt/at);
391 ///////////////////////////////////////////////////////////////////////////
392 void AliCollider::MakeEvent(Int_t npt,Int_t mlist,Int_t medit)
394 // Generate one event.
395 // In case of a nucleus-nucleus interaction, the argument 'npt' denotes
396 // the number of participant nucleons.
397 // Normally also the spectator tracks will be stored into the event structure.
398 // The spectator tracks have a negative user Id to distinguish them from the
399 // ordinary generated tracks.
400 // In case the user has selected the creation of vertex structures, the spectator
401 // tracks will be linked to the primary vertex.
402 // However, specification of npt<0 will suppress the storage of spectator tracks.
403 // In the latter case abs(npt) will be taken as the number of participants.
404 // In case of a standard Pythia run for 'elementary' particle interactions,
405 // the value of npt is totally irrelevant.
407 // The argument 'mlist' denotes the list mode used for Pylist().
408 // Note : mlist<0 suppresses the invokation of Pylist().
409 // By default, no listing is produced (i.e. mlist=-1).
411 // The argument 'medit' denotes the edit mode used for Pyedit().
412 // Note : medit<0 suppresses the invokation of Pyedit().
413 // By default, only 'stable' final particles are kept (i.e. medit=1).
415 // In the case of a standard Pythia run concerning 'elementary' particle
416 // interactions, the projectile and target particle ID's for the created
417 // event structure are set to the corresponding Pythia KF codes.
418 // All the A and Z values are in that case set to zero.
419 // In case of a nucleus-nucleus interaction, the proper A and Z values for
420 // the projectile and target particles are set in the event structure.
421 // However, in this case both particle ID's are set to zero.
423 // Note : Only in case an event passed the selection criteria as specified
424 // via SelectEvent(), the event will appear on the output file.
435 // Counters for the various (proj,targ) combinations : p+p, n+p, p+n and n+n
436 Int_t ncols[4]={0,0,0,0};
446 if (npt<1 || npt>(fAproj+fAtarg))
448 cout << " *AliCollider::MakeEvent* Invalid input value. npt = " << npt
449 << " Aproj = " << fAproj << " Atarg = " << fAtarg << endl;
453 // Determine the number of nucleon-nucleon collisions
455 if (npt%2 && fRan.Uniform()>0.5) ncol+=1;
457 // Determine the number of the various types of N+N interactions
462 Int_t maxa=2; // Indicator whether proj (1) or target (2) has maximal A left
464 Float_t* rans=new Float_t[ncol];
465 fRan.Uniform(rans,ncol);
467 for (Int_t i=0; i<ncol; i++)
469 GetFractions(zp,ap,zt,at);
471 if (rndm<=fFracpp) // p+p interaction
485 if (rndm>fFracpp && rndm<=(fFracpp+fFracnp)) // n+p interaction
498 if (rndm>(fFracpp+fFracnp) && rndm<=(fFracpp+fFracnp+fFracpn)) // p+n interaction
511 if (rndm>(fFracpp+fFracnp+fFracpn)) // n+n interaction
527 if (!(fEventnum%fPrintfreq))
529 cout << " *AliCollider::MakeEvent* Run : " << fRunnum << " Event : " << fEventnum
533 cout << " npart = " << npt << " ncol = " << ncol
534 << " ncolpp = " << ncols[0] << " ncolnp = " << ncols[1]
535 << " ncolpn = " << ncols[2] << " ncolnn = " << ncols[3] << endl;
541 fEvent=new AliEvent();
543 fEvent->SetName(GetName());
544 fEvent->SetTitle(GetTitle());
548 fEvent->SetRunNumber(fRunnum);
549 fEvent->SetEventNumber(fEventnum);
556 Ali3Vector pproj,ptarg;
560 // Make sure the primary vertex gets correct location and Id=1
564 r.SetPosition(v,"car");
568 r.SetPositionErrors(v,"car");
571 vert.SetTrackCopy(0);
572 vert.SetVertexCopy(0);
574 fEvent->AddVertex(vert,0);
578 Float_t charge=0,mass=0;
583 // Singular settings for a normal Pythia elementary particle interation
590 // Generate all the various collisions
591 fSelect=0; // Flag to indicate whether the total event is selected or not
592 Int_t select=0; // Flag to indicate whether the sub-event is selected or not
593 Int_t first=1; // Flag to indicate the first collision process
597 for (Int_t itype=0; itype<ntypes; itype++)
601 if (itype==0 && ncols[itype]) Initialize(fFrame,"p","p",fWin);
602 if (itype==1 && ncols[itype]) Initialize(fFrame,"n","p",fWin);
603 if (itype==2 && ncols[itype]) Initialize(fFrame,"p","n",fWin);
604 if (itype==3 && ncols[itype]) Initialize(fFrame,"n","n",fWin);
606 for (Int_t jcol=0; jcol<ncols[itype]; jcol++)
611 if (select) fSelect=1;
613 if (first) // Store projectile and target information in the event structure
620 pproj.SetVector(v,"car");
621 pnucl=pproj.GetNorm();
622 fEvent->SetProjectile(fAproj,fZproj,pnucl);
626 ptarg.SetVector(v,"car");
627 pnucl=ptarg.GetNorm();
628 fEvent->SetTarget(fAtarg,fZtarg,pnucl);
635 pnucl=sqrt(v[0]*v[0]+v[1]*v[1]+v[2]*v[2]);
637 fEvent->SetProjectile(0,0,pnucl,kf);
641 pnucl=sqrt(v[0]*v[0]+v[1]*v[1]+v[2]*v[2]);
643 fEvent->SetTarget(0,0,pnucl,kf);
648 if (medit >= 0) Pyedit(medit); // Define which particles are to be kept
650 if (mlist>=0 && select) Pylist(mlist);
653 for (Int_t jpart=1; jpart<=npart; jpart++)
656 charge=Pychge(kf)/3.;
660 // 3-momentum in GeV/c
664 p.SetVector(v,"car");
666 // Production location in cm.
667 v[0]=GetV(jpart,1)/10;
668 v[1]=GetV(jpart,2)/10;
669 v[2]=GetV(jpart,3)/10;
670 r.SetPosition(v,"car");
676 t.SetParticleCode(kf);
677 t.SetName(name.Data());
685 // Build the vertex structures if requested
688 // Check if track belongs within the resolution to an existing vertex
690 for (Int_t jv=1; jv<=fEvent->GetNvertices(); jv++)
692 AliVertex* vx=fEvent->GetVertex(jv);
695 rx=vx->GetPosition();
696 dist=rx.GetDistance(r);
697 if (dist < fResolution)
699 AliTrack* tx=fEvent->GetIdTrack(ntk);
707 if (add) break; // No need to look further for vertex candidates
710 // If track was not close enough to an existing vertex
711 // a new secondary vertex is created
712 if (!add && fVertexmode>1)
714 AliTrack* tx=fEvent->GetIdTrack(ntk);
720 r.SetPositionErrors(v,"car");
722 vert.SetTrackCopy(0);
723 vert.SetVertexCopy(0);
724 vert.SetId((fEvent->GetNvertices())+1);
727 fEvent->AddVertex(vert,0);
731 } // End of loop over the produced particles for each collision
732 } // End of loop over number of collisions for each type
733 } // End of loop over collision types
735 // Link sec. vertices to primary if requested
736 // Note that also the connecting tracks are automatically created
739 AliVertex* vp=fEvent->GetIdVertex(1); // Primary vertex
742 for (Int_t i=2; i<=fEvent->GetNvertices(); i++)
744 AliVertex* vx=fEvent->GetVertex(i);
747 if (vx->GetId() != 1) vp->AddVertex(vx);
753 // Include the spectator tracks in the event structure.
754 if (fNucl && specmode)
759 r.SetPosition(v,"car");
761 zp=fZproj-(ncols[0]+ncols[2]);
763 ap=fAproj-(ncols[0]+ncols[1]+ncols[2]+ncols[3]);
765 zt=fZtarg-(ncols[0]+ncols[1]);
767 at=fAtarg-(ncols[0]+ncols[1]+ncols[2]+ncols[3]);
772 if (pproj.GetNorm() > fSpecpmin)
774 kf=2212; // Projectile spectator protons
775 charge=Pychge(kf)/3.;
776 mass=GetPMAS(Pycomp(kf),1);
778 for (Int_t iprojp=1; iprojp<=zp; iprojp++)
783 t.SetParticleCode(kf);
784 t.SetName(name.Data());
785 t.SetTitle("Projectile spectator proton");
788 t.Set3Momentum(pproj);
794 kf=2112; // Projectile spectator neutrons
795 charge=Pychge(kf)/3.;
796 mass=GetPMAS(Pycomp(kf),1);
798 for (Int_t iprojn=1; iprojn<=(ap-zp); iprojn++)
803 t.SetParticleCode(kf);
804 t.SetName(name.Data());
805 t.SetTitle("Projectile spectator neutron");
808 t.Set3Momentum(pproj);
815 if (ptarg.GetNorm() > fSpecpmin)
817 kf=2212; // Target spectator protons
818 charge=Pychge(kf)/3.;
819 mass=GetPMAS(Pycomp(kf),1);
821 for (Int_t itargp=1; itargp<=zt; itargp++)
826 t.SetParticleCode(kf);
827 t.SetName(name.Data());
828 t.SetTitle("Target spectator proton");
831 t.Set3Momentum(ptarg);
837 kf=2112; // Target spectator neutrons
838 charge=Pychge(kf)/3.;
839 mass=GetPMAS(Pycomp(kf),1);
841 for (Int_t itargn=1; itargn<=(at-zt); itargn++)
846 t.SetParticleCode(kf);
847 t.SetName(name.Data());
848 t.SetTitle("Target spectator neutron");
851 t.Set3Momentum(ptarg);
858 // Link the spectator tracks to the primary vertex.
861 AliVertex* vp=fEvent->GetIdVertex(1);
864 for (Int_t ispec=1; ispec<=nspec; ispec++)
866 AliTrack* tx=fEvent->GetIdTrack(-ispec);
867 if (tx) vp->AddTrack(tx);
873 if (mlist && !(fEventnum%fPrintfreq)) cout << endl; // Create empty output line after the event
875 if (fOutTree && fSelect) fOutTree->Fill();
877 ///////////////////////////////////////////////////////////////////////////
878 AliEvent* AliCollider::GetEvent(Int_t select) const
880 // Provide pointer to the generated event structure.
882 // select = 0 : Always return the pointer to the generated event.
883 // 1 : Only return the pointer to the generated event in case
884 // the event passed the selection criteria as specified via
885 // SelectEvent(). Otherwise the value 0 will be returned.
887 // By invoking GetEvent() the default of select=0 will be used.
889 if (!select || fSelect)
898 ///////////////////////////////////////////////////////////////////////////
899 void AliCollider::EndRun()
901 // Properly close the output file (if needed).
906 cout << " *AliCollider::EndRun* Output file correctly closed." << endl;
909 ///////////////////////////////////////////////////////////////////////////
910 void AliCollider::SetStable(Int_t id,Int_t mode)
912 // Declare whether a particle must be regarded as stable or not.
913 // The parameter "id" indicates the Pythia KF particle code, which
914 // basically is the PDG particle identifier code.
915 // The parameter "mode" indicates the action to be taken.
917 // mode = 0 : Particle will be able to decay
918 // 1 : Particle will be regarded as stable.
920 // In case the user does NOT explicitly invoke this function, the standard
921 // Pythia settings for the decay tables are used.
923 // When this function is invoked without the "mode" argument, then the
924 // default of mode=1 will be used for the specified particle.
928 // 1) This function should be invoked after the initialisation call
929 // to AliCollider::Init.
930 // 2) Due to the internals of Pythia, there is no need to specify particles
931 // and their corresponding anti-particles separately as (un)stable.
932 // Once a particle has been declared (un)stable, the corresponding
933 // anti-particle will be treated in the same way.
935 if (mode==0 || mode==1)
945 cout << " *AliCollider::SetStable* Unknown particle code. id = " << id << endl;
950 cout << " *AliCollider::SetStable* Invalid parameter. mode = " << mode << endl;
953 ///////////////////////////////////////////////////////////////////////////
954 void AliCollider::SelectEvent(Int_t id)
956 // Add a particle to the event selection list.
957 // The parameter "id" indicates the Pythia KF particle code, which
958 // basically is the PDG particle identifier code.
959 // In case the user has built a selection list via this procedure, only the
960 // events in which one of the particles specified in the list was generated
962 // The investigation of the generated particles takes place when the complete
963 // event is in memory, including all (shortlived) mother particles and resonances.
964 // So, the settings of the various particle decay modes have no influence on
965 // the event selection described here.
967 // If no list has been specified, all events will be accepted.
969 // Note : id=0 will delete the selection list.
971 // Be aware of the fact that severe selection criteria (i.e. selecting only
972 // rare events) may result in long runtimes before an event sample has been
988 fSelections=new TArrayI(1);
989 fSelections->AddAt(kc,0);
994 Int_t size=fSelections->GetSize();
995 for (Int_t i=0; i<size; i++)
997 if (kc==fSelections->At(i))
1006 fSelections->Set(size+1);
1007 fSelections->AddAt(kc,size);
1012 ///////////////////////////////////////////////////////////////////////////
1013 Int_t AliCollider::GetSelectionFlag() const
1015 // Return the value of the selection flag for the total event.
1016 // When the event passed the selection criteria as specified via
1017 // SelectEvent() the value 1 is returned, otherwise the value 0 is returned.
1020 ///////////////////////////////////////////////////////////////////////////
1021 Int_t AliCollider::IsSelected()
1023 // Check whether the generated (sub)event contains one of the particles
1024 // specified in the selection list via SelectEvent().
1025 // If this is the case or when no selection list is present, the value 1
1026 // will be returned, indicating the event is selected to be kept.
1027 // Otherwise the value 0 will be returned.
1029 if (!fSelections) return 1;
1031 Int_t nsel=fSelections->GetSize();
1036 for (Int_t jpart=1; jpart<=npart; jpart++)
1040 for (Int_t i=0; i<nsel; i++)
1042 if (kc==fSelections->At(i))
1052 ///////////////////////////////////////////////////////////////////////////
1053 void AliCollider::SetSpectatorPmin(Float_t pmin)
1055 // Set minimal momentum in GeV/c for spectator tracks to be stored.
1056 // Spectator tracks with a momentum below this threshold will not be stored
1057 // in the (output) event structure.
1058 // This facility allows to minimise the output file size.
1059 // Note that when the user wants to boost the event into another reference
1060 // frame these spectator tracks might have got momenta above the threshold.
1061 // However, when the spectator tracks were not stored in the event structure
1062 // in the original frame, there is no way to retreive them anymore.
1065 ///////////////////////////////////////////////////////////////////////////
1066 Float_t AliCollider::GetSpectatorPmin() const
1068 // Provide the minimal spectator momentum in GeV/c.
1071 ///////////////////////////////////////////////////////////////////////////
1072 TString AliCollider::GetPyname(Int_t kf)
1074 // Provide the correctly truncated Pythia particle name for PGD code kf
1076 // The TPythia6::Pyname returned name is copied into a TString and truncated
1077 // at the first blank to prevent funny trailing characters due to incorrect
1078 // stripping of empty characters in TPythia6::Pyname.
1079 // The truncation at the first blank is allowed due to the Pythia convention
1080 // that particle names never contain blanks.
1085 for (Int_t i=1; i<16; i++)
1087 if (name[i]==' ') break;
1088 sname=sname+name[i];
1092 ///////////////////////////////////////////////////////////////////////////