]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - PYTHIA6/AliGenPythia.h
1. Adding the new histograms TPC z vertex correlation in order
[u/mrichter/AliRoot.git] / PYTHIA6 / AliGenPythia.h
index e4dfabb3f963c5fc289d272b9da123600b37000b..6cbdb0377de89dd77289b6fbb48a3197d04abf7f 100644 (file)
@@ -25,6 +25,7 @@ class AliGenPythiaEventHeader;
 class AliGenEventHeader;
 class AliStack;
 class AliRunLoader;
+class TObjArray; 
 
 class AliGenPythia : public AliGenMC
 {
@@ -39,87 +40,120 @@ class AliGenPythia : public AliGenMC
     virtual ~AliGenPythia();
     virtual void    Generate();
     virtual void    Init();
-    // set a cut on the Z coord. of the primary vertex (cm)
-    //
+    // Range of events to be printed
     virtual void    SetEventListRange(Int_t eventFirst=-1, Int_t eventLast=-1);
-    // select process type
+    // Select process type
     virtual void    SetProcess(Process_t proc = kPyCharm) {fProcess = proc;}
-    // select structure function
+    virtual void    SetTune(Int_t itune) {fItune = itune;}
+
+    // Select structure function
     virtual void    SetStrucFunc(StrucFunc_t func =  kCTEQ5L) {fStrucFunc = func;}
-    // select pt of hard scattering 
+    // Select pt of hard scattering 
     virtual void    SetPtHard(Float_t ptmin = 0, Float_t ptmax = 1.e10)
        {fPtHardMin = ptmin; fPtHardMax = ptmax; }
+    // y of hard scattering
     virtual void    SetYHard(Float_t ymin = -1.e10, Float_t ymax = 1.e10)
        {fYHardMin = ymin; fYHardMax = ymax; }
     // Set initial and final state gluon radiation
     virtual void    SetGluonRadiation(Int_t iIn, Int_t iFin)
        {fGinit = iIn; fGfinal = iFin;}
+    // Intrinsic kT
     virtual void    SetPtKick(Float_t kt = 1.)
        {fPtKick = kt;}
     // Use the Pythia 6.3 new multiple interations scenario
     virtual void    UseNewMultipleInteractionsScenario() {fNewMIS = kTRUE;}
     // Switch off heavy flavors
     virtual void    SwitchHFOff() {fHFoff = kTRUE;}
-    // set centre of mass energy
+    // Set centre of mass energy
     virtual void    SetEnergyCMS(Float_t energy = 5500) {fEnergyCMS = energy;}
-    // treat protons as inside nuclei
-    virtual void    SetNuclei(Int_t a1, Int_t a2);
+    // Treat protons as inside nuclei with mass numbers a1 and a2
+    virtual void    SetNuclei(Int_t a1, Int_t a2, Int_t pdfset = 0);
+    //
+    // Trigger options
+    //
+    // Energy range for jet trigger
     virtual void    SetJetEtRange(Float_t etmin = 0., Float_t etmax = 1.e4)
        {fEtMinJet = etmin; fEtMaxJet = etmax;}
+    // Eta range for jet trigger
     virtual void    SetJetEtaRange(Float_t etamin = -20., Float_t etamax = 20.)
        {fEtaMinJet = etamin; fEtaMaxJet = etamax;}
-    virtual void    SetJetReconstructionMode(Int_t mode = kCell) {fJetReconstruction = mode;}
+    // Phi range for jet trigger
     virtual void    SetJetPhiRange(Float_t phimin = 0., Float_t phimax = 360.)
        {fPhiMinJet = TMath::Pi()*phimin/180.; fPhiMaxJet = TMath::Pi()*phimax/180.;}
+    // Jet reconstruction mode; default is cone algorithm
+    virtual void    SetJetReconstructionMode(Int_t mode = kCell) {fJetReconstruction = mode;}
+    // Eta range for gamma trigger 
     virtual void    SetGammaEtaRange(Float_t etamin = -20., Float_t etamax = 20.)
        {fEtaMinGamma = etamin; fEtaMaxGamma = etamax;}
+    // Phi range for gamma trigger
     virtual void    SetGammaPhiRange(Float_t phimin = 0., Float_t phimax = 360.)
        {fPhiMinGamma = TMath::Pi()*phimin/180.; fPhiMaxGamma = TMath::Pi()*phimax/180.;}
-
    // Select jets with fragmentation photon or pi0 going to PHOS or EMCAL
+    virtual void  SetFragPhotonInCalo(Bool_t b)  {fFragPhotonInCalo = b;}
+    virtual void  SetPi0InCalo       (Bool_t b)  {fPi0InCalo    = b;}
+    virtual void  SetPhotonInCalo(Bool_t b)      {fPhotonInCalo = b;}
+    virtual void  SetCheckPHOS (Bool_t b)        {fCheckPHOS    = b;}
+    virtual void  SetCheckEMCAL(Bool_t b)        {fCheckEMCAL   = b;}
+    virtual void  SetFragPhotonInEMCAL(Bool_t b) {fCheckEMCAL   = b; fFragPhotonInCalo = b;}
+    virtual void  SetFragPhotonInPHOS(Bool_t b)  {fCheckPHOS    = b; fFragPhotonInCalo = b;}
+    virtual void  SetPi0InEMCAL(Bool_t b)        {fCheckEMCAL   = b; fPi0InCalo        = b;}
+    virtual void  SetPi0InPHOS(Bool_t b)         {fCheckPHOS    = b; fPi0InCalo        = b;}
+    virtual void  SetPhotonInEMCAL(Bool_t b)     {fCheckEMCAL   = b; fPhotonInCalo     = b;}
+    virtual void  SetElectronInEMCAL(Bool_t b)   {fEleInEMCAL   = b;}
+    virtual void  SetPhotonInPHOS(Bool_t b)      {fCheckPHOS    = b; fPhotonInCalo     = b;}
 
-    virtual void  SetFragPhotonInCalo(Bool_t b) {fFragPhotonInCalo = b;}
-    virtual void  SetPi0InCalo       (Bool_t b)      {fPi0InCalo = b;}
-    virtual void  SetCheckPHOS (Bool_t b)        {fCheckPHOS = b;}
-    virtual void  SetCheckEMCAL(Bool_t b)       {fCheckEMCAL = b;}
-
-    virtual void  SetFragPhotonInEMCAL(Bool_t b) 
-    {fCheckEMCAL = b;  fFragPhotonInCalo = b ;}
-    virtual void  SetFragPhotonInPHOS(Bool_t b) 
-    {fCheckPHOS = b ;   fFragPhotonInCalo = b ;}
-    virtual void  SetPi0InEMCAL(Bool_t b) 
-    {fCheckEMCAL = b ;  fPi0InCalo = b ;}
-    virtual void  SetPi0InPHOS(Bool_t b) 
-    {fCheckPHOS = b ;   fPi0InCalo = b ;}
-    
-    virtual void  SetFragPhotonOrPi0MinPt(Float_t pt)      {fFragPhotonOrPi0MinPt=pt;}
-  
-    // Trigger on a particle
+    // Trigger on a minimum multiplicity
+    virtual void  SetTriggerChargedMultiplicity(Int_t multiplicity, Float_t etamax = 0, Float_t ptmin = -1.) 
+    {fTriggerMultiplicity = multiplicity; fTriggerMultiplicityEta = etamax; 
+      fTriggerMultiplicityPtMin = ptmin;}
+       
+    virtual void  SetPhotonInPHOSeta(Bool_t b)   {fCheckPHOSeta = b; fPhotonInCalo     = b;}
+    virtual void  SetFragPhotonOrPi0MinPt(Float_t pt)      {fFragPhotonOrPi0MinPt = pt;}
+    virtual void  SetPhotonMinPt(Float_t pt)     {fPhotonMinPt = pt;}
+    virtual void  SetElectronMinPt(Float_t pt)     {fElectronMinPt = pt;}
+    // Trigger and rotate event 
+    void RotatePhi(Int_t iphcand, Bool_t& okdd);
+    // Trigger on a single particle
     virtual void    SetTriggerParticle(Int_t particle = 0, Float_t etamax = 0.9) 
        {fTriggerParticle = particle; fTriggerEta = etamax;}
+    //
+    // Heavy flavor options
+    //
     // Set option for feed down from higher family
     virtual void SetFeedDownHigherFamily(Bool_t opt) {
-      fFeedDownOpt = opt;
+       fFeedDownOpt = opt;
     }
     // Set option for selecting particles kept in stack according to flavor
     // or to parent selection
     virtual void SetStackFillOpt(StackFillOpt_t opt) {
-      fStackFillOpt = opt;
+       fStackFillOpt = opt;
     }
     // Set fragmentation option
     virtual void SetFragmentation(Bool_t opt) {
-      fFragmentation = opt;
+       fFragmentation = opt;
     }
     // Set counting mode
     virtual void SetCountMode(CountMode_t mode) {
-      fCountMode = mode;
+       fCountMode = mode;
     }
-    // Set quenching mode 0 = no, 1 = AM, 2 = IL
+    //
+    // Quenching
+    //
+    // Set quenching mode 0 = no, 1 = AM, 2 = IL,  3 = NA, 4 = ACS
     virtual void SetQuench(Int_t flag = 0) {fQuench = flag;}
+    // Set transport coefficient.
+    void SetQhat(Float_t qhat) {fQhat = qhat;}
+    //Set initial medium length.
+    void SetLength(Float_t length) {fLength = length;}
+
     virtual void SetHadronisation(Int_t flag = 1) {fHadronisation = flag;}
-    virtual void SetReadFromFile(const Text_t *filname) {fFileName = filname;  fReadFromFile = 1;}    
-  // Get interaction rate for pileup studies
+    virtual void SetPatchOmegaDalitz(Int_t flag = 1) {fPatchOmegaDalitz = flag;}
+    virtual void SetReadFromFile(const Text_t *filname) {fkFileName = filname;  fReadFromFile = 1;}    
+
+    //
+    // Pile-up
+    //
+    // Get interaction rate for pileup studies
     virtual void    SetInteractionRate(Float_t rate,Float_t timewindow = 90.e-6);
     virtual Float_t GetInteractionRate() const {return fInteractionRate;}
     // get cross section of process
@@ -131,13 +165,13 @@ class AliGenPythia : public AliGenMC
                             Float_t thresh = 0., Float_t etseed = 4.,
                             Float_t minet = 10., Float_t r = 1.);
     
-    void LoadEvent(AliStack* stack, Int_t flag = 0, Int_t reHadr = 0);
+  void LoadEvent(AliStack* stack, Int_t flag = 0, Int_t reHadr = 0);
+  void LoadEvent(TObjArray* stack, Int_t flag = 0, Int_t reHadr = 0);
     // Getters
     virtual Process_t    GetProcess() const {return fProcess;}
     virtual StrucFunc_t  GetStrucFunc() const {return fStrucFunc;}
     virtual void         GetPtHard(Float_t& ptmin, Float_t& ptmax) const
        {ptmin = fPtHardMin; ptmax = fPtHardMax;}
-    virtual Float_t      GetEnergyCMS() const {return fEnergyCMS;}
     virtual void         GetNuclei(Int_t&  a1, Int_t& a2) const
        {a1 = fAProjectile; a2 = fATarget;}
     virtual void         GetJetEtRange(Float_t& etamin, Float_t& etamax) const
@@ -154,11 +188,8 @@ class AliGenPythia : public AliGenMC
     //
     virtual void FinishRun();
     Bool_t CheckTrigger(TParticle* jet1, TParticle* jet2);
-
     //Used in some processes to selected child properties
     Bool_t CheckKinematicsOnChild();
-    
-
     void     GetSubEventTime();
 
  protected:
@@ -168,8 +199,8 @@ class AliGenPythia : public AliGenMC
     void     MakeHeader();    
     void     GeneratePileup();
     Process_t   fProcess;           //Process type
+    Int_t       fItune;             // Pythia tune > 6.4
     StrucFunc_t fStrucFunc;         //Structure Function
-    Float_t     fEnergyCMS;         //Centre of mass energy
     Float_t     fKineBias;          //!Bias from kinematic selection
     Int_t       fTrials;            //!Number of trials for current event
     Int_t       fTrialsRun;         //!Number of trials for run
@@ -192,9 +223,13 @@ class AliGenPythia : public AliGenMC
     Int_t       fGinit;             //initial state gluon radiation
     Int_t       fGfinal;            //final state gluon radiation
     Int_t       fHadronisation;     //hadronisation
+    Bool_t      fPatchOmegaDalitz;  //flag for omega dalitz decay patch
     Int_t       fNpartons;          //Number of partons before hadronisation
     Int_t       fReadFromFile;      //read partons from file
     Int_t       fQuench;            //Flag for quenching
+    Float_t     fQhat;              //Transport coefficient (GeV^2/fm)
+    Float_t     fLength;            //Medium length (fm)
+    Float_t     fImpact;            //Impact parameter for quenching simulation (q-pythia) 
     Float_t     fPtKick;            //Transverse momentum kick
     Bool_t      fFullEvent;         //!Write Full event if true
     AliDecayer  *fDecayer;          //!Pointer to the decayer instance
@@ -223,43 +258,52 @@ class AliGenPythia : public AliGenMC
                                     // parents and their decays
     Bool_t fFeedDownOpt;            // Option to set feed down from higher
                                     // quark families (e.g. b->c)
-    Bool_t  fFragmentation;          // Option to activate fragmentation by Pythia
-    Bool_t  fSetNuclei;              // Flag indicating that SetNuclei has been called
-    Bool_t  fNewMIS;                 // Flag for the new multipple interactions scenario
-    Bool_t  fHFoff;                  // Flag for switching heafy flavor production off
-    Int_t   fTriggerParticle;        // Trigger on this particle ...
-    Float_t fTriggerEta;             // .. within |eta| < fTriggerEta
-    //
-    CountMode_t fCountMode;            // Options for counting when the event will be finished.
-    AliGenPythiaEventHeader* fHeader;  //! Event header
-    AliRunLoader*            fRL;      //! Run Loader
-    const Text_t* fFileName;           //! Name of file to read from
-     
+    Bool_t  fFragmentation;         // Option to activate fragmentation by Pythia
+    Bool_t  fSetNuclei;             // Flag indicating that SetNuclei has been called
+    Bool_t  fNewMIS;                // Flag for the new multipple interactions scenario
+    Bool_t  fHFoff;                 // Flag for switching heafy flavor production off
+    Int_t   fNucPdf;                // Nuclear pdf 0: EKS98 1: EPS08
+    Int_t   fTriggerParticle;       // Trigger on this particle ...
+    Float_t fTriggerEta;            // .. within |eta| < fTriggerEta
+    Int_t       fTriggerMultiplicity;       // Trigger on events with a minimum charged multiplicity
+    Float_t     fTriggerMultiplicityEta;    // in a given eta range
+    Float_t     fTriggerMultiplicityPtMin;  // above this pT 
+    CountMode_t fCountMode;         // Options for counting when the event will be finished.     
     // fCountMode = kCountAll         --> All particles that end up in the
     //                                    stack are counted
     // fCountMode = kCountParents     --> Only selected parents are counted
     // fCountMode = kCountTrackabless --> Only particles flagged for tracking
     //                                     are counted
     //
+    //
+
+    AliGenPythiaEventHeader* fHeader;  //! Event header
+    AliRunLoader*            fRL;      //! Run Loader
+    const Text_t* fkFileName;          //! Name of file to read from
 
-    Bool_t fFragPhotonInCalo ; // Option ask for Fragmentation Photon in calorimeters acceptance
-    Bool_t fPi0InCalo ; // Option ask for Pi0 in calorimeters acceptance
-    Bool_t fCheckEMCAL ; //  Option ask for FragPhoton or Pi0 in calorimeters EMCAL acceptance
-    Bool_t fCheckPHOS ; //  Option ask for FragPhoton or Pi0 in calorimeters PHOS acceptance
-    Float_t fFragPhotonOrPi0MinPt ; // Minimum momentum of Fragmentation Photon or Pi0 
-    //Calorimeters eta-phi acceptance 
-    Float_t fPHOSMinPhi ;
-    Float_t fPHOSMaxPhi ;
-    Float_t fPHOSEta ;
-    Float_t fEMCALMinPhi ;
-    Float_t fEMCALMaxPhi ;
-    Float_t fEMCALEta ;
 
+    Bool_t fFragPhotonInCalo; // Option to ask for Fragmentation Photon in calorimeters acceptance
+    Bool_t fPi0InCalo;        // Option to ask for Pi0 in calorimeters acceptance
+    Bool_t fPhotonInCalo;     // Option to ask for Decay Photon in calorimeter acceptance
+    Bool_t fEleInEMCAL;       // Option to ask for Electron in EMCAL acceptance
+    Bool_t fCheckEMCAL;       // Option to ask for FragPhoton or Pi0 in calorimeters EMCAL acceptance
+    Bool_t fCheckPHOS;        // Option to ask for FragPhoton or Pi0 in calorimeters PHOS acceptance
+    Bool_t fCheckPHOSeta;     // Option to ask for PHOS eta acceptance
+    Float_t fFragPhotonOrPi0MinPt; // Minimum momentum of Fragmentation Photon or Pi0
+    Float_t fPhotonMinPt;          // Minimum momentum of Photon 
+    Float_t fElectronMinPt;        // Minimum momentum of Electron 
+    //Calorimeters eta-phi acceptance 
+    Float_t fPHOSMinPhi;           // Minimum phi PHOS
+    Float_t fPHOSMaxPhi;           // Maximum phi PHOS
+    Float_t fPHOSEta;              // Minimum eta PHOS
+    Float_t fEMCALMinPhi;          // Minimum phi EMCAL
+    Float_t fEMCALMaxPhi;          // Maximum phi EMCAL
+    Float_t fEMCALEta;             // Maximum eta EMCAL
  private:
     AliGenPythia(const AliGenPythia &Pythia);
     AliGenPythia & operator=(const AliGenPythia & rhs);
 
-    ClassDef(AliGenPythia,5) // AliGenerator interface to Pythia
+    ClassDef(AliGenPythia, 10) // AliGenerator interface to Pythia
 };
 #endif