/************************************************************************** * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. * * * * Author: The ALICE Off-line Project. * * Contributors are mentioned in the code where appropriate. * * * * Permission to use, copy, modify and distribute this software and its * * documentation strictly for non-commercial purposes is hereby granted * * without fee, provided that the above copyright notice appears in all * * copies and that both the copyright notice and this permission notice * * appear in the supporting documentation. The authors make no claims * * about the suitability of this software for any purpose. It is * * provided "as is" without express or implied warranty. * **************************************************************************/ /* $Id$ */ // // Generator using the TPythia interface (via AliPythia) // to generate pp collisions. // Using SetNuclei() also nuclear modifications to the structure functions // can be taken into account. This makes, of course, only sense for the // generation of the products of hard processes (heavy flavor, jets ...) // // andreas.morsch@cern.ch // #include #include #include #include #include #include #include #include "AliConst.h" #include "AliDecayerPythia.h" #include "AliGenPythia.h" #include "AliFastGlauber.h" #include "AliHeader.h" #include "AliGenPythiaEventHeader.h" #include "AliPythia.h" #include "AliPythiaRndm.h" #include "AliRun.h" #include "AliStack.h" #include "AliRunLoader.h" #include "AliMC.h" #include "PyquenCommon.h" ClassImp(AliGenPythia) AliGenPythia::AliGenPythia(): AliGenMC(), fProcess(kPyCharm), fItune(-1), fStrucFunc(kCTEQ5L), fKineBias(0.), fTrials(0), fTrialsRun(0), fQ(0.), fX1(0.), fX2(0.), fEventTime(0.), fInteractionRate(0.), fTimeWindow(0.), fCurSubEvent(0), fEventsTime(0), fNev(0), fFlavorSelect(0), fXsection(0.), fPythia(0), fPtHardMin(0.), fPtHardMax(1.e4), fYHardMin(-1.e10), fYHardMax(1.e10), fGinit(1), fGfinal(1), fHadronisation(1), fNpartons(0), fReadFromFile(0), fQuench(0), fQhat(0.), fLength(0.), fImpact(0.), fPtKick(1.), fFullEvent(kTRUE), fDecayer(new AliDecayerPythia()), fDebugEventFirst(-1), fDebugEventLast(-1), fEtMinJet(0.), fEtMaxJet(1.e4), fEtaMinJet(-20.), fEtaMaxJet(20.), fPhiMinJet(0.), fPhiMaxJet(2.* TMath::Pi()), fJetReconstruction(kCell), fEtaMinGamma(-20.), fEtaMaxGamma(20.), fPhiMinGamma(0.), fPhiMaxGamma(2. * TMath::Pi()), fPycellEtaMax(2.), fPycellNEta(274), fPycellNPhi(432), fPycellThreshold(0.), fPycellEtSeed(4.), fPycellMinEtJet(10.), fPycellMaxRadius(1.), fStackFillOpt(kFlavorSelection), fFeedDownOpt(kTRUE), fFragmentation(kTRUE), fSetNuclei(kFALSE), fNewMIS(kFALSE), fHFoff(kFALSE), fNucPdf(0), fTriggerParticle(0), fTriggerEta(0.9), fTriggerMultiplicity(0), fTriggerMultiplicityEta(0), fTriggerMultiplicityPtMin(0), fCountMode(kCountAll), fHeader(0), fRL(0), fFileName(0), fFragPhotonInCalo(kFALSE), fPi0InCalo(kFALSE) , fPhotonInCalo(kFALSE), fEleInEMCAL(kFALSE), fCheckEMCAL(kFALSE), fCheckPHOS(kFALSE), fCheckPHOSeta(kFALSE), fFragPhotonOrPi0MinPt(0), fPhotonMinPt(0), fElectronMinPt(0), fPHOSMinPhi(219.), fPHOSMaxPhi(321.), fPHOSEta(0.13), fEMCALMinPhi(79.), fEMCALMaxPhi(191.), fEMCALEta(0.71) { // Default Constructor fEnergyCMS = 5500.; if (!AliPythiaRndm::GetPythiaRandom()) AliPythiaRndm::SetPythiaRandom(GetRandom()); } AliGenPythia::AliGenPythia(Int_t npart) :AliGenMC(npart), fProcess(kPyCharm), fItune(-1), fStrucFunc(kCTEQ5L), fKineBias(0.), fTrials(0), fTrialsRun(0), fQ(0.), fX1(0.), fX2(0.), fEventTime(0.), fInteractionRate(0.), fTimeWindow(0.), fCurSubEvent(0), fEventsTime(0), fNev(0), fFlavorSelect(0), fXsection(0.), fPythia(0), fPtHardMin(0.), fPtHardMax(1.e4), fYHardMin(-1.e10), fYHardMax(1.e10), fGinit(kTRUE), fGfinal(kTRUE), fHadronisation(kTRUE), fNpartons(0), fReadFromFile(kFALSE), fQuench(kFALSE), fQhat(0.), fLength(0.), fImpact(0.), fPtKick(1.), fFullEvent(kTRUE), fDecayer(new AliDecayerPythia()), fDebugEventFirst(-1), fDebugEventLast(-1), fEtMinJet(0.), fEtMaxJet(1.e4), fEtaMinJet(-20.), fEtaMaxJet(20.), fPhiMinJet(0.), fPhiMaxJet(2.* TMath::Pi()), fJetReconstruction(kCell), fEtaMinGamma(-20.), fEtaMaxGamma(20.), fPhiMinGamma(0.), fPhiMaxGamma(2. * TMath::Pi()), fPycellEtaMax(2.), fPycellNEta(274), fPycellNPhi(432), fPycellThreshold(0.), fPycellEtSeed(4.), fPycellMinEtJet(10.), fPycellMaxRadius(1.), fStackFillOpt(kFlavorSelection), fFeedDownOpt(kTRUE), fFragmentation(kTRUE), fSetNuclei(kFALSE), fNewMIS(kFALSE), fHFoff(kFALSE), fNucPdf(0), fTriggerParticle(0), fTriggerEta(0.9), fTriggerMultiplicity(0), fTriggerMultiplicityEta(0), fTriggerMultiplicityPtMin(0), fCountMode(kCountAll), fHeader(0), fRL(0), fFileName(0), fFragPhotonInCalo(kFALSE), fPi0InCalo(kFALSE) , fPhotonInCalo(kFALSE), fEleInEMCAL(kFALSE), fCheckEMCAL(kFALSE), fCheckPHOS(kFALSE), fCheckPHOSeta(kFALSE), fFragPhotonOrPi0MinPt(0), fPhotonMinPt(0), fElectronMinPt(0), fPHOSMinPhi(219.), fPHOSMaxPhi(321.), fPHOSEta(0.13), fEMCALMinPhi(79.), fEMCALMaxPhi(191.), fEMCALEta(0.71) { // default charm production at 5. 5 TeV // semimuonic decay // structure function GRVHO // fEnergyCMS = 5500.; fName = "Pythia"; fTitle= "Particle Generator using PYTHIA"; SetForceDecay(); // Set random number generator if (!AliPythiaRndm::GetPythiaRandom()) AliPythiaRndm::SetPythiaRandom(GetRandom()); } AliGenPythia::~AliGenPythia() { // Destructor if(fEventsTime) delete fEventsTime; } void AliGenPythia::SetInteractionRate(Float_t rate,Float_t timewindow) { // Generate pileup using user specified rate fInteractionRate = rate; fTimeWindow = timewindow; GeneratePileup(); } void AliGenPythia::GeneratePileup() { // Generate sub events time for pileup fEventsTime = 0; if(fInteractionRate == 0.) { Warning("GeneratePileup","Zero interaction specified. Skipping pileup generation.\n"); return; } Int_t npart = NumberParticles(); if(npart < 0) { Warning("GeneratePileup","Negative number of particles. Skipping pileup generation.\n"); return; } if(fEventsTime) delete fEventsTime; fEventsTime = new TArrayF(npart); TArrayF &array = *fEventsTime; for(Int_t ipart = 0; ipart < npart; ipart++) array[ipart] = 0.; Float_t eventtime = 0.; while(1) { eventtime += (AliPythiaRndm::GetPythiaRandom())->Exp(1./fInteractionRate); if(eventtime > fTimeWindow) break; array.Set(array.GetSize()+1); array[array.GetSize()-1] = eventtime; } eventtime = 0.; while(1) { eventtime -= (AliPythiaRndm::GetPythiaRandom())->Exp(1./fInteractionRate); if(TMath::Abs(eventtime) > fTimeWindow) break; array.Set(array.GetSize()+1); array[array.GetSize()-1] = eventtime; } SetNumberParticles(fEventsTime->GetSize()); } void AliGenPythia::SetPycellParameters(Float_t etamax, Int_t neta, Int_t nphi, Float_t thresh, Float_t etseed, Float_t minet, Float_t r) { // Set pycell parameters fPycellEtaMax = etamax; fPycellNEta = neta; fPycellNPhi = nphi; fPycellThreshold = thresh; fPycellEtSeed = etseed; fPycellMinEtJet = minet; fPycellMaxRadius = r; } void AliGenPythia::SetEventListRange(Int_t eventFirst, Int_t eventLast) { // Set a range of event numbers, for which a table // of generated particle will be printed fDebugEventFirst = eventFirst; fDebugEventLast = eventLast; if (fDebugEventLast==-1) fDebugEventLast=fDebugEventFirst; } void AliGenPythia::Init() { // Initialisation SetMC(AliPythia::Instance()); fPythia=(AliPythia*) fMCEvGen; // fParentWeight=1./Float_t(fNpart); // fPythia->SetCKIN(3,fPtHardMin); fPythia->SetCKIN(4,fPtHardMax); fPythia->SetCKIN(7,fYHardMin); fPythia->SetCKIN(8,fYHardMax); if (fAProjectile > 0 && fATarget > 0) fPythia->SetNuclei(fAProjectile, fATarget, fNucPdf); // Fragmentation? if (fFragmentation) { fPythia->SetMSTP(111,1); } else { fPythia->SetMSTP(111,0); } // initial state radiation fPythia->SetMSTP(61,fGinit); // final state radiation fPythia->SetMSTP(71,fGfinal); // pt - kick if (fPtKick > 0.) { fPythia->SetMSTP(91,1); fPythia->SetPARP(91,fPtKick); fPythia->SetPARP(93, 4. * fPtKick); } else { fPythia->SetMSTP(91,0); } if (fReadFromFile) { fRL = AliRunLoader::Open(fFileName, "Partons"); fRL->LoadKinematics(); fRL->LoadHeader(); } else { fRL = 0x0; } // fPythia->ProcInit(fProcess,fEnergyCMS,fStrucFunc, fItune); // Forward Paramters to the AliPythia object fDecayer->SetForceDecay(fForceDecay); // Switch off Heavy Flavors on request if (fHFoff) { // Maximum number of quark flavours used in pdf fPythia->SetMSTP(58, 3); // Maximum number of flavors that can be used in showers fPythia->SetMSTJ(45, 3); // Switch off g->QQbar splitting in decay table ((AliDecayerPythia*) fDecayer)->HeavyFlavourOff(); } fDecayer->Init(); // Parent and Children Selection switch (fProcess) { case kPyOldUEQ2ordered: case kPyOldUEQ2ordered2: case kPyOldPopcorn: break; case kPyCharm: case kPyCharmUnforced: case kPyCharmPbPbMNR: case kPyCharmpPbMNR: case kPyCharmppMNR: case kPyCharmppMNRwmi: fParentSelect[0] = 411; fParentSelect[1] = 421; fParentSelect[2] = 431; fParentSelect[3] = 4122; fParentSelect[4] = 4232; fParentSelect[5] = 4132; fParentSelect[6] = 4332; fFlavorSelect = 4; break; case kPyD0PbPbMNR: case kPyD0pPbMNR: case kPyD0ppMNR: fParentSelect[0] = 421; fFlavorSelect = 4; break; case kPyDPlusPbPbMNR: case kPyDPluspPbMNR: case kPyDPlusppMNR: fParentSelect[0] = 411; fFlavorSelect = 4; break; case kPyDPlusStrangePbPbMNR: case kPyDPlusStrangepPbMNR: case kPyDPlusStrangeppMNR: fParentSelect[0] = 431; fFlavorSelect = 4; break; case kPyBeauty: case kPyBeautyJets: case kPyBeautyPbPbMNR: case kPyBeautypPbMNR: case kPyBeautyppMNR: case kPyBeautyppMNRwmi: fParentSelect[0]= 511; fParentSelect[1]= 521; fParentSelect[2]= 531; fParentSelect[3]= 5122; fParentSelect[4]= 5132; fParentSelect[5]= 5232; fParentSelect[6]= 5332; fFlavorSelect = 5; break; case kPyBeautyUnforced: fParentSelect[0] = 511; fParentSelect[1] = 521; fParentSelect[2] = 531; fParentSelect[3] = 5122; fParentSelect[4] = 5132; fParentSelect[5] = 5232; fParentSelect[6] = 5332; fFlavorSelect = 5; break; case kPyJpsiChi: case kPyJpsi: fParentSelect[0] = 443; break; case kPyMbDefault: case kPyMbAtlasTuneMC09: case kPyMb: case kPyMbWithDirectPhoton: case kPyMbNonDiffr: case kPyMbMSEL1: case kPyJets: case kPyDirectGamma: case kPyLhwgMb: break; case kPyW: case kPyZ: break; } // // // JetFinder for Trigger // // Configure detector (EMCAL like) // fPythia->SetPARU(51, fPycellEtaMax); fPythia->SetMSTU(51, fPycellNEta); fPythia->SetMSTU(52, fPycellNPhi); // // Configure Jet Finder // fPythia->SetPARU(58, fPycellThreshold); fPythia->SetPARU(52, fPycellEtSeed); fPythia->SetPARU(53, fPycellMinEtJet); fPythia->SetPARU(54, fPycellMaxRadius); fPythia->SetMSTU(54, 2); // // This counts the total number of calls to Pyevnt() per run. fTrialsRun = 0; fQ = 0.; fX1 = 0.; fX2 = 0.; fNev = 0 ; // // // AliGenMC::Init(); // // // if (fSetNuclei) { fDyBoost = 0; Warning("Init","SetNuclei used. Use SetProjectile + SetTarget instead. fDyBoost has been reset to 0\n"); } fPythia->SetPARJ(200, 0.0); fPythia->SetPARJ(199, 0.0); fPythia->SetPARJ(198, 0.0); fPythia->SetPARJ(197, 0.0); if (fQuench == 1) { fPythia->InitQuenching(0., 0.1, 0.6e6, 0); } if (fQuench == 3) { // Nestor's change of the splittings fPythia->SetPARJ(200, 0.8); fPythia->SetMSTJ(41, 1); // QCD radiation only fPythia->SetMSTJ(42, 2); // angular ordering fPythia->SetMSTJ(44, 2); // option to run alpha_s fPythia->SetMSTJ(47, 0); // No correction back to hard scattering element fPythia->SetMSTJ(50, 0); // No coherence in first branching fPythia->SetPARJ(82, 1.); // Cut off for parton showers } else if (fQuench == 4) { // Armesto-Cunqueiro-Salgado change of the splittings. AliFastGlauber* glauber = AliFastGlauber::Instance(); glauber->Init(2); //read and store transverse almonds corresponding to differnt //impact parameters. glauber->SetCentralityClass(0.,0.1); fPythia->SetPARJ(200, 1.); fPythia->SetPARJ(198, fQhat); fPythia->SetPARJ(199, fLength); fPythia->SetMSTJ(42, 2); // angular ordering fPythia->SetMSTJ(44, 2); // option to run alpha_s fPythia->SetPARJ(82, 1.); // Cut off for parton showers } } void AliGenPythia::Generate() { // Generate one event if (!fPythia) fPythia=(AliPythia*) fMCEvGen; fDecayer->ForceDecay(); Float_t polar[3] = {0,0,0}; Float_t origin[3] = {0,0,0}; Float_t p[4]; // converts from mm/c to s const Float_t kconv=0.001/2.999792458e8; // Int_t nt=0; Int_t jev=0; Int_t j, kf; fTrials=0; fEventTime = 0.; // Set collision vertex position if (fVertexSmear == kPerEvent) Vertex(); // event loop while(1) { // // Produce event // // // Switch hadronisation off // fPythia->SetMSTJ(1, 0); if (fQuench ==4){ Double_t bimp; // Quenching comes through medium-modified splitting functions. AliFastGlauber::Instance()->GetRandomBHard(bimp); fPythia->SetPARJ(197, bimp); fImpact = bimp; } // // Either produce new event or read partons from file // if (!fReadFromFile) { if (!fNewMIS) { fPythia->Pyevnt(); } else { fPythia->Pyevnw(); } fNpartons = fPythia->GetN(); } else { printf("Loading Event %d\n",AliRunLoader::Instance()->GetEventNumber()); fRL->GetEvent(AliRunLoader::Instance()->GetEventNumber()); fPythia->SetN(0); LoadEvent(fRL->Stack(), 0 , 1); fPythia->Pyedit(21); } // // Run quenching routine // if (fQuench == 1) { fPythia->Quench(); } else if (fQuench == 2){ fPythia->Pyquen(208., 0, 0.); } else if (fQuench == 3) { // Quenching is via multiplicative correction of the splittings } // // Switch hadronisation on // if (fHadronisation) { fPythia->SetMSTJ(1, 1); // // .. and perform hadronisation // printf("Calling hadronisation %d\n", fPythia->GetN()); fPythia->Pyexec(); } fTrials++; fPythia->ImportParticles(&fParticles,"All"); if (TMath::Abs(fDyBoost) > 1.e-4) Boost(); // // // Int_t i; fNprimaries = 0; Int_t np = fParticles.GetEntriesFast(); if (np == 0) continue; // // Int_t* pParent = new Int_t[np]; Int_t* pSelected = new Int_t[np]; Int_t* trackIt = new Int_t[np]; for (i = 0; i < np; i++) { pParent[i] = -1; pSelected[i] = 0; trackIt[i] = 0; } Int_t nc = 0; // Total n. of selected particles Int_t nParents = 0; // Selected parents Int_t nTkbles = 0; // Trackable particles if (fProcess != kPyMbDefault && fProcess != kPyMb && fProcess != kPyMbAtlasTuneMC09 && fProcess != kPyMbWithDirectPhoton && fProcess != kPyJets && fProcess != kPyDirectGamma && fProcess != kPyMbNonDiffr && fProcess != kPyMbMSEL1 && fProcess != kPyW && fProcess != kPyZ && fProcess != kPyCharmppMNRwmi && fProcess != kPyBeautyppMNRwmi && fProcess != kPyBeautyJets) { for (i = 0; i < np; i++) { TParticle* iparticle = (TParticle *) fParticles.At(i); Int_t ks = iparticle->GetStatusCode(); kf = CheckPDGCode(iparticle->GetPdgCode()); // No initial state partons if (ks==21) continue; // // Heavy Flavor Selection // // quark ? kf = TMath::Abs(kf); Int_t kfl = kf; // Resonance if (kfl > 100000) kfl %= 100000; if (kfl > 10000) kfl %= 10000; // meson ? if (kfl > 10) kfl/=100; // baryon if (kfl > 10) kfl/=10; Int_t ipa = iparticle->GetFirstMother()-1; Int_t kfMo = 0; // // Establish mother daughter relation between heavy quarks and mesons // if (kf >= fFlavorSelect && kf <= 6) { Int_t idau = iparticle->GetFirstDaughter() - 1; if (idau > -1) { TParticle* daughter = (TParticle *) fParticles.At(idau); Int_t pdgD = daughter->GetPdgCode(); if (pdgD == 91 || pdgD == 92) { Int_t jmin = daughter->GetFirstDaughter() - 1; Int_t jmax = daughter->GetLastDaughter() - 1; for (Int_t jp = jmin; jp <= jmax; jp++) ((TParticle *) fParticles.At(jp))->SetFirstMother(i+1); } // is string or cluster } // has daughter } // heavy quark if (ipa > -1) { TParticle * mother = (TParticle *) fParticles.At(ipa); kfMo = TMath::Abs(mother->GetPdgCode()); } // What to keep in Stack? Bool_t flavorOK = kFALSE; Bool_t selectOK = kFALSE; if (fFeedDownOpt) { if (kfl >= fFlavorSelect) flavorOK = kTRUE; } else { if (kfl > fFlavorSelect) { nc = -1; break; } if (kfl == fFlavorSelect) flavorOK = kTRUE; } switch (fStackFillOpt) { case kFlavorSelection: selectOK = kTRUE; break; case kParentSelection: if (ParentSelected(kf) || kf <= 10) selectOK = kTRUE; break; } if (flavorOK && selectOK) { // // Heavy flavor hadron or quark // // Kinematic seletion on final state heavy flavor mesons if (ParentSelected(kf) && !KinematicSelection(iparticle, 0)) { continue; } pSelected[i] = 1; if (ParentSelected(kf)) ++nParents; // Update parent count // printf("\n particle (HF) %d %d %d", i, pSelected[i], kf); } else { // Kinematic seletion on decay products if (fCutOnChild && ParentSelected(kfMo) && ChildSelected(kf) && !KinematicSelection(iparticle, 1)) { continue; } // // Decay products // Select if mother was selected and is not tracked if (pSelected[ipa] && !trackIt[ipa] && // mother will be tracked ? kfMo != 5 && // mother is b-quark, don't store fragments kfMo != 4 && // mother is c-quark, don't store fragments kf != 92) // don't store string { // // Semi-stable or de-selected: diselect decay products: // // if (pSelected[i] == -1 || fDecayer->GetLifetime(kf) > fMaxLifeTime) { Int_t ipF = iparticle->GetFirstDaughter(); Int_t ipL = iparticle->GetLastDaughter(); if (ipF > 0) for (j = ipF-1; j < ipL; j++) pSelected[j] = -1; } // printf("\n particle (decay) %d %d %d", i, pSelected[i], kf); pSelected[i] = (pSelected[i] == -1) ? 0 : 1; } } if (pSelected[i] == -1) pSelected[i] = 0; if (!pSelected[i]) continue; // Count quarks only if you did not include fragmentation if (fFragmentation && kf <= 10) continue; nc++; // Decision on tracking trackIt[i] = 0; // // Track final state particle if (ks == 1) trackIt[i] = 1; // Track semi-stable particles if ((ks == 1) || (fDecayer->GetLifetime(kf) > fMaxLifeTime)) trackIt[i] = 1; // Track particles selected by process if undecayed. if (fForceDecay == kNoDecay) { if (ParentSelected(kf)) trackIt[i] = 1; } else { if (ParentSelected(kf)) trackIt[i] = 0; } if (trackIt[i] == 1) ++nTkbles; // Update trackable counter // // } // particle selection loop if (nc > 0) { for (i = 0; iGetPdgCode()); Int_t ks = iparticle->GetStatusCode(); p[0] = iparticle->Px(); p[1] = iparticle->Py(); p[2] = iparticle->Pz(); p[3] = iparticle->Energy(); origin[0] = fVertex[0]+iparticle->Vx()/10; // [cm] origin[1] = fVertex[1]+iparticle->Vy()/10; // [cm] origin[2] = fVertex[2]+iparticle->Vz()/10; // [cm] Float_t tof = kconv*iparticle->T(); Int_t ipa = iparticle->GetFirstMother()-1; Int_t iparent = (ipa > -1) ? pParent[ipa] : -1; PushTrack(fTrackIt*trackIt[i], iparent, kf, p[0], p[1], p[2], p[3], origin[0], origin[1], origin[2], tof, polar[0], polar[1], polar[2], kPPrimary, nt, 1., ks); pParent[i] = nt; KeepTrack(nt); fNprimaries++; } // PushTrack loop } } else { nc = GenerateMB(); } // mb ? GetSubEventTime(); delete[] pParent; delete[] pSelected; delete[] trackIt; if (nc > 0) { switch (fCountMode) { case kCountAll: // printf(" Count all \n"); jev += nc; break; case kCountParents: // printf(" Count parents \n"); jev += nParents; break; case kCountTrackables: // printf(" Count trackable \n"); jev += nTkbles; break; } if (jev >= fNpart || fNpart == -1) { fKineBias=Float_t(fNpart)/Float_t(fTrials); fQ += fPythia->GetVINT(51); fX1 += fPythia->GetVINT(41); fX2 += fPythia->GetVINT(42); fTrialsRun += fTrials; fNev++; MakeHeader(); break; } } } // event loop SetHighWaterMark(nt); // adjust weight due to kinematic selection // AdjustWeights(); // get cross-section fXsection=fPythia->GetPARI(1); } Int_t AliGenPythia::GenerateMB() { // // Min Bias selection and other global selections // Int_t i, kf, nt, iparent; Int_t nc = 0; Float_t p[4]; Float_t polar[3] = {0,0,0}; Float_t origin[3] = {0,0,0}; // converts from mm/c to s const Float_t kconv=0.001/2.999792458e8; Int_t np = (fHadronisation) ? fParticles.GetEntriesFast() : fNpartons; Int_t* pParent = new Int_t[np]; for (i=0; i< np; i++) pParent[i] = -1; if (fProcess == kPyJets || fProcess == kPyDirectGamma || fProcess == kPyBeautyJets || fProcess == kPyBeautyppMNRwmi) { TParticle* jet1 = (TParticle *) fParticles.At(6); TParticle* jet2 = (TParticle *) fParticles.At(7); if (!CheckTrigger(jet1, jet2)) { delete [] pParent; return 0; } } // Select jets with fragmentation photon or pi0 going to PHOS or EMCAL if (fProcess == kPyJets && (fFragPhotonInCalo || fPi0InCalo) ) { Bool_t ok = kFALSE; Int_t pdg = 0; if (fFragPhotonInCalo) pdg = 22 ; // Photon else if (fPi0InCalo) pdg = 111 ; // Pi0 for (i=0; i< np; i++) { TParticle* iparticle = (TParticle *) fParticles.At(i); if(iparticle->GetStatusCode()==1 && iparticle->GetPdgCode()==pdg && iparticle->Pt() > fFragPhotonOrPi0MinPt){ Int_t imother = iparticle->GetFirstMother() - 1; TParticle* pmother = (TParticle *) fParticles.At(imother); if(pdg == 111 || (pdg == 22 && pmother->GetStatusCode() != 11)) //No photon from hadron decay { Float_t phi = iparticle->Phi()*180./TMath::Pi(); //Convert to degrees Float_t eta =TMath::Abs(iparticle->Eta()); //in calos etamin=-etamax if((fCheckEMCAL && IsInEMCAL(phi,eta)) || (fCheckPHOS && IsInPHOS(phi,eta)) ) ok =kTRUE; } } } if(!ok) return 0; } // Select beauty jets with electron in EMCAL if (fProcess == kPyBeautyJets && fEleInEMCAL) { Bool_t ok = kFALSE; Int_t pdg = 11; //electron Float_t pt = 0.; Float_t eta = 0.; Float_t phi = 0.; for (i=0; i< np; i++) { TParticle* iparticle = (TParticle *) fParticles.At(i); if(iparticle->GetStatusCode()==1 && TMath::Abs(iparticle->GetPdgCode())==pdg && iparticle->Pt() > fElectronMinPt){ pt = iparticle->Pt(); phi = iparticle->Phi()*180./TMath::Pi(); //Convert to degrees eta =TMath::Abs(iparticle->Eta()); //in calos etamin=-etamax if(IsInEMCAL(phi,eta)) ok =kTRUE; } } if(!ok) return 0; AliDebug(5,Form("Found an electron jet (pt,eta,phi) = (%f,%f,%f)",pt,eta,phi)); } // Check for minimum multiplicity if (fTriggerMultiplicity > 0) { Int_t multiplicity = 0; for (i = 0; i < np; i++) { TParticle * iparticle = (TParticle *) fParticles.At(i); Int_t statusCode = iparticle->GetStatusCode(); // Initial state particle if (statusCode != 1) continue; // eta cut if (fTriggerMultiplicityEta > 0 && TMath::Abs(iparticle->Eta()) > fTriggerMultiplicityEta) continue; // pt cut if (iparticle->Pt() < fTriggerMultiplicityPtMin) continue; TParticlePDG* pdgPart = iparticle->GetPDG(); if (pdgPart && pdgPart->Charge() == 0) continue; ++multiplicity; } if (multiplicity < fTriggerMultiplicity) { delete [] pParent; return 0; } Printf("Triggered on event with multiplicity of %d >= %d", multiplicity, fTriggerMultiplicity); } // Select events with a photon pt > min pt going to PHOS eta acceptance or exactly PHOS eta phi if ((fProcess == kPyJets || fProcess == kPyDirectGamma) && fPhotonInCalo && (fCheckPHOSeta || fCheckPHOS)){ Bool_t okd = kFALSE; Int_t pdg = 22; Int_t iphcand = -1; for (i=0; i< np; i++) { TParticle* iparticle = (TParticle *) fParticles.At(i); Float_t phi = iparticle->Phi()*180./TMath::Pi(); //Convert to degrees Float_t eta =TMath::Abs(iparticle->Eta());//in calos etamin=-etamax if(iparticle->GetStatusCode() == 1 && iparticle->GetPdgCode() == pdg && iparticle->Pt() > fPhotonMinPt && eta < fPHOSEta){ // first check if the photon is in PHOS phi if(IsInPHOS(phi,eta)){ okd = kTRUE; break; } if(fCheckPHOSeta) iphcand = i; // candiate photon to rotate in phi } } if(!okd && iphcand != -1) // execute rotation in phi RotatePhi(iphcand,okd); if(!okd) return 0; } if (fTriggerParticle) { Bool_t triggered = kFALSE; for (i = 0; i < np; i++) { TParticle * iparticle = (TParticle *) fParticles.At(i); kf = CheckPDGCode(iparticle->GetPdgCode()); if (kf != fTriggerParticle) continue; if (iparticle->Pt() == 0.) continue; if (TMath::Abs(iparticle->Eta()) > fTriggerEta) continue; triggered = kTRUE; break; } if (!triggered) { delete [] pParent; return 0; } } // Check if there is a ccbar or bbbar pair with at least one of the two // in fYMin < y < fYMax if (fProcess == kPyCharmppMNRwmi || fProcess == kPyBeautyppMNRwmi || fProcess == kPyBeautyJets) { TParticle *partCheck; TParticle *mother; Bool_t theQ=kFALSE,theQbar=kFALSE,inYcut=kFALSE; Bool_t theChild=kFALSE; Float_t y; Int_t pdg,mpdg,mpdgUpperFamily; for(i=0; iGetPdgCode(); if(TMath::Abs(pdg) == fFlavorSelect) { // quark if(pdg>0) { theQ=kTRUE; } else { theQbar=kTRUE; } y = 0.5*TMath::Log((partCheck->Energy()+partCheck->Pz()+1.e-13)/ (partCheck->Energy()-partCheck->Pz()+1.e-13)); if(y>fYMin && yGetFirstMother() - 1; if(mi<0) continue; mother = (TParticle*)fParticles.At(mi); mpdg=TMath::Abs(mother->GetPdgCode()); mpdgUpperFamily=(mpdg>1000 ? mpdg+1000 : mpdg+100); // keep e from c from b if ( ParentSelected(mpdg) || (fFlavorSelect==5 && ParentSelected(mpdgUpperFamily))) { if (KinematicSelection(partCheck,1)) { theChild=kTRUE; } } } } if (!theQ || !theQbar || !inYcut) { // one of the c/b conditions not satisfied delete[] pParent; return 0; } if (fCutOnChild && !theChild) { // one of the child conditions not satisfied delete[] pParent; return 0; } } //Introducing child cuts in case kPyW, kPyZ, kPyMb, and kPyMbNonDiff if ( (fProcess == kPyW || fProcess == kPyZ || fProcess == kPyMbDefault || fProcess == kPyMb || fProcess == kPyMbAtlasTuneMC09 || fProcess == kPyMbWithDirectPhoton || fProcess == kPyMbNonDiffr) && (fCutOnChild == 1) ) { if ( !CheckKinematicsOnChild() ) { delete[] pParent; return 0; } } for (i = 0; i < np; i++) { Int_t trackIt = 0; TParticle * iparticle = (TParticle *) fParticles.At(i); kf = CheckPDGCode(iparticle->GetPdgCode()); Int_t ks = iparticle->GetStatusCode(); Int_t km = iparticle->GetFirstMother(); if ((ks == 1 && kf!=0 && KinematicSelection(iparticle, 0)) || (ks != 1) || ((fProcess == kPyJets || fProcess == kPyBeautyJets) && ks == 21 && km == 0 && i>1)) { nc++; if (ks == 1) trackIt = 1; Int_t ipa = iparticle->GetFirstMother()-1; iparent = (ipa > -1) ? pParent[ipa] : -1; // // store track information p[0] = iparticle->Px(); p[1] = iparticle->Py(); p[2] = iparticle->Pz(); p[3] = iparticle->Energy(); origin[0] = fVertex[0]+iparticle->Vx()/10; // [cm] origin[1] = fVertex[1]+iparticle->Vy()/10; // [cm] origin[2] = fVertex[2]+iparticle->Vz()/10; // [cm] Float_t tof = fEventTime + kconv * iparticle->T(); PushTrack(fTrackIt*trackIt, iparent, kf, p[0], p[1], p[2], p[3], origin[0], origin[1], origin[2], tof, polar[0], polar[1], polar[2], kPPrimary, nt, 1., ks); fNprimaries++; KeepTrack(nt); pParent[i] = nt; SetHighWaterMark(nt); } // select particle } // particle loop delete[] pParent; return 1; } void AliGenPythia::FinishRun() { // Print x-section summary fPythia->Pystat(1); if (fNev > 0.) { fQ /= fNev; fX1 /= fNev; fX2 /= fNev; } printf("\nTotal number of Pyevnt() calls %d\n", fTrialsRun); printf("\nMean Q, x1, x2: %f %f %f\n", fQ, fX1, fX2); } void AliGenPythia::AdjustWeights() const { // Adjust the weights after generation of all events // if (gAlice) { TParticle *part; Int_t ntrack=gAlice->GetMCApp()->GetNtrack(); for (Int_t i=0; iGetMCApp()->Particle(i); part->SetWeight(part->GetWeight()*fKineBias); } } } void AliGenPythia::SetNuclei(Int_t a1, Int_t a2, Int_t pdfset) { // Treat protons as inside nuclei with mass numbers a1 and a2 fAProjectile = a1; fATarget = a2; fNucPdf = pdfset; // 0 EKS98 1 EPS08 fSetNuclei = kTRUE; } void AliGenPythia::MakeHeader() { // // Make header for the simulated event // if (gAlice) { if (gAlice->GetEvNumber()>=fDebugEventFirst && gAlice->GetEvNumber()<=fDebugEventLast) fPythia->Pylist(2); } // Builds the event header, to be called after each event if (fHeader) delete fHeader; fHeader = new AliGenPythiaEventHeader("Pythia"); // // Event type ((AliGenPythiaEventHeader*) fHeader)->SetProcessType(fPythia->GetMSTI(1)); // // Number of trials ((AliGenPythiaEventHeader*) fHeader)->SetTrials(fTrials); // // Event Vertex fHeader->SetPrimaryVertex(fVertex); // // Number of primaries fHeader->SetNProduced(fNprimaries); // // Jets that have triggered //Need to store jets for b-jet studies too! if (fProcess == kPyJets || fProcess == kPyDirectGamma || fProcess == kPyBeautyJets || fProcess == kPyBeautyppMNRwmi) { Int_t ntrig, njet; Float_t jets[4][10]; GetJets(njet, ntrig, jets); for (Int_t i = 0; i < ntrig; i++) { ((AliGenPythiaEventHeader*) fHeader)->AddJet(jets[0][i], jets[1][i], jets[2][i], jets[3][i]); } } // // Copy relevant information from external header, if present. // Float_t uqJet[4]; if (fRL) { AliGenPythiaEventHeader* exHeader = (AliGenPythiaEventHeader*) (fRL->GetHeader()->GenEventHeader()); for (Int_t i = 0; i < exHeader->NTriggerJets(); i++) { printf("Adding Jet %d %d \n", i, exHeader->NTriggerJets()); exHeader->TriggerJet(i, uqJet); ((AliGenPythiaEventHeader*) fHeader)->AddUQJet(uqJet[0], uqJet[1], uqJet[2], uqJet[3]); } } // // Store quenching parameters // if (fQuench){ Double_t z[4]; Double_t xp, yp; if (fQuench == 1) { // Pythia::Quench() fPythia->GetQuenchingParameters(xp, yp, z); } else if (fQuench == 2){ // Pyquen Double_t r1 = PARIMP.rb1; Double_t r2 = PARIMP.rb2; Double_t b = PARIMP.b1; Double_t r = 0.5 * TMath::Sqrt(2. * (r1 * r1 + r2 * r2) - b * b); Double_t phi = PARIMP.psib1; xp = r * TMath::Cos(phi); yp = r * TMath::Sin(phi); } else if (fQuench == 4) { // QPythia Double_t xy[2]; Double_t i0i1[2]; AliFastGlauber::Instance()->GetSavedXY(xy); AliFastGlauber::Instance()->GetSavedI0I1(i0i1); xp = xy[0]; yp = xy[1]; ((AliGenPythiaEventHeader*) fHeader)->SetImpactParameter(fImpact); } ((AliGenPythiaEventHeader*) fHeader)->SetXYJet(xp, yp); ((AliGenPythiaEventHeader*) fHeader)->SetZQuench(z); } // // Store pt^hard ((AliGenPythiaEventHeader*) fHeader)->SetPtHard(fPythia->GetVINT(47)); // // Pass header // AddHeader(fHeader); fHeader = 0x0; } Bool_t AliGenPythia::CheckTrigger(TParticle* jet1, TParticle* jet2) { // Check the kinematic trigger condition // Double_t eta[2]; eta[0] = jet1->Eta(); eta[1] = jet2->Eta(); Double_t phi[2]; phi[0] = jet1->Phi(); phi[1] = jet2->Phi(); Int_t pdg[2]; pdg[0] = jet1->GetPdgCode(); pdg[1] = jet2->GetPdgCode(); Bool_t triggered = kFALSE; if (fProcess == kPyJets || fProcess == kPyBeautyJets || fProcess == kPyBeautyppMNRwmi) { Int_t njets = 0; Int_t ntrig = 0; Float_t jets[4][10]; // // Use Pythia clustering on parton level to determine jet axis // GetJets(njets, ntrig, jets); if (ntrig || fEtMinJet == 0.) triggered = kTRUE; // } else { Int_t ij = 0; Int_t ig = 1; if (pdg[0] == kGamma) { ij = 1; ig = 0; } //Check eta range first... if ((eta[ij] < fEtaMaxJet && eta[ij] > fEtaMinJet) && (eta[ig] < fEtaMaxGamma && eta[ig] > fEtaMinGamma)) { //Eta is okay, now check phi range if ((phi[ij] < fPhiMaxJet && phi[ij] > fPhiMinJet) && (phi[ig] < fPhiMaxGamma && phi[ig] > fPhiMinGamma)) { triggered = kTRUE; } } } return triggered; } Bool_t AliGenPythia::CheckKinematicsOnChild(){ // //Checking Kinematics on Child (status code 1, particle code ?, kin cuts // Bool_t checking = kFALSE; Int_t j, kcode, ks, km; Int_t nPartAcc = 0; //number of particles in the acceptance range Int_t numberOfAcceptedParticles = 1; if (fNumberOfAcceptedParticles != 0) { numberOfAcceptedParticles = fNumberOfAcceptedParticles; } Int_t npart = fParticles.GetEntriesFast(); for (j = 0; jGetPdgCode()) ); ks = jparticle->GetStatusCode(); km = jparticle->GetFirstMother(); if( (ks == 1) && (kcode == fPdgCodeParticleforAcceptanceCut) && (KinematicSelection(jparticle,1)) ){ nPartAcc++; } if( numberOfAcceptedParticles <= nPartAcc){ checking = kTRUE; break; } } return checking; } void AliGenPythia::LoadEvent(AliStack* stack, Int_t flag, Int_t reHadr) { // // Load event into Pythia Common Block // Int_t npart = stack -> GetNprimary(); Int_t n0 = 0; if (!flag) { (fPythia->GetPyjets())->N = npart; } else { n0 = (fPythia->GetPyjets())->N; (fPythia->GetPyjets())->N = n0 + npart; } for (Int_t part = 0; part < npart; part++) { TParticle *mPart = stack->Particle(part); Int_t kf = mPart->GetPdgCode(); Int_t ks = mPart->GetStatusCode(); Int_t idf = mPart->GetFirstDaughter(); Int_t idl = mPart->GetLastDaughter(); if (reHadr) { if (ks == 11 || ks == 12) { ks -= 10; idf = -1; idl = -1; } } Float_t px = mPart->Px(); Float_t py = mPart->Py(); Float_t pz = mPart->Pz(); Float_t e = mPart->Energy(); Float_t m = mPart->GetCalcMass(); (fPythia->GetPyjets())->P[0][part+n0] = px; (fPythia->GetPyjets())->P[1][part+n0] = py; (fPythia->GetPyjets())->P[2][part+n0] = pz; (fPythia->GetPyjets())->P[3][part+n0] = e; (fPythia->GetPyjets())->P[4][part+n0] = m; (fPythia->GetPyjets())->K[1][part+n0] = kf; (fPythia->GetPyjets())->K[0][part+n0] = ks; (fPythia->GetPyjets())->K[3][part+n0] = idf + 1; (fPythia->GetPyjets())->K[4][part+n0] = idl + 1; (fPythia->GetPyjets())->K[2][part+n0] = mPart->GetFirstMother() + 1; } } void AliGenPythia::LoadEvent(TObjArray* stack, Int_t flag, Int_t reHadr) { // // Load event into Pythia Common Block // Int_t npart = stack -> GetEntries(); Int_t n0 = 0; if (!flag) { (fPythia->GetPyjets())->N = npart; } else { n0 = (fPythia->GetPyjets())->N; (fPythia->GetPyjets())->N = n0 + npart; } for (Int_t part = 0; part < npart; part++) { TParticle *mPart = dynamic_cast(stack->At(part)); Int_t kf = mPart->GetPdgCode(); Int_t ks = mPart->GetStatusCode(); Int_t idf = mPart->GetFirstDaughter(); Int_t idl = mPart->GetLastDaughter(); if (reHadr) { if (ks == 11 || ks == 12) { ks -= 10; idf = -1; idl = -1; } } Float_t px = mPart->Px(); Float_t py = mPart->Py(); Float_t pz = mPart->Pz(); Float_t e = mPart->Energy(); Float_t m = mPart->GetCalcMass(); (fPythia->GetPyjets())->P[0][part+n0] = px; (fPythia->GetPyjets())->P[1][part+n0] = py; (fPythia->GetPyjets())->P[2][part+n0] = pz; (fPythia->GetPyjets())->P[3][part+n0] = e; (fPythia->GetPyjets())->P[4][part+n0] = m; (fPythia->GetPyjets())->K[1][part+n0] = kf; (fPythia->GetPyjets())->K[0][part+n0] = ks; (fPythia->GetPyjets())->K[3][part+n0] = idf + 1; (fPythia->GetPyjets())->K[4][part+n0] = idl + 1; (fPythia->GetPyjets())->K[2][part+n0] = mPart->GetFirstMother() + 1; } } void AliGenPythia::RecJetsUA1(Int_t& njets, Float_t jets [4][50]) { // // Calls the Pythia jet finding algorithm to find jets in the current event // // // // Save jets Int_t n = fPythia->GetN(); // // Run Jet Finder fPythia->Pycell(njets); Int_t i; for (i = 0; i < njets; i++) { Float_t px = (fPythia->GetPyjets())->P[0][n+i]; Float_t py = (fPythia->GetPyjets())->P[1][n+i]; Float_t pz = (fPythia->GetPyjets())->P[2][n+i]; Float_t e = (fPythia->GetPyjets())->P[3][n+i]; jets[0][i] = px; jets[1][i] = py; jets[2][i] = pz; jets[3][i] = e; } } void AliGenPythia::GetJets(Int_t& nJets, Int_t& nJetsTrig, Float_t jets[4][10]) { // // Calls the Pythia clustering algorithm to find jets in the current event // Int_t n = fPythia->GetN(); nJets = 0; nJetsTrig = 0; if (fJetReconstruction == kCluster) { // // Configure cluster algorithm // fPythia->SetPARU(43, 2.); fPythia->SetMSTU(41, 1); // // Call cluster algorithm // fPythia->Pyclus(nJets); // // Loading jets from common block // } else { // // Run Jet Finder fPythia->Pycell(nJets); } Int_t i; for (i = 0; i < nJets; i++) { Float_t px = (fPythia->GetPyjets())->P[0][n+i]; Float_t py = (fPythia->GetPyjets())->P[1][n+i]; Float_t pz = (fPythia->GetPyjets())->P[2][n+i]; Float_t e = (fPythia->GetPyjets())->P[3][n+i]; Float_t pt = TMath::Sqrt(px * px + py * py); Float_t phi = TMath::Pi() + TMath::ATan2(-py, -px); Float_t theta = TMath::ATan2(pt,pz); Float_t et = e * TMath::Sin(theta); Float_t eta = -TMath::Log(TMath::Tan(theta / 2.)); if ( eta > fEtaMinJet && eta < fEtaMaxJet && phi > fPhiMinJet && phi < fPhiMaxJet && et > fEtMinJet && et < fEtMaxJet ) { jets[0][nJetsTrig] = px; jets[1][nJetsTrig] = py; jets[2][nJetsTrig] = pz; jets[3][nJetsTrig] = e; nJetsTrig++; // printf("\n........-Jet #%d: %10.3f %10.3f %10.3f %10.3f \n", i, pt, et, eta, phi * kRaddeg); } else { // printf("\n........-Jet #%d: %10.3f %10.3f %10.3f %10.3f \n", i, pt, et, eta, phi * kRaddeg); } } } void AliGenPythia::GetSubEventTime() { // Calculates time of the next subevent fEventTime = 0.; if (fEventsTime) { TArrayF &array = *fEventsTime; fEventTime = array[fCurSubEvent++]; } // printf(" Event time: %d %f %p",fCurSubEvent,fEventTime,fEventsTime); return; } Bool_t AliGenPythia::IsInEMCAL(Float_t phi, Float_t eta) { // Is particle in EMCAL acceptance? // phi in degrees, etamin=-etamax if(phi > fEMCALMinPhi && phi < fEMCALMaxPhi && eta < fEMCALEta ) return kTRUE; else return kFALSE; } Bool_t AliGenPythia::IsInPHOS(Float_t phi, Float_t eta) { // Is particle in PHOS acceptance? // Acceptance slightly larger considered. // phi in degrees, etamin=-etamax if(phi > fPHOSMinPhi && phi < fPHOSMaxPhi && eta < fPHOSEta ) return kTRUE; else return kFALSE; } void AliGenPythia::RotatePhi(Int_t iphcand, Bool_t& okdd) { //calculate the new position random between fPHOSMinPhi and fPHOSMaxPhi Double_t phiPHOSmin = TMath::Pi()*fPHOSMinPhi/180; Double_t phiPHOSmax = TMath::Pi()*fPHOSMaxPhi/180; Double_t phiPHOS = gRandom->Uniform(phiPHOSmin,phiPHOSmax); //calculate deltaphi TParticle* ph = (TParticle *) fParticles.At(iphcand); Double_t phphi = ph->Phi(); Double_t deltaphi = phiPHOS - phphi; //loop for all particles and produce the phi rotation Int_t np = (fHadronisation) ? fParticles.GetEntriesFast() : fNpartons; Double_t oldphi, newphi; Double_t newVx, newVy, R, Vz, time; Double_t newPx, newPy, pt, Pz, e; for(Int_t i=0; i< np; i++) { TParticle* iparticle = (TParticle *) fParticles.At(i); oldphi = iparticle->Phi(); newphi = oldphi + deltaphi; if(newphi < 0) newphi = 2*TMath::Pi() + newphi; // correct angle if(newphi > 2*TMath::Pi()) newphi = newphi - 2*TMath::Pi(); // correct angle R = iparticle->R(); newVx = R*TMath::Cos(newphi); newVy = R*TMath::Sin(newphi); Vz = iparticle->Vz(); // don't transform time = iparticle->T(); // don't transform pt = iparticle->Pt(); newPx = pt*TMath::Cos(newphi); newPy = pt*TMath::Sin(newphi); Pz = iparticle->Pz(); // don't transform e = iparticle->Energy(); // don't transform // apply rotation iparticle->SetProductionVertex(newVx, newVy, Vz, time); iparticle->SetMomentum(newPx, newPy, Pz, e); } //end particle loop // now let's check that we put correctly the candidate photon in PHOS Float_t phi = ph->Phi()*180./TMath::Pi(); //Convert to degrees Float_t eta =TMath::Abs(ph->Eta());//in calos etamin=-etamax if(IsInPHOS(phi,eta)) okdd = kTRUE; } #ifdef never void AliGenPythia::Streamer(TBuffer &R__b) { // Stream an object of class AliGenPythia. if (R__b.IsReading()) { Version_t R__v = R__b.ReadVersion(); if (R__v) { } AliGenerator::Streamer(R__b); R__b >> (Int_t&)fProcess; R__b >> (Int_t&)fStrucFunc; R__b >> (Int_t&)fForceDecay; R__b >> fEnergyCMS; R__b >> fKineBias; R__b >> fTrials; fParentSelect.Streamer(R__b); fChildSelect.Streamer(R__b); R__b >> fXsection; // (AliPythia::Instance())->Streamer(R__b); R__b >> fPtHardMin; R__b >> fPtHardMax; // if (fDecayer) fDecayer->Streamer(R__b); } else { R__b.WriteVersion(AliGenPythia::IsA()); AliGenerator::Streamer(R__b); R__b << (Int_t)fProcess; R__b << (Int_t)fStrucFunc; R__b << (Int_t)fForceDecay; R__b << fEnergyCMS; R__b << fKineBias; R__b << fTrials; fParentSelect.Streamer(R__b); fChildSelect.Streamer(R__b); R__b << fXsection; // R__b << fPythia; R__b << fPtHardMin; R__b << fPtHardMax; // fDecayer->Streamer(R__b); } } #endif