/************************************************************************** * 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 "AliConst.h" #include "AliDecayerPythia.h" #include "AliGenPythia.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" ClassImp(AliGenPythia) AliGenPythia::AliGenPythia() :AliGenMC() { // Default Constructor fParticles = 0; fPythia = 0; fHeader = 0; fReadFromFile = 0; fDecayer = new AliDecayerPythia(); SetEventListRange(); SetJetPhiRange(); SetJetEtaRange(); SetJetEtRange(); SetGammaPhiRange(); SetGammaEtaRange(); SetPtKick(); SetQuench(); SetHadronisation(); fSetNuclei = kFALSE; if (!AliPythiaRndm::GetPythiaRandom()) AliPythiaRndm::SetPythiaRandom(GetRandom()); } AliGenPythia::AliGenPythia(Int_t npart) :AliGenMC(npart) { // default charm production at 5. 5 TeV // semimuonic decay // structure function GRVHO // fName = "Pythia"; fTitle= "Particle Generator using PYTHIA"; fXsection = 0.; fReadFromFile = 0; SetProcess(); SetStrucFunc(); SetForceDecay(); SetPtHard(); SetYHard(); SetEnergyCMS(); fDecayer = new AliDecayerPythia(); // Set random number generator if (!AliPythiaRndm::GetPythiaRandom()) AliPythiaRndm::SetPythiaRandom(GetRandom()); fFlavorSelect = 0; // Produced particles fParticles = new TClonesArray("TParticle",1000); fHeader = 0; SetEventListRange(); SetJetPhiRange(); SetJetEtaRange(); SetJetEtRange(); SetGammaPhiRange(); SetGammaEtaRange(); SetJetReconstructionMode(); SetQuench(); SetHadronisation(); SetPtKick(); // Options determining what to keep in the stack (Heavy flavour generation) fStackFillOpt = kFlavorSelection; // Keep particle with selected flavor fFeedDownOpt = kTRUE; // allow feed down from higher family // Fragmentation on/off fFragmentation = kTRUE; // Default counting mode fCountMode = kCountAll; // Pycel SetPycellParameters(); fSetNuclei = kFALSE; } AliGenPythia::AliGenPythia(const AliGenPythia & Pythia) :AliGenMC(Pythia) { // copy constructor Pythia.Copy(*this); } AliGenPythia::~AliGenPythia() { // Destructor } 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); // // Forward Paramters to the AliPythia object fDecayer->SetForceDecay(fForceDecay); fDecayer->Init(); 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); // 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); } else { fPythia->SetMSTP(91,0); } if (fReadFromFile) { fRL = AliRunLoader::Open(fFileName, "Partons"); fRL->LoadKinematics(); fRL->LoadHeader(); } else { fRL = 0x0; } // fPythia->ProcInit(fProcess,fEnergyCMS,fStrucFunc); // Parent and Children Selection switch (fProcess) { case kPyCharm: case kPyCharmUnforced: case kPyCharmPbPbMNR: case kPyCharmppMNR: case kPyCharmpPbMNR: fParentSelect[0] = 411; fParentSelect[1] = 421; fParentSelect[2] = 431; fParentSelect[3] = 4122; 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 kPyBeauty: case kPyBeautyPbPbMNR: case kPyBeautypPbMNR: case kPyBeautyppMNR: 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 kPyMb: case kPyMbNonDiffr: case kPyJets: case kPyDirectGamma: 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"); } if (fQuench) { fPythia->InitQuenching(0., 0.1, 0.6e6, 0); } } void AliGenPythia::Generate() { // Generate one event 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; // Set collision vertex position if (fVertexSmear == kPerEvent) Vertex(); // event loop while(1) { // // Produce event // // // Switch hadronisation off // fPythia->SetMSTJ(1, 0); // // Either produce new event or read partons from file // if (!fReadFromFile) { fPythia->Pyevnt(); fNpartons = fPythia->GetN(); } else { printf("Loading Event %d\n",AliRunLoader::GetRunLoader()->GetEventNumber()); fRL->GetEvent(AliRunLoader::GetRunLoader()->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.); } // // Switch hadronisation on // fPythia->SetMSTJ(1, 1); // // .. and perform hadronisation // printf("Calling hadronisation %d\n", fPythia->GetN()); fPythia->Pyexec(); if (gAlice) { if (gAlice->GetEvNumber()>=fDebugEventFirst && gAlice->GetEvNumber()<=fDebugEventLast) fPythia->Pylist(1); } fTrials++; fPythia->ImportParticles(fParticles,"All"); Boost(); // // // Int_t i; 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 != kPyMb && fProcess != kPyJets && fProcess != kPyDirectGamma && fProcess != kPyMbNonDiffr) { 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; // meson ? if (kfl > 10) kfl/=100; // baryon if (kfl > 10) kfl/=10; if (kfl > 10) kfl/=10; Int_t ipa = iparticle->GetFirstMother()-1; Int_t kfMo = 0; 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; iAt(i); kf = CheckPDGCode(iparticle->GetPdgCode()); 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); } // PushTrack loop } } else { nc = GenerateMB(); } // mb ? if (pParent) delete[] pParent; if (pSelected) delete[] pSelected; if (trackIt) 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); printf("\n Trials: %i %i %i\n",fTrials, fNpart, jev); 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) { TParticle* jet1 = (TParticle *) fParticles->At(6); TParticle* jet2 = (TParticle *) fParticles->At(7); if (!CheckTrigger(jet1, jet2)) return 0; } for (i = 0; iAt(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 && 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=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); // // Special Treatment to store color-flow // if (ks == 3 || ks == 13 || ks == 14) { TParticle* particle = 0; if (fStack) { particle = fStack->Particle(nt); } else { particle = gAlice->Stack()->Particle(nt); } particle->SetFirstDaughter(fPythia->GetK(2, i)); particle->SetLastDaughter(fPythia->GetK(3, i)); } KeepTrack(nt); pParent[i] = nt; } // select particle } // particle loop if (pParent) delete[] pParent; printf("\n I've put %i particles on the stack \n",nc); return nc; } void AliGenPythia::FinishRun() { // Print x-section summary fPythia->Pystat(1); 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() { // 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) { // Treat protons as inside nuclei with mass numbers a1 and a2 fAProjectile = a1; fATarget = a2; fSetNuclei = kTRUE; } void AliGenPythia::MakeHeader() { // 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); // // Jets that have triggered if (fProcess == kPyJets) { 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; fPythia->GetQuenchingParameters(xp, yp, z); ((AliGenPythiaEventHeader*) fHeader)->SetXYJet(xp, yp); ((AliGenPythiaEventHeader*) fHeader)->SetZQuench(z); } // // Pass header to RunLoader // AliRunLoader::GetRunLoader()->GetHeader()->SetGenEventHeader(fHeader); } 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) { 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) 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; } AliGenPythia& AliGenPythia::operator=(const AliGenPythia& rhs) { // Assignment operator rhs.Copy(*this); return *this; } 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::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); } } } #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