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