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 **************************************************************************/
19 // Generator using the TPythia interface (via AliPythia)
20 // to generate pp collisions.
21 // Using SetNuclei() also nuclear modifications to the structure functions
22 // can be taken into account. This makes, of course, only sense for the
23 // generation of the products of hard processes (heavy flavor, jets ...)
25 // andreas.morsch@cern.ch
28 #include <TDatabasePDG.h>
29 #include <TParticle.h>
34 #include "AliDecayerPythia.h"
35 #include "AliGenPythia.h"
36 #include "AliHeader.h"
37 #include "AliGenPythiaEventHeader.h"
38 #include "AliPythia.h"
39 #include "AliPythiaRndm.h"
42 #include "AliRunLoader.h"
44 #include "pyquenCommon.h"
46 ClassImp(AliGenPythia)
48 AliGenPythia::AliGenPythia()
51 // Default Constructor
57 fInteractionRate = 0.;
61 fDecayer = new AliDecayerPythia();
72 if (!AliPythiaRndm::GetPythiaRandom())
73 AliPythiaRndm::SetPythiaRandom(GetRandom());
76 AliGenPythia::AliGenPythia(Int_t npart)
79 // default charm production at 5. 5 TeV
81 // structure function GRVHO
84 fTitle= "Particle Generator using PYTHIA";
88 fInteractionRate = 0.;
98 fDecayer = new AliDecayerPythia();
99 // Set random number generator
100 if (!AliPythiaRndm::GetPythiaRandom())
101 AliPythiaRndm::SetPythiaRandom(GetRandom());
103 // Produced particles
104 fParticles = new TClonesArray("TParticle",1000);
112 SetJetReconstructionMode();
116 // Options determining what to keep in the stack (Heavy flavour generation)
117 fStackFillOpt = kFlavorSelection; // Keep particle with selected flavor
118 fFeedDownOpt = kTRUE; // allow feed down from higher family
119 // Fragmentation on/off
120 fFragmentation = kTRUE;
121 // Default counting mode
122 fCountMode = kCountAll;
124 SetPycellParameters();
128 AliGenPythia::AliGenPythia(const AliGenPythia & Pythia)
135 AliGenPythia::~AliGenPythia()
138 if(fEventsTime) delete fEventsTime;
141 void AliGenPythia::SetInteractionRate(Float_t rate,Float_t timewindow)
143 // Generate pileup using user specified rate
144 fInteractionRate = rate;
145 fTimeWindow = timewindow;
149 void AliGenPythia::GeneratePileup()
151 // Generate sub events time for pileup
153 if(fInteractionRate == 0.) {
154 Warning("GeneratePileup","Zero interaction specified. Skipping pileup generation.\n");
158 Int_t npart = NumberParticles();
160 Warning("GeneratePileup","Negative number of particles. Skipping pileup generation.\n");
164 if(fEventsTime) delete fEventsTime;
165 fEventsTime = new TArrayF(npart);
166 TArrayF &array = *fEventsTime;
167 for(Int_t ipart = 0; ipart < npart; ipart++)
170 Float_t eventtime = 0.;
173 eventtime += (AliPythiaRndm::GetPythiaRandom())->Exp(1./fInteractionRate);
174 if(eventtime > fTimeWindow) break;
175 array.Set(array.GetSize()+1);
176 array[array.GetSize()-1] = eventtime;
182 eventtime -= (AliPythiaRndm::GetPythiaRandom())->Exp(1./fInteractionRate);
183 if(TMath::Abs(eventtime) > fTimeWindow) break;
184 array.Set(array.GetSize()+1);
185 array[array.GetSize()-1] = eventtime;
188 SetNumberParticles(fEventsTime->GetSize());
191 void AliGenPythia::SetPycellParameters(Float_t etamax, Int_t neta, Int_t nphi,
192 Float_t thresh, Float_t etseed, Float_t minet, Float_t r)
194 // Set pycell parameters
195 fPycellEtaMax = etamax;
198 fPycellThreshold = thresh;
199 fPycellEtSeed = etseed;
200 fPycellMinEtJet = minet;
201 fPycellMaxRadius = r;
206 void AliGenPythia::SetEventListRange(Int_t eventFirst, Int_t eventLast)
208 // Set a range of event numbers, for which a table
209 // of generated particle will be printed
210 fDebugEventFirst = eventFirst;
211 fDebugEventLast = eventLast;
212 if (fDebugEventLast==-1) fDebugEventLast=fDebugEventFirst;
215 void AliGenPythia::Init()
219 SetMC(AliPythia::Instance());
220 fPythia=(AliPythia*) fMCEvGen;
223 fParentWeight=1./Float_t(fNpart);
225 // Forward Paramters to the AliPythia object
226 fDecayer->SetForceDecay(fForceDecay);
230 fPythia->SetCKIN(3,fPtHardMin);
231 fPythia->SetCKIN(4,fPtHardMax);
232 fPythia->SetCKIN(7,fYHardMin);
233 fPythia->SetCKIN(8,fYHardMax);
235 if (fAProjectile > 0 && fATarget > 0) fPythia->SetNuclei(fAProjectile, fATarget);
237 if (fFragmentation) {
238 fPythia->SetMSTP(111,1);
240 fPythia->SetMSTP(111,0);
244 // initial state radiation
245 fPythia->SetMSTP(61,fGinit);
246 // final state radiation
247 fPythia->SetMSTP(71,fGfinal);
250 fPythia->SetMSTP(91,1);
251 fPythia->SetPARP(91,fPtKick);
253 fPythia->SetMSTP(91,0);
258 fRL = AliRunLoader::Open(fFileName, "Partons");
259 fRL->LoadKinematics();
267 fPythia->ProcInit(fProcess,fEnergyCMS,fStrucFunc);
269 // Parent and Children Selection
273 case kPyCharmUnforced:
274 case kPyCharmPbPbMNR:
277 fParentSelect[0] = 411;
278 fParentSelect[1] = 421;
279 fParentSelect[2] = 431;
280 fParentSelect[3] = 4122;
286 fParentSelect[0] = 421;
289 case kPyDPlusPbPbMNR:
292 fParentSelect[0] = 411;
296 case kPyBeautyPbPbMNR:
297 case kPyBeautypPbMNR:
299 fParentSelect[0]= 511;
300 fParentSelect[1]= 521;
301 fParentSelect[2]= 531;
302 fParentSelect[3]= 5122;
303 fParentSelect[4]= 5132;
304 fParentSelect[5]= 5232;
305 fParentSelect[6]= 5332;
308 case kPyBeautyUnforced:
309 fParentSelect[0] = 511;
310 fParentSelect[1] = 521;
311 fParentSelect[2] = 531;
312 fParentSelect[3] = 5122;
313 fParentSelect[4] = 5132;
314 fParentSelect[5] = 5232;
315 fParentSelect[6] = 5332;
320 fParentSelect[0] = 443;
332 // JetFinder for Trigger
334 // Configure detector (EMCAL like)
336 fPythia->SetPARU(51, fPycellEtaMax);
337 fPythia->SetMSTU(51, fPycellNEta);
338 fPythia->SetMSTU(52, fPycellNPhi);
340 // Configure Jet Finder
342 fPythia->SetPARU(58, fPycellThreshold);
343 fPythia->SetPARU(52, fPycellEtSeed);
344 fPythia->SetPARU(53, fPycellMinEtJet);
345 fPythia->SetPARU(54, fPycellMaxRadius);
346 fPythia->SetMSTU(54, 2);
348 // This counts the total number of calls to Pyevnt() per run.
363 Warning("Init","SetNuclei used. Use SetProjectile + SetTarget instead. fDyBoost has been reset to 0\n");
367 fPythia->InitQuenching(0., 0.1, 0.6e6, 0);
372 void AliGenPythia::Generate()
374 // Generate one event
376 fDecayer->ForceDecay();
378 Float_t polar[3] = {0,0,0};
379 Float_t origin[3] = {0,0,0};
381 // converts from mm/c to s
382 const Float_t kconv=0.001/2.999792458e8;
392 // Set collision vertex position
393 if (fVertexSmear == kPerEvent) Vertex();
402 // Switch hadronisation off
404 fPythia->SetMSTJ(1, 0);
406 // Either produce new event or read partons from file
408 if (!fReadFromFile) {
410 fNpartons = fPythia->GetN();
412 printf("Loading Event %d\n",AliRunLoader::GetRunLoader()->GetEventNumber());
413 fRL->GetEvent(AliRunLoader::GetRunLoader()->GetEventNumber());
415 LoadEvent(fRL->Stack(), 0 , 1);
420 // Run quenching routine
424 } else if (fQuench == 2){
425 fPythia->Pyquen(208., 0, 0.);
428 // Switch hadronisation on
430 fPythia->SetMSTJ(1, 1);
432 // .. and perform hadronisation
433 // printf("Calling hadronisation %d\n", fPythia->GetN());
436 fPythia->ImportParticles(fParticles,"All");
444 Int_t np = fParticles->GetEntriesFast();
446 if (np == 0) continue;
450 Int_t* pParent = new Int_t[np];
451 Int_t* pSelected = new Int_t[np];
452 Int_t* trackIt = new Int_t[np];
453 for (i = 0; i < np; i++) {
459 Int_t nc = 0; // Total n. of selected particles
460 Int_t nParents = 0; // Selected parents
461 Int_t nTkbles = 0; // Trackable particles
462 if (fProcess != kPyMb && fProcess != kPyJets &&
463 fProcess != kPyDirectGamma &&
464 fProcess != kPyMbNonDiffr &&
467 for (i = 0; i < np; i++) {
468 TParticle* iparticle = (TParticle *) fParticles->At(i);
469 Int_t ks = iparticle->GetStatusCode();
470 kf = CheckPDGCode(iparticle->GetPdgCode());
471 // No initial state partons
472 if (ks==21) continue;
474 // Heavy Flavor Selection
481 if (kfl > 100000) kfl %= 100000;
482 if (kfl > 10000) kfl %= 10000;
484 if (kfl > 10) kfl/=100;
486 if (kfl > 10) kfl/=10;
487 Int_t ipa = iparticle->GetFirstMother()-1;
490 // Establish mother daughter relation between heavy quarks and mesons
492 if (kf >= fFlavorSelect && kf <= 6) {
493 Int_t idau = iparticle->GetFirstDaughter() - 1;
495 TParticle* daughter = (TParticle *) fParticles->At(idau);
496 Int_t pdgD = daughter->GetPdgCode();
497 if (pdgD == 91 || pdgD == 92) {
498 Int_t jmin = daughter->GetFirstDaughter() - 1;
499 Int_t jmax = daughter->GetLastDaughter() - 1;
500 for (Int_t j = jmin; j <= jmax; j++)
501 ((TParticle *) fParticles->At(j))->SetFirstMother(i+1);
502 } // is string or cluster
508 TParticle * mother = (TParticle *) fParticles->At(ipa);
509 kfMo = TMath::Abs(mother->GetPdgCode());
512 // What to keep in Stack?
513 Bool_t flavorOK = kFALSE;
514 Bool_t selectOK = kFALSE;
516 if (kfl >= fFlavorSelect) flavorOK = kTRUE;
518 if (kfl > fFlavorSelect) {
522 if (kfl == fFlavorSelect) flavorOK = kTRUE;
524 switch (fStackFillOpt) {
525 case kFlavorSelection:
528 case kParentSelection:
529 if (ParentSelected(kf) || kf <= 10) selectOK = kTRUE;
532 if (flavorOK && selectOK) {
534 // Heavy flavor hadron or quark
536 // Kinematic seletion on final state heavy flavor mesons
537 if (ParentSelected(kf) && !KinematicSelection(iparticle, 0))
542 if (ParentSelected(kf)) ++nParents; // Update parent count
543 // printf("\n particle (HF) %d %d %d", i, pSelected[i], kf);
545 // Kinematic seletion on decay products
546 if (fCutOnChild && ParentSelected(kfMo) && ChildSelected(kf)
547 && !KinematicSelection(iparticle, 1))
553 // Select if mother was selected and is not tracked
555 if (pSelected[ipa] &&
556 !trackIt[ipa] && // mother will be tracked ?
557 kfMo != 5 && // mother is b-quark, don't store fragments
558 kfMo != 4 && // mother is c-quark, don't store fragments
559 kf != 92) // don't store string
562 // Semi-stable or de-selected: diselect decay products:
565 if (pSelected[i] == -1 || fDecayer->GetLifetime(kf) > fMaxLifeTime)
567 Int_t ipF = iparticle->GetFirstDaughter();
568 Int_t ipL = iparticle->GetLastDaughter();
569 if (ipF > 0) for (j = ipF-1; j < ipL; j++) pSelected[j] = -1;
571 // printf("\n particle (decay) %d %d %d", i, pSelected[i], kf);
572 pSelected[i] = (pSelected[i] == -1) ? 0 : 1;
575 if (pSelected[i] == -1) pSelected[i] = 0;
576 if (!pSelected[i]) continue;
577 // Count quarks only if you did not include fragmentation
578 if (fFragmentation && kf <= 10) continue;
581 // Decision on tracking
584 // Track final state particle
585 if (ks == 1) trackIt[i] = 1;
586 // Track semi-stable particles
587 if ((ks == 1) || (fDecayer->GetLifetime(kf) > fMaxLifeTime)) trackIt[i] = 1;
588 // Track particles selected by process if undecayed.
589 if (fForceDecay == kNoDecay) {
590 if (ParentSelected(kf)) trackIt[i] = 1;
592 if (ParentSelected(kf)) trackIt[i] = 0;
594 if (trackIt[i] == 1) ++nTkbles; // Update trackable counter
598 } // particle selection loop
600 for (i = 0; i<np; i++) {
601 if (!pSelected[i]) continue;
602 TParticle * iparticle = (TParticle *) fParticles->At(i);
603 kf = CheckPDGCode(iparticle->GetPdgCode());
604 Int_t ks = iparticle->GetStatusCode();
605 p[0] = iparticle->Px();
606 p[1] = iparticle->Py();
607 p[2] = iparticle->Pz();
608 p[3] = iparticle->Energy();
610 origin[0] = fVertex[0]+iparticle->Vx()/10; // [cm]
611 origin[1] = fVertex[1]+iparticle->Vy()/10; // [cm]
612 origin[2] = fVertex[2]+iparticle->Vz()/10; // [cm]
614 Float_t tof = kconv*iparticle->T();
615 Int_t ipa = iparticle->GetFirstMother()-1;
616 Int_t iparent = (ipa > -1) ? pParent[ipa] : -1;
618 PushTrack(fTrackIt*trackIt[i], iparent, kf,
619 p[0], p[1], p[2], p[3],
620 origin[0], origin[1], origin[2], tof,
621 polar[0], polar[1], polar[2],
622 kPPrimary, nt, 1., ks);
633 if (pParent) delete[] pParent;
634 if (pSelected) delete[] pSelected;
635 if (trackIt) delete[] trackIt;
638 switch (fCountMode) {
640 // printf(" Count all \n");
644 // printf(" Count parents \n");
647 case kCountTrackables:
648 // printf(" Count trackable \n");
652 if (jev >= fNpart || fNpart == -1) {
653 fKineBias=Float_t(fNpart)/Float_t(fTrials);
654 printf("\n Trials: %i %i %i\n",fTrials, fNpart, jev);
656 fQ += fPythia->GetVINT(51);
657 fX1 += fPythia->GetVINT(41);
658 fX2 += fPythia->GetVINT(42);
659 fTrialsRun += fTrials;
666 SetHighWaterMark(nt);
667 // adjust weight due to kinematic selection
670 fXsection=fPythia->GetPARI(1);
673 Int_t AliGenPythia::GenerateMB()
676 // Min Bias selection and other global selections
678 Int_t i, kf, nt, iparent;
681 Float_t polar[3] = {0,0,0};
682 Float_t origin[3] = {0,0,0};
683 // converts from mm/c to s
684 const Float_t kconv=0.001/2.999792458e8;
688 Int_t np = (fHadronisation) ? fParticles->GetEntriesFast() : fNpartons;
691 Int_t* pParent = new Int_t[np];
692 for (i=0; i< np; i++) pParent[i] = -1;
693 if (fProcess == kPyJets || fProcess == kPyDirectGamma) {
694 TParticle* jet1 = (TParticle *) fParticles->At(6);
695 TParticle* jet2 = (TParticle *) fParticles->At(7);
696 if (!CheckTrigger(jet1, jet2)) return 0;
700 //Introducing child cuts in case kPyW
701 if ( (fProcess == kPyW) && (fCutOnChild == 1) ) {
702 if ( !CheckKinematicsOnChild() ) return 0;
706 for (i = 0; i < np; i++) {
708 TParticle * iparticle = (TParticle *) fParticles->At(i);
709 kf = CheckPDGCode(iparticle->GetPdgCode());
710 Int_t ks = iparticle->GetStatusCode();
711 Int_t km = iparticle->GetFirstMother();
712 if ((ks == 1 && kf!=0 && KinematicSelection(iparticle, 0)) ||
714 (fProcess == kPyJets && ks == 21 && km == 0 && i>1)) {
716 if (ks == 1) trackIt = 1;
717 Int_t ipa = iparticle->GetFirstMother()-1;
719 iparent = (ipa > -1) ? pParent[ipa] : -1;
722 // store track information
723 p[0] = iparticle->Px();
724 p[1] = iparticle->Py();
725 p[2] = iparticle->Pz();
726 p[3] = iparticle->Energy();
729 origin[0] = fVertex[0]+iparticle->Vx()/10; // [cm]
730 origin[1] = fVertex[1]+iparticle->Vy()/10; // [cm]
731 origin[2] = fVertex[2]+iparticle->Vz()/10; // [cm]
733 Float_t tof = fEventTime + kconv * iparticle->T();
735 PushTrack(fTrackIt*trackIt, iparent, kf,
736 p[0], p[1], p[2], p[3],
737 origin[0], origin[1], origin[2], tof,
738 polar[0], polar[1], polar[2],
739 kPPrimary, nt, 1., ks);
741 // Special Treatment to store color-flow
743 if (ks == 3 || ks == 13 || ks == 14) {
744 TParticle* particle = 0;
746 particle = fStack->Particle(nt);
748 particle = gAlice->Stack()->Particle(nt);
750 particle->SetFirstDaughter(fPythia->GetK(2, i));
751 particle->SetLastDaughter(fPythia->GetK(3, i));
756 SetHighWaterMark(nt);
761 if (pParent) delete[] pParent;
763 printf("\n I've put %i particles on the stack \n",nc);
768 void AliGenPythia::FinishRun()
770 // Print x-section summary
779 printf("\nTotal number of Pyevnt() calls %d\n", fTrialsRun);
780 printf("\nMean Q, x1, x2: %f %f %f\n", fQ, fX1, fX2);
783 void AliGenPythia::AdjustWeights()
785 // Adjust the weights after generation of all events
789 Int_t ntrack=gAlice->GetMCApp()->GetNtrack();
790 for (Int_t i=0; i<ntrack; i++) {
791 part= gAlice->GetMCApp()->Particle(i);
792 part->SetWeight(part->GetWeight()*fKineBias);
797 void AliGenPythia::SetNuclei(Int_t a1, Int_t a2)
799 // Treat protons as inside nuclei with mass numbers a1 and a2
807 void AliGenPythia::MakeHeader()
810 if (gAlice->GetEvNumber()>=fDebugEventFirst &&
811 gAlice->GetEvNumber()<=fDebugEventLast) fPythia->Pylist(2);
814 // Builds the event header, to be called after each event
815 if (fHeader) delete fHeader;
816 fHeader = new AliGenPythiaEventHeader("Pythia");
819 ((AliGenPythiaEventHeader*) fHeader)->SetProcessType(fPythia->GetMSTI(1));
822 ((AliGenPythiaEventHeader*) fHeader)->SetTrials(fTrials);
825 fHeader->SetPrimaryVertex(fVertex);
827 // Jets that have triggered
829 if (fProcess == kPyJets)
833 GetJets(njet, ntrig, jets);
836 for (Int_t i = 0; i < ntrig; i++) {
837 ((AliGenPythiaEventHeader*) fHeader)->AddJet(jets[0][i], jets[1][i], jets[2][i],
842 // Copy relevant information from external header, if present.
847 AliGenPythiaEventHeader* exHeader = (AliGenPythiaEventHeader*) (fRL->GetHeader()->GenEventHeader());
848 for (Int_t i = 0; i < exHeader->NTriggerJets(); i++)
850 printf("Adding Jet %d %d \n", i, exHeader->NTriggerJets());
853 exHeader->TriggerJet(i, uqJet);
854 ((AliGenPythiaEventHeader*) fHeader)->AddUQJet(uqJet[0], uqJet[1], uqJet[2], uqJet[3]);
858 // Store quenching parameters
865 fPythia->GetQuenchingParameters(xp, yp, z);
868 Double_t r1 = PARIMP.rb1;
869 Double_t r2 = PARIMP.rb2;
870 Double_t b = PARIMP.b1;
871 Double_t r = 0.5 * TMath::Sqrt(2. * (r1 * r1 + r2 * r2) - b * b);
872 Double_t phi = PARIMP.psib1;
873 xp = r * TMath::Cos(phi);
874 yp = r * TMath::Sin(phi);
877 ((AliGenPythiaEventHeader*) fHeader)->SetXYJet(xp, yp);
878 ((AliGenPythiaEventHeader*) fHeader)->SetZQuench(z);
887 void AliGenPythia::AddHeader(AliGenEventHeader* header)
889 // Add header to container or runloader
891 fContainer->AddHeader(header);
893 AliRunLoader::GetRunLoader()->GetHeader()->SetGenEventHeader(header);
898 Bool_t AliGenPythia::CheckTrigger(TParticle* jet1, TParticle* jet2)
900 // Check the kinematic trigger condition
903 eta[0] = jet1->Eta();
904 eta[1] = jet2->Eta();
906 phi[0] = jet1->Phi();
907 phi[1] = jet2->Phi();
909 pdg[0] = jet1->GetPdgCode();
910 pdg[1] = jet2->GetPdgCode();
911 Bool_t triggered = kFALSE;
913 if (fProcess == kPyJets) {
918 // Use Pythia clustering on parton level to determine jet axis
920 GetJets(njets, ntrig, jets);
922 if (ntrig) triggered = kTRUE;
927 if (pdg[0] == kGamma) {
931 //Check eta range first...
932 if ((eta[ij] < fEtaMaxJet && eta[ij] > fEtaMinJet) &&
933 (eta[ig] < fEtaMaxGamma && eta[ig] > fEtaMinGamma))
935 //Eta is okay, now check phi range
936 if ((phi[ij] < fPhiMaxJet && phi[ij] > fPhiMinJet) &&
937 (phi[ig] < fPhiMaxGamma && phi[ig] > fPhiMinGamma))
947 //Checking Kinematics on Child (status code 1, particle code ?, kin cuts
948 Bool_t AliGenPythia::CheckKinematicsOnChild(){
950 Bool_t checking = kFALSE;
951 Int_t j, kcode, ks, km;
952 Int_t nPartAcc = 0; //number of particles in the acceptance range
953 Int_t numberOfAcceptedParticles = 1;
954 if (fNumberOfAcceptedParticles != 0) { numberOfAcceptedParticles = fNumberOfAcceptedParticles; }
955 Int_t npart = fParticles->GetEntriesFast();
957 for (j = 0; j<npart; j++) {
958 TParticle * jparticle = (TParticle *) fParticles->At(j);
959 kcode = TMath::Abs( CheckPDGCode(jparticle->GetPdgCode()) );
960 ks = jparticle->GetStatusCode();
961 km = jparticle->GetFirstMother();
963 if( (ks == 1) && (kcode == fPdgCodeParticleforAcceptanceCut) && (KinematicSelection(jparticle,1)) ){
968 if( numberOfAcceptedParticles <= nPartAcc){
977 AliGenPythia& AliGenPythia::operator=(const AliGenPythia& rhs)
979 // Assignment operator
984 void AliGenPythia::LoadEvent(AliStack* stack, Int_t flag, Int_t reHadr)
987 // Load event into Pythia Common Block
990 Int_t npart = stack -> GetNprimary();
994 (fPythia->GetPyjets())->N = npart;
996 n0 = (fPythia->GetPyjets())->N;
997 (fPythia->GetPyjets())->N = n0 + npart;
1001 for (Int_t part = 0; part < npart; part++) {
1002 TParticle *MPart = stack->Particle(part);
1004 Int_t kf = MPart->GetPdgCode();
1005 Int_t ks = MPart->GetStatusCode();
1006 Int_t idf = MPart->GetFirstDaughter();
1007 Int_t idl = MPart->GetLastDaughter();
1010 if (ks == 11 || ks == 12) {
1017 Float_t px = MPart->Px();
1018 Float_t py = MPart->Py();
1019 Float_t pz = MPart->Pz();
1020 Float_t e = MPart->Energy();
1021 Float_t m = MPart->GetCalcMass();
1024 (fPythia->GetPyjets())->P[0][part+n0] = px;
1025 (fPythia->GetPyjets())->P[1][part+n0] = py;
1026 (fPythia->GetPyjets())->P[2][part+n0] = pz;
1027 (fPythia->GetPyjets())->P[3][part+n0] = e;
1028 (fPythia->GetPyjets())->P[4][part+n0] = m;
1030 (fPythia->GetPyjets())->K[1][part+n0] = kf;
1031 (fPythia->GetPyjets())->K[0][part+n0] = ks;
1032 (fPythia->GetPyjets())->K[3][part+n0] = idf + 1;
1033 (fPythia->GetPyjets())->K[4][part+n0] = idl + 1;
1034 (fPythia->GetPyjets())->K[2][part+n0] = MPart->GetFirstMother() + 1;
1039 void AliGenPythia::RecJetsUA1(Int_t& njets, Float_t jets [4][50])
1042 // Calls the Pythia jet finding algorithm to find jets in the current event
1047 Int_t n = fPythia->GetN();
1051 fPythia->Pycell(njets);
1053 for (i = 0; i < njets; i++) {
1054 Float_t px = (fPythia->GetPyjets())->P[0][n+i];
1055 Float_t py = (fPythia->GetPyjets())->P[1][n+i];
1056 Float_t pz = (fPythia->GetPyjets())->P[2][n+i];
1057 Float_t e = (fPythia->GetPyjets())->P[3][n+i];
1068 void AliGenPythia::GetJets(Int_t& nJets, Int_t& nJetsTrig, Float_t jets[4][10])
1071 // Calls the Pythia clustering algorithm to find jets in the current event
1073 Int_t n = fPythia->GetN();
1076 if (fJetReconstruction == kCluster) {
1078 // Configure cluster algorithm
1080 fPythia->SetPARU(43, 2.);
1081 fPythia->SetMSTU(41, 1);
1083 // Call cluster algorithm
1085 fPythia->Pyclus(nJets);
1087 // Loading jets from common block
1093 fPythia->Pycell(nJets);
1097 for (i = 0; i < nJets; i++) {
1098 Float_t px = (fPythia->GetPyjets())->P[0][n+i];
1099 Float_t py = (fPythia->GetPyjets())->P[1][n+i];
1100 Float_t pz = (fPythia->GetPyjets())->P[2][n+i];
1101 Float_t e = (fPythia->GetPyjets())->P[3][n+i];
1102 Float_t pt = TMath::Sqrt(px * px + py * py);
1103 Float_t phi = TMath::Pi() + TMath::ATan2(-py, -px);
1104 Float_t theta = TMath::ATan2(pt,pz);
1105 Float_t et = e * TMath::Sin(theta);
1106 Float_t eta = -TMath::Log(TMath::Tan(theta / 2.));
1108 eta > fEtaMinJet && eta < fEtaMaxJet &&
1109 phi > fPhiMinJet && phi < fPhiMaxJet &&
1110 et > fEtMinJet && et < fEtMaxJet
1113 jets[0][nJetsTrig] = px;
1114 jets[1][nJetsTrig] = py;
1115 jets[2][nJetsTrig] = pz;
1116 jets[3][nJetsTrig] = e;
1118 // printf("\n........-Jet #%d: %10.3f %10.3f %10.3f %10.3f \n", i, pt, et, eta, phi * kRaddeg);
1120 // printf("\n........-Jet #%d: %10.3f %10.3f %10.3f %10.3f \n", i, pt, et, eta, phi * kRaddeg);
1125 void AliGenPythia::GetSubEventTime()
1127 // Calculates time of the next subevent
1130 TArrayF &array = *fEventsTime;
1131 fEventTime = array[fCurSubEvent++];
1133 // printf(" Event time: %d %f %p",fCurSubEvent,fEventTime,fEventsTime);
1138 void AliGenPythia::Streamer(TBuffer &R__b)
1140 // Stream an object of class AliGenPythia.
1142 if (R__b.IsReading()) {
1143 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
1144 AliGenerator::Streamer(R__b);
1145 R__b >> (Int_t&)fProcess;
1146 R__b >> (Int_t&)fStrucFunc;
1147 R__b >> (Int_t&)fForceDecay;
1151 fParentSelect.Streamer(R__b);
1152 fChildSelect.Streamer(R__b);
1154 // (AliPythia::Instance())->Streamer(R__b);
1157 // if (fDecayer) fDecayer->Streamer(R__b);
1159 R__b.WriteVersion(AliGenPythia::IsA());
1160 AliGenerator::Streamer(R__b);
1161 R__b << (Int_t)fProcess;
1162 R__b << (Int_t)fStrucFunc;
1163 R__b << (Int_t)fForceDecay;
1167 fParentSelect.Streamer(R__b);
1168 fChildSelect.Streamer(R__b);
1173 // fDecayer->Streamer(R__b);