]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - PYTHIA6/AliGenPythia.cxx
Additional software trigger posiibilities for pions and gammas (G. Conesa)
[u/mrichter/AliRoot.git] / PYTHIA6 / AliGenPythia.cxx
index 93cdb30fa2ae31e0883dedad3d03ac9c129f0738..3a3fd632bd1872f02b4ede8c78c1317fe1a7714b 100644 (file)
@@ -30,7 +30,6 @@
 #include <TPDGCode.h>
 #include <TSystem.h>
 #include <TTree.h>
-
 #include "AliConst.h"
 #include "AliDecayerPythia.h"
 #include "AliGenPythia.h"
 #include "AliStack.h"
 #include "AliRunLoader.h"
 #include "AliMC.h"
+#include "pyquenCommon.h"
 
 ClassImp(AliGenPythia)
 
-AliGenPythia::AliGenPythia()
-                 :AliGenMC()
+
+AliGenPythia::AliGenPythia():
+    AliGenMC(),
+    fProcess(kPyCharm),          
+    fStrucFunc(kCTEQ5L), 
+    fEnergyCMS(5500.),
+    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),
+    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),    
+    fTriggerParticle(0),
+    fTriggerEta(0.9),     
+    fCountMode(kCountAll),      
+    fHeader(0),  
+    fRL(0),      
+    fFileName(0),
+    fBremssInCalo(kFALSE),
+    fPi0InCalo(kFALSE) ,
+    fBremssPi0InEMCAL(kFALSE),
+    fBremssPi0InPHOS(kFALSE),
+    fBremssPi0MinPt(0), 
+    fPHOSMinPhi(220.),
+    fPHOSMaxPhi(320.),
+    fPHOSEta(0.12),
+    fEMCALMinPhi(80.),
+    fEMCALMaxPhi(190.),
+    fEMCALEta(0.7)
+
 {
 // Default Constructor
-  fParticles = 0;
-  fPythia    = 0;
-  fHeader = 0;
-  fReadFromFile = 0;
-  fDecayer   = new AliDecayerPythia();
-  SetEventListRange();
-  SetJetPhiRange();
-  SetJetEtaRange();
-  SetJetEtRange();
-  SetGammaPhiRange();
-  SetGammaEtaRange();
-  SetPtKick();
-  SetQuench();
-  SetHadronisation();  
-  fSetNuclei = kFALSE;
+  SetNuclei(0,0);
   if (!AliPythiaRndm::GetPythiaRandom()) 
-    AliPythiaRndm::SetPythiaRandom(GetRandom());
+      AliPythiaRndm::SetPythiaRandom(GetRandom());
 }
 
 AliGenPythia::AliGenPythia(Int_t npart)
-                 :AliGenMC(npart)
+    :AliGenMC(npart),
+     fProcess(kPyCharm),          
+     fStrucFunc(kCTEQ5L), 
+     fEnergyCMS(5500.),
+     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),
+     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),    
+     fTriggerParticle(0),
+     fTriggerEta(0.9),     
+     fCountMode(kCountAll),      
+     fHeader(0),  
+     fRL(0),      
+     fFileName(0),
+     fBremssInCalo(kFALSE),
+     fPi0InCalo(kFALSE) ,
+     fBremssPi0InEMCAL(kFALSE),
+     fBremssPi0InPHOS(kFALSE),
+     fBremssPi0MinPt(0),
+     fPHOSMinPhi(220.),
+     fPHOSMaxPhi(320.),
+     fPHOSEta(0.12),
+     fEMCALMinPhi(80.),
+     fEMCALMaxPhi(190.),
+     fEMCALEta(0.7)
 {
 // default charm production at 5. 5 TeV
 // semimuonic decay
@@ -77,46 +213,90 @@ AliGenPythia::AliGenPythia(Int_t npart)
 //
     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;
-}
+    SetNuclei(0,0);
+ }
 
 AliGenPythia::AliGenPythia(const AliGenPythia & Pythia)
-    :AliGenMC(Pythia)
+    :AliGenMC(Pythia),
+     fProcess(kPyCharm),          
+     fStrucFunc(kCTEQ5L), 
+     fEnergyCMS(5500.),
+     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),
+     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),    
+     fTriggerParticle(0),
+     fTriggerEta(0.9),     
+     fCountMode(kCountAll),      
+     fHeader(0),  
+     fRL(0),      
+     fFileName(0),
+     fBremssInCalo(kFALSE),
+     fPi0InCalo(kFALSE) ,
+     fBremssPi0InEMCAL(kFALSE),
+     fBremssPi0InPHOS(kFALSE),
+     fBremssPi0MinPt(0),
+     fPHOSMinPhi(220.),
+     fPHOSMaxPhi(320.),
+     fPHOSEta(0.12),
+     fEMCALMinPhi(80.),
+     fEMCALMaxPhi(190.),
+     fEMCALEta(0.7)   
 {
 // copy constructor
     Pythia.Copy(*this);
@@ -125,6 +305,57 @@ AliGenPythia::AliGenPythia(const AliGenPythia & Pythia)
 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,
@@ -200,19 +431,28 @@ void AliGenPythia::Init()
     } else {
        fRL = 0x0;
     }
-    
-    
+// Switch off Heavy Flavors on request  
+    if (fHFoff) {
+       fPythia->SetMSTP(58, 3);
+       fPythia->SetMSTJ(45, 3);        
+       for (Int_t i = 156; i <= 160; i++) fPythia->SetMDME(i, 1, 0);
+    }
  //
     fPythia->ProcInit(fProcess,fEnergyCMS,fStrucFunc);
 
 //  Parent and Children Selection
     switch (fProcess) 
     {
+    case kPyOldUEQ2ordered:
+    case kPyOldUEQ2ordered2:
+    case kPyOldPopcorn:
+      break;
     case kPyCharm:
     case kPyCharmUnforced:
     case kPyCharmPbPbMNR:
-    case kPyCharmppMNR:
     case kPyCharmpPbMNR:
+    case kPyCharmppMNR:
+    case kPyCharmppMNRwmi:
        fParentSelect[0] =   411;
        fParentSelect[1] =   421;
        fParentSelect[2] =   431;
@@ -231,10 +471,17 @@ void AliGenPythia::Init()
        fParentSelect[0] =   411;
        fFlavorSelect    =   4; 
        break;
+    case kPyDPlusStrangePbPbMNR:
+    case kPyDPlusStrangepPbMNR:
+    case kPyDPlusStrangeppMNR:
+       fParentSelect[0] =   431;
+       fFlavorSelect    =   4; 
+       break;
     case kPyBeauty:
     case kPyBeautyPbPbMNR:
     case kPyBeautypPbMNR:
     case kPyBeautyppMNR:
+    case kPyBeautyppMNRwmi:
        fParentSelect[0]=  511;
        fParentSelect[1]=  521;
        fParentSelect[2]=  531;
@@ -260,9 +507,14 @@ void AliGenPythia::Init()
        break;
     case kPyMb:
     case kPyMbNonDiffr:
+    case kPyMbMSEL1:
     case kPyJets:
     case kPyDirectGamma:
+    case kPyLhwgMb:    
        break;
+    case kPyW:
+    case kPyZ:
+        break;
     }
 //
 //
@@ -270,17 +522,17 @@ void AliGenPythia::Init()
 //
 //  Configure detector (EMCAL like)
 //
-       fPythia->SetPARU(51, fPycellEtaMax);
-       fPythia->SetMSTU(51, fPycellNEta);
-       fPythia->SetMSTU(52, fPycellNPhi);
+    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);
+    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;
@@ -299,11 +551,22 @@ void AliGenPythia::Init()
        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);
     }
-    
+    fPythia->SetPARJ(200, 0.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
+    }
 }
 
 void AliGenPythia::Generate()
@@ -322,6 +585,8 @@ void AliGenPythia::Generate()
     Int_t jev=0;
     Int_t j, kf;
     fTrials=0;
+    fEventTime = 0.;
+    
     
 
     //  Set collision vertex position 
@@ -341,7 +606,11 @@ void AliGenPythia::Generate()
 // Either produce new event or read partons from file
 //     
        if (!fReadFromFile) {
-           fPythia->Pyevnt();
+           if (!fNewMIS) {
+               fPythia->Pyevnt();
+           } else {
+               fPythia->Pyevnw();
+           }
            fNpartons = fPythia->GetN();
        } else {
            printf("Loading Event %d\n",AliRunLoader::GetRunLoader()->GetEventNumber());
@@ -358,7 +627,10 @@ void AliGenPythia::Generate()
            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
 //
@@ -367,12 +639,6 @@ void AliGenPythia::Generate()
 // .. 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();
@@ -384,7 +650,7 @@ void AliGenPythia::Generate()
 
        Int_t np = fParticles->GetEntriesFast();
        
-       if (np == 0 ) continue;
+       if (np == 0) continue;
 //
        
 //
@@ -402,7 +668,10 @@ void AliGenPythia::Generate()
        Int_t nTkbles = 0;   // Trackable particles
        if (fProcess != kPyMb && fProcess != kPyJets && 
            fProcess != kPyDirectGamma &&
-           fProcess != kPyMbNonDiffr) {
+           fProcess != kPyMbNonDiffr  &&
+           fProcess != kPyMbMSEL1     &&
+           fProcess != kPyW && fProcess != kPyZ &&
+           fProcess != kPyCharmppMNRwmi && fProcess != kPyBeautyppMNRwmi) {
            
            for (i = 0; i < np; i++) {
                TParticle* iparticle = (TParticle *) fParticles->At(i);
@@ -416,19 +685,39 @@ void AliGenPythia::Generate()
                // 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;
-               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 j = jmin; j <= jmax; j++)
+                               ((TParticle *) fParticles->At(j))->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;
@@ -456,7 +745,7 @@ void AliGenPythia::Generate()
 // Kinematic seletion on final state heavy flavor mesons
                    if (ParentSelected(kf) && !KinematicSelection(iparticle, 0)) 
                    {
-                     continue;
+                       continue;
                    }
                    pSelected[i] = 1;
                    if (ParentSelected(kf)) ++nParents; // Update parent count
@@ -464,9 +753,9 @@ void AliGenPythia::Generate()
                } else {
 // Kinematic seletion on decay products
                    if (fCutOnChild && ParentSelected(kfMo) && ChildSelected(kf) 
-                       && !KinematicSelection(iparticle, 1))
+                       && !KinematicSelection(iparticle, 1)) 
                    {
-                     continue;
+                       continue;
                    }
 //
 // Decay products 
@@ -496,6 +785,7 @@ void AliGenPythia::Generate()
                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;
@@ -546,10 +836,12 @@ void AliGenPythia::Generate()
        } else {
            nc = GenerateMB();
        } // mb ?
+       
+       GetSubEventTime();
 
-       if (pParent)   delete[] pParent;
-       if (pSelected) delete[] pSelected;
-       if (trackIt)   delete[] trackIt;
+       delete[] pParent;
+       delete[] pSelected;
+       delete[] trackIt;
 
        if (nc > 0) {
          switch (fCountMode) {
@@ -568,8 +860,7 @@ void AliGenPythia::Generate()
          }
            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);
@@ -605,15 +896,93 @@ Int_t  AliGenPythia::GenerateMB()
     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;
+       if (!CheckTrigger(jet1, jet2)) {
+         delete [] pParent;
+         return 0;
+       }
+    }
+
+    // Select jets with bremsstrahlung or pi0 going to PHOS or EMCAL
+    if (fProcess == kPyJets && (fBremssInCalo || fPi0InCalo) ) {
+
+      Bool_t ok = kFALSE;
+
+      Int_t pdg  = 0; 
+      if (fBremssInCalo) 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() > fBremssPi0MinPt){
+         Float_t phi = iparticle->Phi()*180./TMath::Pi(); //Convert to degrees
+         Float_t eta =TMath::Abs(iparticle->Eta());//in calos etamin=-etamax     
+         if((fBremssPi0InEMCAL && IsInEMCAL(phi,eta)) ||
+             fBremssPi0InPHOS    && IsInPHOS(phi,eta)) 
+           ok =kTRUE;  
+       }
+      }
+      if(!ok)
+       return 0;
     }
     
-    for (i = 0; i<np; i++) {
+    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) {
+      TParticle *hvq;
+      Bool_t  theQ=kFALSE,theQbar=kFALSE,inYcut=kFALSE;
+      Float_t yQ;  
+      Int_t   pdgQ;
+      for(i=0; i<np; i++) {
+       hvq = (TParticle*)fParticles->At(i);
+       pdgQ = hvq->GetPdgCode();  
+       if(TMath::Abs(pdgQ) != fFlavorSelect) continue; 
+       if(pdgQ>0) { theQ=kTRUE; } else { theQbar=kTRUE; }
+       yQ = 0.5*TMath::Log((hvq->Energy()+hvq->Pz()+1.e-13)/
+                           (hvq->Energy()-hvq->Pz()+1.e-13));
+       if(yQ>fYMin && yQ<fYMax) inYcut=kTRUE;
+      }
+      if (!theQ || !theQbar || !inYcut) {
+       delete[] pParent;
+       return 0;
+      }
+    }
+
+    //Introducing child cuts in case kPyW, kPyZ, kPyMb, and kPyMbNonDiff
+    if ( (fProcess == kPyW || fProcess == kPyZ || fProcess == kPyMb || 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());
@@ -640,7 +1009,7 @@ Int_t  AliGenPythia::GenerateMB()
            origin[1] = fVertex[1]+iparticle->Vy()/10; // [cm]
            origin[2] = fVertex[2]+iparticle->Vz()/10; // [cm]
            
-           Float_t tof=kconv*iparticle->T();
+           Float_t tof = fEventTime + kconv * iparticle->T();
 
            PushTrack(fTrackIt*trackIt, iparent, kf, 
                      p[0], p[1], p[2], p[3], 
@@ -663,13 +1032,14 @@ Int_t  AliGenPythia::GenerateMB()
            
            KeepTrack(nt);
            pParent[i] = nt;
+           SetHighWaterMark(nt);
+           
        } // select particle
     } // particle loop 
 
-    if (pParent) delete[] pParent;
+    delete[] pParent;
     
-    printf("\n I've put %i particles on the stack \n",nc);
-    return nc;
+    return 1;
 }
 
 
@@ -677,16 +1047,18 @@ void AliGenPythia::FinishRun()
 {
 // Print x-section summary
     fPythia->Pystat(1);
-    fQ  /= fNev;
-    fX1 /= fNev;
-    fX2 /= fNev;    
+
+    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()
+void AliGenPythia::AdjustWeights() const
 {
 // Adjust the weights after generation of all events
 //
@@ -712,6 +1084,14 @@ void AliGenPythia::SetNuclei(Int_t a1, Int_t a2)
 
 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");
@@ -726,11 +1106,13 @@ void AliGenPythia::MakeHeader()
     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], 
@@ -759,19 +1141,42 @@ void AliGenPythia::MakeHeader()
     if (fQuench){
        Double_t z[4];
        Double_t xp, yp;
-       
-       fPythia->GetQuenchingParameters(xp, yp, z);
-       
-       ((AliGenPythiaEventHeader*) fHeader)->SetXYJet(xp, yp);
-       ((AliGenPythiaEventHeader*) fHeader)->SetZQuench(z);
-    }
-    
+       if (fQuench == 1) {
+           // Pythia::Quench()
+           fPythia->GetQuenchingParameters(xp, yp, z);
+       } else {
+           // 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);
+           
+       }
+           ((AliGenPythiaEventHeader*) fHeader)->SetXYJet(xp, yp);
+           ((AliGenPythiaEventHeader*) fHeader)->SetZQuench(z);
+       }
+//
+// Store pt^hard 
+    ((AliGenPythiaEventHeader*) fHeader)->SetPtHard(fPythia->GetVINT(47));
 //
-//  Pass header to RunLoader
+//  Pass header
 //
-    AliRunLoader::GetRunLoader()->GetHeader()->SetGenEventHeader(fHeader);   
+    AddHeader(fHeader);
 }
-       
+
+void AliGenPythia::AddHeader(AliGenEventHeader* header)
+{
+    // Add header to container or runloader
+    if (fContainer) {
+       fContainer->AddHeader(header);
+    } else {
+       AliRunLoader::GetRunLoader()->GetHeader()->SetGenEventHeader(header);   
+    }
+}
+
 
 Bool_t AliGenPythia::CheckTrigger(TParticle* jet1, TParticle* jet2)
 {
@@ -797,7 +1202,7 @@ Bool_t AliGenPythia::CheckTrigger(TParticle* jet1, TParticle* jet2)
 //
        GetJets(njets, ntrig, jets);
        
-       if (ntrig) triggered = kTRUE;
+       if (ntrig || fEtMinJet == 0.) triggered = kTRUE;
 //
     } else {
        Int_t ij = 0;
@@ -820,6 +1225,38 @@ Bool_t AliGenPythia::CheckTrigger(TParticle* jet1, TParticle* jet2)
     }
     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; j<npart; j++) {
+       TParticle *  jparticle = (TParticle *) fParticles->At(j);
+       kcode = TMath::Abs( CheckPDGCode(jparticle->GetPdgCode()) );
+       ks = jparticle->GetStatusCode();
+       km = jparticle->GetFirstMother(); 
+       
+       if( (ks == 1)  &&  (kcode == fPdgCodeParticleforAcceptanceCut)  &&  (KinematicSelection(jparticle,1)) ){
+           nPartAcc++;
+       }
+       if( numberOfAcceptedParticles <= nPartAcc){
+         checking = kTRUE;
+         break;
+       }
+    }
+
+    return checking;
+}
+
          
 AliGenPythia& AliGenPythia::operator=(const  AliGenPythia& rhs)
 {
@@ -846,12 +1283,12 @@ void  AliGenPythia::LoadEvent(AliStack* stack, Int_t flag, Int_t reHadr)
     
     
     for (Int_t part = 0; part < npart; part++) {
-       TParticle *MPart = stack->Particle(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();
+       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) {
@@ -861,11 +1298,11 @@ void  AliGenPythia::LoadEvent(AliStack* stack, Int_t flag, Int_t reHadr)
            }
        }
        
-       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();
+       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;
@@ -878,7 +1315,7 @@ void  AliGenPythia::LoadEvent(AliStack* stack, Int_t flag, Int_t reHadr)
        (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;
+       (fPythia->GetPyjets())->K[2][part+n0] = mPart->GetFirstMother() + 1;
     }
 }
 
@@ -951,7 +1388,6 @@ void  AliGenPythia::GetJets(Int_t& nJets, Int_t& nJetsTrig, Float_t jets[4][10])
        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 &&
@@ -970,6 +1406,43 @@ void  AliGenPythia::GetJets(Int_t& nJets, Int_t& nJetsTrig, Float_t jets[4][10])
     }
 }
 
+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? 
+  // Acceptance slightly larger considered.
+  // phi in degrees, etamin=-etamax
+  if(phi > fEMCALMinPhi-0.1  && phi < fEMCALMaxPhi+0.1 && 
+     eta < fEMCALEta+0.01   ) 
+    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-0.1  && phi < fPHOSMaxPhi+0.1 && 
+     eta < fPHOSEta+0.01   ) 
+    return kTRUE;
+  else 
+    return kFALSE;
+}
+
+
 
 #ifdef never
 void AliGenPythia::Streamer(TBuffer &R__b)
@@ -1013,3 +1486,4 @@ void AliGenPythia::Streamer(TBuffer &R__b)
 #endif
 
 
+