]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PYTHIA6/AliGenPythia.h
Trigger condition with eta range
[u/mrichter/AliRoot.git] / PYTHIA6 / AliGenPythia.h
1 #ifndef ALIGENPYTHIA_H
2 #define ALIGENPYTHIA_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice                               */
5
6
7 /* $Id$ */
8
9 //
10 // Generator using the TPythia interface (via AliPythia)
11 // to generate pp collisions.
12 // Using SetNuclei() also nuclear modifications to the structure functions
13 // can be taken into account. This makes, of course, only sense for the
14 // generation of the products of hard processes (heavy flavor, jets ...)
15 //
16 // andreas.morsch@cern.ch
17 //
18
19 #include "AliGenMC.h"
20 #include "AliPythia.h"
21
22 class AliPythia;
23 class TParticle;
24 class AliGenPythiaEventHeader;
25 class AliGenEventHeader;
26 class AliStack;
27 class AliRunLoader;
28 class TObjArray; 
29
30 class AliGenPythia : public AliGenMC
31 {
32  public:
33
34     typedef enum {kFlavorSelection, kParentSelection, kHeavyFlavor} StackFillOpt_t;
35     typedef enum {kCountAll, kCountParents, kCountTrackables} CountMode_t;
36     typedef enum {kCluster, kCell} JetRecMode_t;
37           
38     AliGenPythia();
39     AliGenPythia(Int_t npart);
40     virtual ~AliGenPythia();
41     virtual void    Generate();
42     virtual void    Init();
43     // Range of events to be printed
44     virtual void    SetEventListRange(Int_t eventFirst=-1, Int_t eventLast=-1);
45     // Select process type
46     virtual void    SetProcess(Process_t proc = kPyCharm) {fProcess = proc;}
47     virtual void    SetTune(Int_t itune) {fItune = itune;}
48
49     // Select structure function
50     virtual void    SetStrucFunc(StrucFunc_t func =  kCTEQ5L) {fStrucFunc = func;}
51     // Select pt of hard scattering 
52     virtual void    SetPtHard(Float_t ptmin = 0, Float_t ptmax = 1.e10)
53         {fPtHardMin = ptmin; fPtHardMax = ptmax; }
54     // y of hard scattering
55     virtual void    SetYHard(Float_t ymin = -1.e10, Float_t ymax = 1.e10)
56         {fYHardMin = ymin; fYHardMax = ymax; }
57     // Set initial and final state gluon radiation
58     virtual void    SetGluonRadiation(Int_t iIn, Int_t iFin)
59         {fGinit = iIn; fGfinal = iFin;}
60     virtual void SetColorReconnectionOff(Int_t iflag=0){fCRoff=iflag;} 
61    // Intrinsic kT
62     virtual void    SetPtKick(Float_t kt = 1.)
63         {fPtKick = kt;}
64     // Use the Pythia 6.3 new multiple interations scenario
65     virtual void    UseNewMultipleInteractionsScenario() {fNewMIS = kTRUE;}
66     // Switch off heavy flavors
67     virtual void    SwitchHFOff() {fHFoff = kTRUE;}
68     // Set centre of mass energy
69     virtual void    SetEnergyCMS(Float_t energy = 5500) {fEnergyCMS = energy;}
70     // Treat protons as inside nuclei with mass numbers a1 and a2
71     virtual void    SetNuclei(Int_t a1, Int_t a2, Int_t pdfset = 0);
72     // Set colliding nuclei ("p","n",...)
73     virtual void    SetCollisionSystem(TString projectile, TString target) { fProjectile = projectile; fTarget = target; }
74     virtual void    SetNuclearPDF(Int_t pdf) {fNucPdf = pdf;}
75     virtual void    SetUseNuclearPDF(Bool_t val) {fUseNuclearPDF = val;}
76     virtual void    SetUseLorentzBoost(Bool_t val) {fUseLorentzBoost = val;}
77     //
78     // Trigger options
79     //
80     // Energy range for jet trigger
81     virtual void    SetJetEtRange(Float_t etmin = 0., Float_t etmax = 1.e4)
82         {fEtMinJet = etmin; fEtMaxJet = etmax;}
83     // Eta range for jet trigger
84     virtual void    SetJetEtaRange(Float_t etamin = -20., Float_t etamax = 20.)
85         {fEtaMinJet = etamin; fEtaMaxJet = etamax;}
86     // Phi range for jet trigger
87     virtual void    SetJetPhiRange(Float_t phimin = 0., Float_t phimax = 360.)
88         {fPhiMinJet = TMath::Pi()*phimin/180.; fPhiMaxJet = TMath::Pi()*phimax/180.;}
89     // Jet reconstruction mode; default is cone algorithm
90     virtual void    SetJetReconstructionMode(Int_t mode = kCell) {fJetReconstruction = mode;}
91     // Eta range for gamma trigger 
92     virtual void    SetGammaEtaRange(Float_t etamin = -20., Float_t etamax = 20.)
93         {fEtaMinGamma = etamin; fEtaMaxGamma = etamax;}
94     // Phi range for gamma trigger
95     virtual void    SetGammaPhiRange(Float_t phimin = 0., Float_t phimax = 360.)
96         {fPhiMinGamma = TMath::Pi()*phimin/180.; fPhiMaxGamma = TMath::Pi()*phimax/180.;}
97   
98     // Select events with fragmentation photon, decay photon, pi0 or eta going to PHOS or EMCAL and central barrel
99     virtual Bool_t TriggerOnSelectedParticles(Int_t np);
100   
101     virtual void  SetCheckPHOS         (Bool_t b) {fCheckPHOS    = b;}
102     virtual void  SetCheckEMCAL        (Bool_t b) {fCheckEMCAL   = b;}
103     virtual void  SetCheckBarrel       (Bool_t b) {fCheckBarrel  = b;}
104   
105     //virtual void  SetElectronInEMCAL   (Bool_t b) {fEleInEMCAL   = b;}
106     //virtual void  SetPhotonInPHOS      (Bool_t b) {fCheckPHOS    = b; fPhotonInCalo     = b;} // Not in use
107
108     virtual void  SetFragPhotonInCalo  (Bool_t b) {                   fFragPhotonInCalo = b;}
109     virtual void  SetFragPhotonInBarrel(Bool_t b) {fCheckBarrel  = b; fFragPhotonInCalo = b;}
110     virtual void  SetFragPhotonInEMCAL (Bool_t b) {fCheckEMCAL   = b; fFragPhotonInCalo = b;}
111     virtual void  SetFragPhotonInPHOS  (Bool_t b) {fCheckPHOS    = b; fFragPhotonInCalo = b;}
112   
113     virtual void  SetHadronInCalo      (Bool_t b) {                   fHadronInCalo     = b;}
114     virtual void  SetHadronInBarrel    (Bool_t b) {fCheckBarrel  = b; fHadronInCalo     = b;}
115     virtual void  SetHadronInEMCAL     (Bool_t b) {fCheckEMCAL   = b; fHadronInCalo     = b;}
116     virtual void  SetHadronInPHOS      (Bool_t b) {fCheckPHOS    = b; fHadronInCalo     = b;}
117   
118     virtual void  SetElectronInCalo    (Bool_t b) {                   fEleInCalo        = b;}
119     virtual void  SetElectronInBarrel  (Bool_t b) {fCheckBarrel  = b; fEleInCalo        = b;}
120     virtual void  SetElectronInEMCAL   (Bool_t b) {fCheckEMCAL   = b; fEleInCalo        = b;}
121     virtual void  SetElectronInPHOS    (Bool_t b) {fCheckPHOS    = b; fEleInCalo        = b;}
122   
123     virtual void  SetDecayPhotonInCalo (Bool_t d)  {fDecayPhotonInCalo = d;}
124     virtual void  SetDecayPhotonInBarrel(Bool_t d) {fDecayPhotonInCalo = d; fCheckBarrel  = d;}
125     virtual void  SetDecayPhotonInEMCAL(Bool_t d)  {fDecayPhotonInCalo = d; fCheckEMCAL   = d;}
126     virtual void  SetDecayPhotonInPHOS (Bool_t d)  {fDecayPhotonInCalo = d; fCheckPHOS    = d;}
127   
128     virtual void  SetPi0InCalo         (Bool_t b, Bool_t f = kFALSE) {fPi0InCalo = b; fForceNeutralMeson2PhotonDecay = f;}
129     virtual void  SetPi0InBarrel       (Bool_t b, Bool_t f = kFALSE) {fPi0InCalo = b; fForceNeutralMeson2PhotonDecay = f; fCheckBarrel= b; }
130     virtual void  SetPi0InEMCAL        (Bool_t b, Bool_t f = kFALSE) {fPi0InCalo = b; fForceNeutralMeson2PhotonDecay = f; fCheckEMCAL = b; }
131     virtual void  SetPi0InPHOS         (Bool_t b, Bool_t f = kFALSE) {fPi0InCalo = b; fForceNeutralMeson2PhotonDecay = f; fCheckPHOS  = b; }
132  
133     virtual void  SetEtaInCalo         (Bool_t b, Bool_t f = kFALSE) {fEtaInCalo = b; fForceNeutralMeson2PhotonDecay = f;}
134     virtual void  SetEtaInBarrel       (Bool_t b, Bool_t f = kFALSE) {fEtaInCalo = b; fForceNeutralMeson2PhotonDecay = f; fCheckBarrel= b; }
135     virtual void  SetEtaInEMCAL        (Bool_t b, Bool_t f = kFALSE) {fEtaInCalo = b; fForceNeutralMeson2PhotonDecay = f; fCheckEMCAL = b; }
136     virtual void  SetEtaInPHOS         (Bool_t b, Bool_t f = kFALSE) {fEtaInCalo = b; fForceNeutralMeson2PhotonDecay = f; fCheckPHOS  = b; }
137
138     virtual void  SetPi0PhotonDecayInBarrel(Bool_t b, Bool_t f = kFALSE) {fPi0InCalo = b; fDecayPhotonInCalo = b; fForceNeutralMeson2PhotonDecay = f; fCheckBarrel  = b; }
139     virtual void  SetPi0PhotonDecayInEMCAL (Bool_t b, Bool_t f = kFALSE) {fPi0InCalo = b; fDecayPhotonInCalo = b; fForceNeutralMeson2PhotonDecay = f; fCheckEMCAL   = b; }
140     virtual void  SetPi0PhotonDecayInPHOS  (Bool_t b, Bool_t f = kFALSE) {fPi0InCalo = b; fDecayPhotonInCalo = b; fForceNeutralMeson2PhotonDecay = f; fCheckPHOS    = b; }
141
142     virtual void  SetEtaPhotonDecayInBarrel(Bool_t b, Bool_t f = kFALSE) {fEtaInCalo = b; fDecayPhotonInCalo = b; fForceNeutralMeson2PhotonDecay = f; fCheckBarrel  = b; }
143     virtual void  SetEtaPhotonDecayInEMCAL (Bool_t b, Bool_t f = kFALSE) {fEtaInCalo = b; fDecayPhotonInCalo = b; fForceNeutralMeson2PhotonDecay = f; fCheckEMCAL   = b; }
144     virtual void  SetEtaPhotonDecayInPHOS  (Bool_t b, Bool_t f = kFALSE) {fEtaInCalo = b; fDecayPhotonInCalo = b; fForceNeutralMeson2PhotonDecay = f; fCheckPHOS    = b; }
145
146   
147     // Trigger on a minimum multiplicity
148     virtual void  SetTriggerChargedMultiplicity(Int_t multiplicity, Float_t etamax = 0, Float_t ptmin = -1.) 
149     {fTriggerMultiplicity = multiplicity; fTriggerMultiplicityEta = etamax; 
150       fTriggerMultiplicityPtMin = ptmin;}
151     
152     // Trigger on a minimum multiplicity for a given eta range
153     virtual void  SetTriggerMultiplicityEtaRange(Int_t multiplicity, Float_t etamin = 0., Float_t etamax = 0., Float_t ptmin = -1.) 
154     {fTriggerMultiplicity = multiplicity; fTriggerMultiplicityEtaMin = etamin; fTriggerMultiplicityEtaMax = etamax; 
155       fTriggerMultiplicityPtMin = ptmin;}
156     
157     // Calorimeters acceptance
158     // Set Phi in degrees, and Eta coverage, should not be negative
159     virtual void  SetBarrelAcceptance(Float_t deta) {fTriggerEta = deta ;}
160     virtual void  SetEMCALAcceptance (Float_t phimin, Float_t phimax, Float_t deta) {fEMCALMinPhi = phimin ; fEMCALMaxPhi = phimax ; fEMCALEta = deta ; }
161     virtual void  SetPHOSAcceptance  (Float_t phimin, Float_t phimax, Float_t deta) {fPHOSMinPhi  = phimin ; fPHOSMaxPhi  = phimax ; fPHOSEta  = deta ; }
162     virtual void  SetRotateParticleInPHOSeta(Bool_t b) {fCheckPHOSeta = b;}
163   
164     virtual void  SetTriggerParticleMinPt(Float_t pt) {fTriggerParticleMinPt = pt;}
165 //    virtual void  SetPhotonMinPt(Float_t pt)          {fPhotonMinPt = pt;}
166 //    virtual void  SetElectronMinPt(Float_t pt)        {fElectronMinPt = pt;}
167     // Trigger and rotate event 
168     void RotatePhi(Bool_t& okdd);
169   
170     // Trigger on a single particle (not related to calorimeter trigger above)
171     virtual void    SetTriggerParticle(Int_t particle = 0, Float_t etamax = 0.9, Float_t ptmin = -1, Float_t ptmax = 1000) 
172         {fTriggerParticle = particle; fTriggerEta = etamax; fTriggerMinPt = ptmin; fTriggerMaxPt = ptmax;}
173
174     //
175     // Heavy flavor options
176     //
177     // Set option for feed down from higher family
178     virtual void SetFeedDownHigherFamily(Bool_t opt) {
179         fFeedDownOpt = opt;
180     }
181     // Set option for selecting particles kept in stack according to flavor
182     // or to parent selection
183     virtual void SetStackFillOpt(StackFillOpt_t opt) {
184         fStackFillOpt = opt;
185     }
186     // Set fragmentation option
187     virtual void SetFragmentation(Bool_t opt) {
188         fFragmentation = opt;
189     }
190     // Set counting mode
191     virtual void SetCountMode(CountMode_t mode) {
192         fCountMode = mode;
193     }
194     //
195     // Quenching
196     //
197     // Set quenching mode 0 = no, 1 = AM, 2 = IL,  3 = NA, 4 = ACS
198     virtual void SetQuench(Int_t flag = 0) {fQuench = flag;}
199     // Set transport coefficient.
200     void SetQhat(Float_t qhat) {fQhat = qhat;}
201     //Set initial medium length.
202     void SetLength(Float_t length) {fLength = length;}
203     //set parameters for pyquen afterburner
204     virtual void SetPyquenPar(Float_t t0=1., Float_t tau0=0.1, Int_t nf=0,Int_t iengl=0, Int_t iangl=3)
205     {fpyquenT = t0; fpyquenTau = tau0; fpyquenNf=nf;fpyquenEloss=iengl;fpyquenAngle=iangl;}
206     virtual void SetHadronisation(Int_t flag = 1) {fHadronisation = flag;}
207     virtual void SetPatchOmegaDalitz(Int_t flag = 1) {fPatchOmegaDalitz = flag;}
208     virtual void SetDecayerExodus(Int_t flag = 1) {fDecayerExodus = flag;}
209     virtual void SetReadFromFile(const Text_t *filname) {fkFileName = filname;  fReadFromFile = 1;}    
210     virtual void SetReadLHEF(const Text_t *filename) {fkNameLHEF = filename; fReadLHEF = 1;}    
211
212     //
213     // Pile-up
214     //
215     // Get interaction rate for pileup studies
216     virtual void    SetInteractionRate(Float_t rate,Float_t timewindow = 90.e-6);
217     virtual Float_t GetInteractionRate() const {return fInteractionRate;}
218     // get cross section of process
219     virtual Float_t GetXsection() const {return fXsection;}
220     // get triggered jets
221     void GetJets(Int_t& njets, Int_t& ntrig, Float_t jets[4][10]);
222     void RecJetsUA1(Int_t& njets, Float_t jets[4][50]);
223     void SetPycellParameters(Float_t etamax = 2., Int_t neta = 274, Int_t nphi = 432,
224                              Float_t thresh = 0., Float_t etseed = 4.,
225                              Float_t minet = 10., Float_t r = 1.);
226     
227   void LoadEvent(AliStack* stack, Int_t flag = 0, Int_t reHadr = 0);
228   void LoadEvent(const TObjArray* stack, Int_t flag = 0, Int_t reHadr = 0);
229     // Getters
230     virtual Process_t    GetProcess() const {return fProcess;}
231     virtual StrucFunc_t  GetStrucFunc() const {return fStrucFunc;}
232     virtual void         GetPtHard(Float_t& ptmin, Float_t& ptmax) const
233         {ptmin = fPtHardMin; ptmax = fPtHardMax;}
234     virtual void         GetNuclei(Int_t&  a1, Int_t& a2) const
235         {a1 = fAProjectile; a2 = fATarget;}
236     virtual void         GetJetEtRange(Float_t& etamin, Float_t& etamax) const
237         {etamin = fEtaMinJet; etamax = fEtaMaxJet;}
238     virtual void         GetJetPhiRange(Float_t& phimin, Float_t& phimax) const
239         {phimin = fPhiMinJet*180./TMath::Pi(); phimax = fPhiMaxJet*180/TMath::Pi();}
240     virtual void         GetGammaEtaRange(Float_t& etamin, Float_t& etamax) const
241         {etamin = fEtaMinGamma; etamax = fEtaMaxGamma;}
242     virtual void         GetGammaPhiRange(Float_t& phimin, Float_t& phimax) const
243         {phimin = fPhiMinGamma*180./TMath::Pi(); phimax = fPhiMaxGamma*180./TMath::Pi();}
244     //
245     Bool_t CheckDetectorAcceptance(Float_t phi, Float_t eta, Int_t iparticle);
246     Bool_t IsInEMCAL (Float_t phi, Float_t eta) const;
247     Bool_t IsInPHOS  (Float_t phi, Float_t eta, Int_t iparticle) ;
248     Bool_t IsInBarrel(Float_t eta) const;
249     Bool_t IsFromHeavyFlavor(Int_t ipart);
250     //
251     virtual void FinishRun();
252     Bool_t CheckTrigger(const TParticle* jet1, const TParticle* jet2);
253     //Used in some processes to selected child properties
254     Bool_t CheckKinematicsOnChild();
255     void     GetSubEventTime();
256
257     void SetTuneForDiff(Bool_t a=kTRUE) {fkTuneForDiff=a;}
258     AliDecayer * GetDecayer(){return fDecayer;}
259
260  protected:
261     // adjust the weight from kinematic cuts
262     void     AdjustWeights() const;
263     Int_t    GenerateMB();
264     void     MakeHeader();    
265     void     GeneratePileup();
266     Process_t   fProcess;           //Process type
267     Int_t       fItune;             // Pythia tune > 6.4
268     StrucFunc_t fStrucFunc;         //Structure Function
269     Float_t     fKineBias;          //!Bias from kinematic selection
270     Int_t       fTrials;            //!Number of trials for current event
271     Int_t       fTrialsRun;         //!Number of trials for run
272     Float_t     fQ;                 //Mean Q
273     Float_t     fX1;                //Mean x1
274     Float_t     fX2;                //Mean x2
275     Float_t     fEventTime;         //Time of the subevent
276     Float_t     fInteractionRate;   //Interaction rate (set by user)
277     Float_t     fTimeWindow;        //Time window for pileup events (set by user)
278     Int_t       fCurSubEvent;       //Index of the current sub-event
279     TArrayF     *fEventsTime;       //Subevents time for pileup
280     Int_t       fNev;               //Number of events 
281     Int_t       fFlavorSelect;      //Heavy Flavor Selection
282     Float_t     fXsection;          //Cross-section
283     AliPythia   *fPythia;           //!Pythia 
284     Float_t     fPtHardMin;         //lower pT-hard cut 
285     Float_t     fPtHardMax;         //higher pT-hard cut
286     Float_t     fYHardMin;          //lower  y-hard cut 
287     Float_t     fYHardMax;          //higher y-hard cut
288     Int_t       fGinit;             //initial state gluon radiation
289     Int_t       fGfinal;            //final state gluon radiation
290     Int_t       fCRoff;             //color reconnection off in the pythia6 annealying model
291     Int_t       fHadronisation;     //hadronisation
292     Bool_t      fPatchOmegaDalitz;  //flag for omega dalitz decay patch
293     Bool_t      fDecayerExodus;     //flag for exodus decayer
294     Int_t       fNpartons;          //Number of partons before hadronisation
295     Int_t       fReadFromFile;      //read partons from file
296     Int_t       fReadLHEF;          //read lhef file
297     Int_t       fQuench;            //Flag for quenching
298     Float_t     fQhat;              //Transport coefficient (GeV^2/fm)
299     Float_t     fLength;            //Medium length (fm)
300     Float_t     fpyquenT;           //Pyquen initial temperature 
301     Float_t     fpyquenTau;         //Pyquen initial proper time
302     Int_t       fpyquenNf;          //Pyquen number of flavours into the game
303     Int_t       fpyquenEloss;       //Pyquen type of energy loss 
304     Int_t       fpyquenAngle;       //Pyquen radiation angle for gluons
305     Float_t     fImpact;            //Impact parameter for quenching simulation (q-pythia) 
306     Float_t     fPtKick;            //Transverse momentum kick
307     Bool_t      fFullEvent;         //!Write Full event if true
308     AliDecayer  *fDecayer;          //!Pointer to the decayer instance
309     Int_t       fDebugEventFirst;   //!First event to debug
310     Int_t       fDebugEventLast;    //!Last  event to debug
311     Float_t     fEtMinJet;          //Minimum et of triggered Jet
312     Float_t     fEtMaxJet;          //Maximum et of triggered Jet
313     Float_t     fEtaMinJet;         //Minimum eta of triggered Jet
314     Float_t     fEtaMaxJet;         //Maximum eta of triggered Jet
315     Float_t     fPhiMinJet;         //Minimum phi of triggered Jet
316     Float_t     fPhiMaxJet;         //Maximum phi of triggered Jet
317     Int_t       fJetReconstruction; //Jet Reconstruction mode 
318     Float_t     fEtaMinGamma;       // Minimum eta of triggered gamma
319     Float_t     fEtaMaxGamma;       // Maximum eta of triggered gamma
320     Float_t     fPhiMinGamma;       // Minimum phi of triggered gamma
321     Float_t     fPhiMaxGamma;       // Maximum phi of triggered gamma
322     Float_t     fPycellEtaMax;      // Max. eta for Pycell 
323     Int_t       fPycellNEta;        // Number of eta bins for Pycell 
324     Int_t       fPycellNPhi;        // Number of phi bins for Pycell
325     Float_t     fPycellThreshold;   // Pycell threshold
326     Float_t     fPycellEtSeed;      // Pycell seed
327     Float_t     fPycellMinEtJet;    // Pycell min. jet et
328     Float_t     fPycellMaxRadius;   // Pycell cone radius
329     StackFillOpt_t fStackFillOpt;   // Stack filling with all particles with
330                                     // that flavour or only with selected
331                                     // parents and their decays
332     Bool_t fFeedDownOpt;            // Option to set feed down from higher
333                                     // quark families (e.g. b->c)
334     Bool_t  fFragmentation;         // Option to activate fragmentation by Pythia
335     Bool_t  fSetNuclei;             // Flag indicating that SetNuclei has been called
336     Bool_t  fUseNuclearPDF;         // flag if nuclear pdf should be applied
337     Bool_t  fUseLorentzBoost;       // flag if lorentz boost should be applied
338     Bool_t  fNewMIS;                // Flag for the new multipple interactions scenario
339     Bool_t  fHFoff;                 // Flag for switching heafy flavor production off
340     Int_t   fNucPdf;                // Nuclear pdf 0: EKS98 1: EPS08
341     Int_t   fTriggerParticle;       // Trigger on this particle ...
342     Float_t fTriggerEta;            // .. within |eta| < fTriggerEta
343     Float_t fTriggerMinPt;          // .. within pt > fTriggerMinPt
344     Float_t fTriggerMaxPt;          // .. within pt < fTriggerMaxPt
345     Int_t       fTriggerMultiplicity;       // Trigger on events with a minimum charged multiplicity
346     Float_t     fTriggerMultiplicityEta;    // in a given eta range
347     Float_t     fTriggerMultiplicityEtaMin;    // in a given eta min
348     Float_t     fTriggerMultiplicityEtaMax;    // in a given eta max
349     Float_t     fTriggerMultiplicityPtMin;  // above this pT 
350     CountMode_t fCountMode;         // Options for counting when the event will be finished.     
351     // fCountMode = kCountAll         --> All particles that end up in the
352     //                                    stack are counted
353     // fCountMode = kCountParents     --> Only selected parents are counted
354     // fCountMode = kCountTrackabless --> Only particles flagged for tracking
355     //                                     are counted
356     //
357     //
358
359     AliGenPythiaEventHeader* fHeader;  //! Event header
360     AliRunLoader*            fRL;      //! Run Loader
361     const Text_t* fkFileName;          //! Name of file to read from
362     const Text_t* fkNameLHEF;          //! Name of lhef file to read from
363     Bool_t fFragPhotonInCalo; // Option to ask for Fragmentation Photon in calorimeters acceptance
364     Bool_t fHadronInCalo;     // Option to ask for hadron (not pi0) in calorimeters acceptance
365     Bool_t fPi0InCalo;        // Option to ask for Pi0 in calorimeters acceptance
366     Bool_t fEtaInCalo;        // Option to ask for Eta in calorimeters acceptance
367     Bool_t fPhotonInCalo;     // Option to ask for Photon in calorimeter acceptance (not in use)
368     Bool_t fDecayPhotonInCalo;// Option to ask for Decay Photon in calorimeter acceptance
369     Bool_t fForceNeutralMeson2PhotonDecay; // Option to ask for Pi0/Eta in calorimeters acceptance when decay into 2 photons
370     Bool_t fEleInCalo;        // Option to ask for Electron in EMCAL acceptance
371     Bool_t fEleInEMCAL;       // Option to ask for Electron in EMCAL acceptance (not in use)
372     Bool_t fCheckBarrel;      // Option to ask for FragPhoton or Pi0 or Eta or gamma decays in central barrel acceptance
373     Bool_t fCheckEMCAL;       // Option to ask for FragPhoton or Pi0 or Eta or gamma decays in calorimeters EMCAL acceptance
374     Bool_t fCheckPHOS;        // Option to ask for FragPhoton or Pi0 or Eta or gamma decays in calorimeters PHOS acceptance
375     Bool_t fCheckPHOSeta;     // Option to ask for rotate event particles in phi to have in PHOS acceptance a requested particle that previously had the good eta
376     Int_t  fPHOSRotateCandidate; // Internal member to select the particle candidate to trigger the event phi rotation, to put it in PHOS phi acceptance
377     Float_t fTriggerParticleMinPt; // Minimum momentum of Fragmentation Photon or Pi0 or other hadron
378     Float_t fPhotonMinPt;          // Minimum momentum of Photon  (not in use)
379     Float_t fElectronMinPt;        // Minimum momentum of Electron (not in use) 
380     //Calorimeters eta-phi acceptance 
381     Float_t fPHOSMinPhi;           // Minimum phi PHOS, degrees
382     Float_t fPHOSMaxPhi;           // Maximum phi PHOS, degrees
383     Float_t fPHOSEta;              // Minimum eta PHOS, coverage delta eta
384     Float_t fEMCALMinPhi;          // Minimum phi EMCAL, degrees
385     Float_t fEMCALMaxPhi;          // Maximum phi EMCAL, degrees
386     Float_t fEMCALEta;             // Maximum eta EMCAL, coverage delta eta
387
388     Bool_t fkTuneForDiff;    // Pythia tune 
389     Int_t  fProcDiff;
390  private:
391     AliGenPythia(const AliGenPythia &Pythia);
392     AliGenPythia & operator=(const AliGenPythia & rhs);
393
394
395     Bool_t CheckDiffraction();
396     Bool_t GetWeightsDiffraction(Double_t M, Double_t &Mmin, Double_t &Mmax, 
397                                                Double_t &wSD, Double_t &wDD, Double_t &wND);
398
399     ClassDef(AliGenPythia, 14) // AliGenerator interface to Pythia
400 };
401 #endif
402
403
404
405
406