]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGDQ/dielectron/AliReducedEvent.h
including switch to set on/off iso-track core removal, cleaning and bug fix
[u/mrichter/AliRoot.git] / PWGDQ / dielectron / AliReducedEvent.h
1 // Classes used for creating a reduced information tree
2 // Author: Ionut-Cristian Arsene (i.c.arsene@gsi.de)
3 // 
4 //  Basic structure:
5 //  1. Event wise information
6 //  2. List of tracks in the event
7 //  3. List of resonance candidates
8
9 #ifndef ALIREDUCEDEVENT_H
10 #define ALIREDUCEDEVENT_H
11
12 #include <TClonesArray.h>
13 #include <TBits.h>
14 #include <TMath.h>
15
16
17 const Int_t fgkNMaxHarmonics = 10;
18 const Float_t gZdcNalpha= 1.0;
19
20 //_____________________________________________________________________
21 class AliReducedTrack : public TObject {
22
23   friend class AliAnalysisTaskReducedTree;  // friend analysis task which fills the object
24
25  public:
26   AliReducedTrack();
27   ~AliReducedTrack();
28
29   // getters
30   UShort_t TrackId()                     const {return fTrackId;}
31   ULong_t  Status()                      const {return fStatus;}
32   Bool_t   CheckTrackStatus(UInt_t flag) const {return (flag<8*sizeof(ULong_t) ? (fStatus&(1<<flag)) : kFALSE);}
33   Int_t    Charge()                      const {return (fGlobalPt>0.0 ? +1 : -1);}
34   Float_t  Px()                          const {return TMath::Abs(fGlobalPt)*TMath::Cos(fGlobalPhi);}
35   Float_t  Py()                          const {return TMath::Abs(fGlobalPt)*TMath::Sin(fGlobalPhi);}
36   Float_t  Pz()                          const {return TMath::Abs(fGlobalPt)*TMath::SinH(fGlobalEta);}
37   Float_t  P()                           const {return TMath::Abs(fGlobalPt)*TMath::CosH(fGlobalEta);};
38   Float_t  Phi()                         const {return fGlobalPhi;}
39   Float_t  Pt()                          const {return TMath::Abs(fGlobalPt);}
40   Float_t  Eta()                         const {return fGlobalEta;}
41   Float_t  Theta()                       const {return TMath::ACos(TMath::TanH(fGlobalEta));}
42   Float_t  PxTPC()                       const {return fTPCPt*TMath::Cos(fTPCPhi);}
43   Float_t  PyTPC()                       const {return fTPCPt*TMath::Sin(fTPCPhi);}
44   Float_t  PzTPC()                       const {return fTPCPt*TMath::SinH(fTPCEta);}
45   Float_t  PTPC()                        const {return fTPCPt*TMath::CosH(fTPCEta);};
46   Float_t  PhiTPC()                      const {return fTPCPhi;}
47   Float_t  PtTPC()                       const {return fTPCPt;}
48   Float_t  EtaTPC()                      const {return fTPCEta;}
49   Float_t  ThetaTPC()                    const {return TMath::ACos(TMath::TanH(fTPCEta));}
50   Float_t  Pin()                         const {return fMomentumInner;}
51   Float_t  DCAxy()                       const {return fDCA[0];}
52   Float_t  DCAz()                        const {return fDCA[1];}
53   Float_t  TrackLength()                 const {return fTrackLength;}
54   
55   UShort_t ITSncls()                const;
56   UChar_t  ITSclusterMap()          const {return fITSclusterMap;}
57   Bool_t   ITSLayerHit(Int_t layer) const {return (layer>=0 && layer<6 ? (fITSclusterMap&(1<<layer)) : kFALSE);};
58   Float_t  ITSsignal()              const {return fITSsignal;}
59   Float_t  ITSnSig(Int_t specie)    const {return (specie>=0 && specie<=3 ? fITSnSig[specie] : -999.);}
60   Float_t  ITSchi2()                const {return fITSchi2;}
61   
62   UChar_t TPCncls()                        const {return fTPCNcls;}
63   UChar_t TPCFindableNcls()                const {return fTPCNclsF;}
64   UChar_t TPCCrossedRows()                 const {return fTPCCrossedRows;}
65   UChar_t TPCnclsIter1()                   const {return fTPCNclsIter1;}
66   UChar_t TPCClusterMap()                  const {return fTPCClusterMap;}
67   Int_t   TPCClusterMapBitsFired()         const;
68   Bool_t  TPCClusterMapBitFired(Int_t bit) const {return (bit>=0 && bit<8 ? (fTPCClusterMap&(1<<bit)) : kFALSE);};
69   Float_t TPCsignal()                      const {return fTPCsignal;}
70   UChar_t TPCsignalN()                     const {return fTPCsignalN;}
71   Float_t TPCnSig(Int_t specie)            const {return (specie>=0 && specie<=3 ? fTPCnSig[specie] : -999.);}
72   Float_t TPCchi2()                        const {return fTPCchi2;}
73   
74   Float_t  TOFbeta()             const {return fTOFbeta;}    
75   Float_t  TOFnSig(Int_t specie) const {return (specie>=0 && specie<=3 ? fTOFnSig[specie] : -999.);}
76   Short_t  TOFdeltaBC()          const {return fTOFdeltaBC;}
77   
78   Int_t    TRDntracklets(Int_t type)  const {return (type==0 || type==1 ? fTRDntracklets[type] : -1);}
79   Float_t  TRDpid(Int_t specie)       const {return (specie>=0 && specie<=1 ? fTRDpid[specie] : -999.);}
80   Float_t  TRDpidLQ1D(Int_t specie)   const {return (specie>=0 && specie<=1 ? fTRDpid[specie] : -999.);}
81   Float_t  TRDpidLQ2D(Int_t specie)   const {return (specie>=0 && specie<=1 ? fTRDpidLQ2D[specie] : -999.);}
82   
83   Int_t    CaloClusterId() const {return fCaloClusterId;}
84     
85   Float_t  BayesPID(Int_t specie) const {return (specie>=0 && specie<=2 ? fBayesPID[specie] : -999.);}
86   
87   Bool_t UsedForQvector()             const {return fFlags&(UShort_t(1)<<0);}
88   Bool_t TestFlag(UShort_t iflag)     const {return (iflag<8*sizeof(UShort_t) ? fFlags&(UShort_t(1)<<iflag) : kFALSE);}
89   Bool_t SetFlag(UShort_t iflag)            {if (iflag>=8*sizeof(UShort_t)) return kFALSE; fFlags|=(UShort_t(1)<<iflag); return kTRUE;}
90   Bool_t IsGammaLeg()                 const {return fFlags&(UShort_t(1)<<1);}
91   Bool_t IsPureGammaLeg()             const {return fFlags&(UShort_t(1)<<8);}
92   Bool_t IsK0sLeg()                   const {return fFlags&(UShort_t(1)<<2);}
93   Bool_t IsPureK0sLeg()               const {return fFlags&(UShort_t(1)<<9);}
94   Bool_t IsLambdaLeg()                const {return fFlags&(UShort_t(1)<<3);}
95   Bool_t IsPureLambdaLeg()            const {return fFlags&(UShort_t(1)<<10);}
96   Bool_t IsALambdaLeg()               const {return fFlags&(UShort_t(1)<<4);}
97   Bool_t IsPureALambdaLeg()           const {return fFlags&(UShort_t(1)<<11);}
98   Bool_t IsKink(Int_t i=0)            const {return (i>=0 && i<3 ? fFlags&(UShort_t(1)<<(5+i)) : kFALSE);}
99   Bool_t TestFlagMore(UShort_t iflag) const {return (iflag<8*sizeof(ULong_t) ? fMoreFlags&(ULong_t(1)<<iflag) : kFALSE);}
100   Bool_t SetFlagMore(UShort_t iflag)        {if(iflag>=8*sizeof(ULong_t)) return kFALSE; fMoreFlags|=(ULong_t(1)<<iflag); return kTRUE;}
101   Bool_t UnsetFlagMore(UShort_t iflag)      {if(iflag>=8*sizeof(ULong_t)) return kFALSE; fMoreFlags^=(ULong_t(1)<<iflag); return kTRUE;}
102   ULong_t GetFlagsMore()               const {return fMoreFlags;}
103   
104  private:
105   UShort_t fTrackId;            // track id 
106   ULong_t fStatus;              // tracking status
107   Float_t fGlobalPhi;           // phi at the vertex from global track, in the [0,2pi) interval
108   Float_t fGlobalPt;            // pt*charge at the vertex from global track
109   Float_t fGlobalEta;           // eta at the vertex from global track
110   Float_t fTPCPhi;              // phi at the vertex from TPC alone tracking , in the [0,2pi) interval
111   Float_t fTPCPt;               // pt at the vertex from TPC alone tracking  
112   Float_t fTPCEta;              // eta at the vertex from TPC alone tracking 
113   Float_t fMomentumInner;       // inner param momentum (only the magnitude)
114   Float_t fDCA[2];              // DCA xy,z
115   Float_t fTrackLength;         // track length
116   
117   // ITS
118   UChar_t  fITSclusterMap;      // ITS cluster map
119   Float_t  fITSsignal;          // ITS signal
120   Float_t  fITSnSig[4];         // 0-electron; 1-pion; 2-kaon; 3-proton
121   Float_t  fITSchi2;            // ITS chi2 / cls
122   
123   // TPC
124   UChar_t fTPCNcls;            // TPC ncls                          
125   UChar_t fTPCCrossedRows;     // TPC crossed rows                  
126   UChar_t fTPCNclsF;           // TPC findable ncls                 
127   UChar_t fTPCNclsIter1;       // TPC no clusters after first iteration
128   UChar_t fTPCClusterMap;      // TPC cluster distribution map
129   Float_t fTPCsignal;          // TPC de/dx
130   UChar_t fTPCsignalN;         // TPC no clusters de/dx
131   Float_t fTPCnSig[4];         // 0-electron; 1-pion; 2-kaon; 3-proton
132   Float_t fTPCchi2;            // TPC chi2 / cls
133     
134   // TOF
135   Float_t fTOFbeta;             // TOF pid info
136   Float_t fTOFnSig[4];          // TOF n-sigma deviation from expected signal
137   Short_t fTOFdeltaBC;          // BC(event) - BC(track) estimated by TOF
138   
139   // TRD
140   UChar_t fTRDntracklets[2];       // 0 - AliESDtrack::GetTRDntracklets(); 1 - AliESDtrack::GetTRDntrackletsPID()   TODO: use only 1 char
141   Float_t fTRDpid[2];              // TRD pid 1D likelihoods, [0]-electron , [1]- pion
142   Float_t fTRDpidLQ2D[2];          // TRD pid 2D likelihoods, [0]-electron , [1]- pion
143   
144   // EMCAL/PHOS
145   Int_t  fCaloClusterId;          // ID for the calorimeter cluster (if any)
146   
147   // Bayesian PID
148   Float_t fBayesPID[3];           // Combined Bayesian PID   pi/K/p
149   
150   UShort_t fFlags;                // BIT0 toggled if track used for TPC event plane
151                                   // BIT1 toggled if track belongs to a gamma conversion
152                                   // BIT2 toggled if track belongs to a K0s
153                                   // BIT3 toggled if track belongs to a Lambda
154                                   // BIT4 toggled if track belongs to an Anti-Lambda
155                                   // BIT5 toggled if the track has kink0 index > 0
156                                   // BIT6 toggled if the track has kink1 index > 0
157                                   // BIT7 toggled if the track has kink2 index > 0
158                                   // BIT8 toggled if track belongs to a pure gamma conversion
159                                   // BIT9 toggled if track belongs to a pure K0s
160                                   // BIT10 toggled if track belongs to a pure Lambda
161                                   // BIT11 toggled if track belongs to a pure ALambda
162   ULong_t  fMoreFlags;            // Space reserved for more information which might be needed later for analysis
163     
164   AliReducedTrack(const AliReducedTrack &c);
165   AliReducedTrack& operator= (const AliReducedTrack &c);
166
167   ClassDef(AliReducedTrack, 4);
168 };
169
170
171 //_____________________________________________________________________
172 class AliReducedPair : public TObject {
173
174   friend class AliAnalysisTaskReducedTree;  // friend analysis task which fills the object
175
176  public:
177   enum CandidateType {
178     kK0sToPiPi=0,
179     kPhiToKK,
180     kLambda0ToPPi,
181     kALambda0ToPPi,
182     kJpsiToEE,
183     kUpsilon,
184     kGammaConv,
185     kNMaxCandidateTypes
186   };
187   AliReducedPair();
188   AliReducedPair(const AliReducedPair &c);
189   ~AliReducedPair();
190
191   // getters
192   Char_t   CandidateId()         const {return fCandidateId;}
193   Char_t   PairType()            const {return fPairType;}
194   Int_t    LegId(Int_t leg)      const {return (leg==0 || leg==1 ? fLegIds[leg] : -1);}
195   Float_t  Mass(Int_t idx=0)     const {return (idx>=0 && idx<4 ? fMass[idx] : -999.);}
196   Float_t  Px()                  const {return fPt*TMath::Cos(fPhi);}
197   Float_t  Py()                  const {return fPt*TMath::Sin(fPhi);}
198   Float_t  Pz()                  const {return fPt*TMath::SinH(fEta);}
199   Float_t  P()                   const {return fPt*TMath::CosH(fEta);}
200   Float_t  Phi()                 const {return fPhi;}
201   Float_t  Pt()                  const {return fPt;}
202   Float_t  Eta()                 const {return fEta;}
203   Float_t  Energy()              const;
204   Float_t  Rapidity()            const;
205   Float_t  Theta()               const {return TMath::ACos(TMath::TanH(fEta));}
206   Float_t  Lxy()                 const {return fLxy;}
207   Float_t  DecayRadius()         const {return fLxy;}
208   Float_t  LxyErr()              const {return fLxyErr;}
209   Float_t  PointingAngle()       const {return fPointingAngle;}
210   Float_t  Chi2()                const {return fChisquare;}
211   Bool_t   IsOnTheFly()          const {return fPairType;}
212   Bool_t   IsPureV0K0s()         const {return (fMCid&(UInt_t(1)<<1));}
213   Bool_t   IsPureV0Lambda()      const {return (fMCid&(UInt_t(1)<<2));}
214   Bool_t   IsPureV0ALambda()     const {return (fMCid&(UInt_t(1)<<3));}
215   Bool_t   IsPureV0Gamma()       const {return (fMCid&(UInt_t(1)<<4));}
216   UInt_t   MCid()                const {return fMCid;}
217   Bool_t   CheckMC(const Int_t flag) const {return (flag<32 ? (fMCid&(1<<flag)) : kFALSE);}
218   
219  private:
220   Char_t  fCandidateId;         // candidate type (K0s, Lambda, J/psi, phi, etc)
221   Char_t  fPairType;            // 0 ++; 1 +-; 2 -- for dielectron pairs; 0- offline, 1- on the fly for V0 candidates
222   UShort_t fLegIds[2];          // leg ids 
223   Float_t fMass[4];             // invariant mass for pairs (3 extra mass values for other V0 pid assumptions)
224                                 // idx=0 -> K0s assumption; idx=1 -> Lambda; idx=2 -> anti-Lambda; idx=3 -> gamma conversion
225   Float_t fPhi;                 // pair phi in the [0,2*pi) interval
226   Float_t fPt;                  // pair pt
227   Float_t fEta;                 // pair eta 
228   Float_t fLxy;                 // pseudo-proper decay length (pair candidates) or radius of the secondary vertex for V0s 
229   Float_t fLxyErr;              // error on Lxy
230   Float_t fPointingAngle;       // angle between the pair momentum vector and the secondary vertex position vector
231   Float_t fChisquare;           // chi2 for the legs matching
232   UInt_t  fMCid;                // Bit map with Monte Carlo info about the pair
233
234   AliReducedPair& operator= (const AliReducedPair &c);
235
236   ClassDef(AliReducedPair, 3);
237 };
238
239
240 //_____________________________________________________________________
241 class AliReducedFMD : public TObject {
242
243   friend class AliAnalysisTaskReducedTree;  // friend analysis task which fills the object
244
245  public:
246   AliReducedFMD();
247   ~AliReducedFMD();
248
249   // getters
250   UShort_t Multiplicity()           const {return fMultiplicity;}
251   //Float_t Eta()                                 const {return fEta;}
252   UShort_t Id()                 const {return fId;}
253
254   //void SetIgnoreSteamer();
255   
256  private:
257
258   UShort_t fMultiplicity;        
259   //Float_t fEta;            
260   UShort_t fId;           
261
262     
263   AliReducedFMD(const AliReducedFMD &c);
264   AliReducedFMD& operator= (const AliReducedFMD &c);
265
266   ClassDef(AliReducedFMD, 1);
267 };
268
269
270 //_________________________________________________________________________
271 class AliReducedEventFriend : public TObject {
272   
273   friend class AliAnalysisTaskReducedTree;    // friend analysis task which fills the object
274   
275  public: 
276   enum EventPlaneStatus {
277     kRaw=0,
278     kCalibrated,
279     kRecentered,
280     kShifted,
281     kNMaxFlowFlags
282   };
283   enum EventPlaneDetector {
284     kTPC=0,       
285     kTPCptWeights,
286     kTPCpos,
287     kTPCneg,
288     kVZEROA,
289     kVZEROC,
290     kFMD,
291     kZDCA,
292     kZDCC,
293     kNdetectors
294   };
295   
296   AliReducedEventFriend();
297   ~AliReducedEventFriend();
298   
299   Double_t Qx(Int_t det, Int_t harmonic)  const {return (det>=0 && det<kNdetectors && harmonic>0 && harmonic<=fgkNMaxHarmonics ? fQvector[det][harmonic-1][0] : -999.);}
300   Double_t Qy(Int_t det, Int_t harmonic)  const {return (det>=0 && det<kNdetectors && harmonic>0 && harmonic<=fgkNMaxHarmonics ? fQvector[det][harmonic-1][1] : -999.);}
301   Double_t EventPlane(Int_t det, Int_t h) const;
302   UChar_t GetEventPlaneStatus(Int_t det, Int_t h) const {return (det>=0 && det<kNdetectors && h>0 && h<=fgkNMaxHarmonics ? fEventPlaneStatus[det][h-1] : UChar_t(255));} 
303   Bool_t  CheckEventPlaneStatus(Int_t det, Int_t h, EventPlaneStatus flag) const;
304   void    CopyEvent(const AliReducedEventFriend* event);
305
306   void SetQx(Int_t det, Int_t harmonic, Float_t qx) { if(det>=0 && det<kNdetectors && harmonic>0 && harmonic<=fgkNMaxHarmonics) fQvector[det][harmonic-1][0]=qx;}
307   void SetQy(Int_t det, Int_t harmonic, Float_t qy) { if(det>=0 && det<kNdetectors && harmonic>0 && harmonic<=fgkNMaxHarmonics) fQvector[det][harmonic-1][1]=qy;}
308   void SetEventPlaneStatus(Int_t det, Int_t harmonic, EventPlaneStatus status) { 
309     if(det>=0 && det<kNdetectors && harmonic>0 && harmonic<=fgkNMaxHarmonics) 
310       fEventPlaneStatus[det][harmonic-1] |= (1<<status);
311   }
312   
313  private:
314   // Q-vectors for the first 10 harmonics from TPC, VZERO, FMD and ZDC detectors
315   Double_t fQvector[kNdetectors][fgkNMaxHarmonics][2];     // Q vector components for all detectors and 6 harmonics
316   UChar_t fEventPlaneStatus[kNdetectors][fgkNMaxHarmonics];  // Bit maps for the event plane status (1 char per detector and per harmonic)
317    
318   void ClearEvent();
319   AliReducedEventFriend(const AliReducedEventFriend &c);
320   AliReducedEventFriend& operator= (const AliReducedEventFriend &c);
321
322   ClassDef(AliReducedEventFriend, 1);
323 };
324
325
326 //_________________________________________________________________________
327 class AliReducedCaloCluster : public TObject {
328   
329   friend class AliAnalysisTaskReducedTree;         // friend analysis task which fills the object
330   
331  public:
332   enum ClusterType {
333     kUndefined=0, kEMCAL, kPHOS  
334   };
335    
336   AliReducedCaloCluster();
337   ~AliReducedCaloCluster();
338   
339   Bool_t  IsEMCAL()    const {return (fType==kEMCAL ? kTRUE : kFALSE);}
340   Bool_t  IsPHOS()     const {return (fType==kPHOS ? kTRUE : kFALSE);}
341   Float_t Energy()     const {return fEnergy;}
342   Float_t Dx()         const {return fTrackDx;}
343   Float_t Dz()         const {return fTrackDz;}
344   Float_t M20()        const {return fM20;}
345   Float_t M02()        const {return fM02;}
346   Float_t Dispersion() const {return fDispersion;}
347   
348  private:
349   Char_t  fType;         // cluster type (EMCAL/PHOS)
350   Float_t fEnergy;       // cluster energy
351   Float_t fTrackDx;      // distance to closest track in phi
352   Float_t fTrackDz;      // distance to closest track in z
353   Float_t fM20;          // short axis
354   Float_t fM02;          // long axis
355   Float_t fDispersion;   // dispersion
356   
357   AliReducedCaloCluster(const AliReducedCaloCluster &c);
358   AliReducedCaloCluster& operator= (const AliReducedCaloCluster &c);
359
360   ClassDef(AliReducedCaloCluster, 2);
361 };
362
363
364 //_________________________________________________________________________
365 class AliReducedEvent : public TObject {
366
367   friend class AliAnalysisTaskReducedTree;     // friend analysis task which fills the object
368
369  public:
370   AliReducedEvent();
371   AliReducedEvent(const Char_t* name);
372   ~AliReducedEvent();
373
374   // getters
375   ULong64_t EventTag()                        const {return fEventTag;}
376   Bool_t    EventTag(UShort_t bit)            const {return (bit<8*sizeof(ULong64_t) ? (fEventTag&(ULong64_t(1)<<bit)) : kFALSE);}
377   Int_t     EventNumberInFile()               const {return fEventNumberInFile;}
378   UInt_t    L0TriggerInputs()                 const {return fL0TriggerInputs;}
379   Bool_t    L0TriggerInput(UShort_t bit)         const {return (bit<8*sizeof(UInt_t) ? (fL0TriggerInputs&(UInt_t(1)<<bit)) : kFALSE);}
380   UInt_t    L1TriggerInputs()                 const {return fL1TriggerInputs;}
381   Bool_t    L1TriggerInput(UShort_t bit)         const {return (bit<8*sizeof(UInt_t) ? (fL1TriggerInputs&(UInt_t(1)<<bit)) : kFALSE);}
382   UShort_t  L2TriggerInputs()                 const {return fL2TriggerInputs;}
383   Bool_t    L2TriggerInput(UShort_t bit)         const {return (bit<8*sizeof(UShort_t) ? (fL2TriggerInputs&(UShort_t(1)<<bit)) : kFALSE);}
384   Int_t     RunNo()                           const {return fRunNo;}
385   UShort_t  BC()                              const {return fBC;}
386   UInt_t    TimeStamp()                       const {return fTimeStamp;}
387   UInt_t    EventType()                       const {return fEventType;}
388   ULong64_t TriggerMask()                     const {return fTriggerMask;}
389   Bool_t    IsPhysicsSelection()              const {return fIsPhysicsSelection;}
390   Bool_t    IsSPDPileup()                     const {return fIsSPDPileup;}
391   Bool_t    IsSPDPileupMultBins()             const {return fIsSPDPileupMultBins;}
392   Int_t     IRIntClosestIntMap(Int_t id)      const {return (id>=0 && id<2 ? fIRIntClosestIntMap[id] : -999);}
393   Bool_t    IsFMDReduced()                    const {return fIsFMDReduced;}
394   Float_t   Vertex(Int_t axis)                const {return (axis>=0 && axis<=2 ? fVtx[axis] : 0);}
395   Int_t     VertexNContributors()             const {return fNVtxContributors;}
396   Float_t   VertexTPC(Int_t axis)             const {return (axis>=0 && axis<=2 ? fVtxTPC[axis] : 0);}
397   Int_t     VertexTPCContributors()           const {return fNVtxTPCContributors;}
398   Float_t   VertexTZERO()                     const {return fT0zVertex;}
399   Int_t     NpileupSPD()                      const {return fNpileupSPD;}
400   Int_t     NpileupTracks()                   const {return fNpileupTracks;}
401   Int_t     NPMDtracks()                      const {return fNPMDtracks;}
402   Int_t     NTRDtracks()                      const {return fNTRDtracks;}
403   Int_t     NTRDtracklets()                   const {return fNTRDtracklets;}
404   Float_t   CentralityVZERO()                 const {return fCentrality[0];}
405   Float_t   CentralitySPD()                   const {return fCentrality[1];}
406   Float_t   CentralityTPC()                   const {return fCentrality[2];}
407   Float_t   CentralityZEMvsZDC()              const {return fCentrality[3];}
408   Int_t     CentralityQuality()               const {return fCentQuality;}
409   Int_t     NV0CandidatesTotal()              const {return fNV0candidates[0];}
410   Int_t     NV0Candidates()                   const {return fNV0candidates[1];}
411   Int_t     NDielectrons()                    const {return fNDielectronCandidates;}
412   Int_t     NTracksTotal()                    const {return fNtracks[0];}
413   Int_t     NTracks()                         const {return fNtracks[1];}
414   Int_t     SPDntracklets()                   const {return fSPDntracklets;}
415   Int_t     SPDntracklets(Int_t bin)          const {return (bin>=0 && bin<32 ? fSPDntrackletsEta[bin] : -999);}
416   Int_t     TracksPerTrackingFlag(Int_t flag) const {return (flag>=0 && flag<32 ? fNtracksPerTrackingFlag[flag] : -999);}
417   UShort_t  NFMDchannels(Int_t det)           const {return (det>=0 && det<5 ? fNFMDchannels[det] : 0);}
418   UShort_t  FMDtotalMult(Int_t det)           const {return (det>=0 && det<5 ? fFMDtotalMult[det] : 0);}
419   
420   Float_t   MultChannelVZERO(Int_t channel)   const {return (channel>=0 && channel<=63 ? fVZEROMult[channel] : -999.);}
421   Float_t   MultVZEROA()                      const;
422   Float_t   MultVZEROC()                      const;
423   Float_t   MultVZERO()                       const;
424   Float_t   MultRingVZEROA(Int_t ring)        const;
425   Float_t   MultRingVZEROC(Int_t ring)        const;
426   
427   Float_t   AmplitudeTZEROA()                       const;
428   Float_t   AmplitudeTZEROC()                       const;
429   Float_t   AmplitudeTZERO()                        const;
430   Float_t   AmplitudeTZEROch(Int_t ch)              const {return (ch>=0 && ch<=25 ? fT0amplitude[ch] : -999.);}
431   Float_t   EventTZEROStartTime()                   const {return fT0start;}
432   Float_t   EventTZEROStartTimeTOFfirst(Int_t side) const {return (side>=0 && side<3 ? fT0TOF[side] : -999.);}
433   Float_t   EventTZEROStartTimeTOFbest(Int_t side)  const {return (side>=0 && side<3 ? fT0TOFbest[side] : -999.);}
434   Bool_t    IsPileupTZERO()                         const {return fT0pileup;}
435   Bool_t    IsSatteliteCollisionTZERO()             const {return fT0sattelite;}
436   
437   Float_t   EnergyZDCnTree(UShort_t channel)  const {return (channel<10 ? fZDCnEnergy[channel] : -999.);};
438   Float_t   EnergyZDCpTree(UShort_t channel)  const {return (channel<10 ? fZDCpEnergy[channel] : -999.);};
439   Float_t   EnergyZDCn(Int_t channel)  const;
440   Float_t   EnergyZDCA() const;
441   Float_t   EnergyZDCC() const;
442   Float_t   EnergyZDC() const;
443   
444   Bool_t    TestEventTag(UShort_t iflag)  const {return (iflag<8*sizeof(ULong64_t) ? fEventTag&(ULong64_t(1)<<iflag) : kFALSE);}
445   Bool_t    SetEventTag(UShort_t iflag)         {if (iflag>=8*sizeof(ULong64_t)) return kFALSE; fEventTag|=(ULong64_t(1)<<iflag); return kTRUE;}
446   
447   AliReducedTrack* GetTrack(Int_t i)         const 
448     {return (i<fNtracks[1] ? (AliReducedTrack*)fTracks->At(i) : 0x0);}
449   AliReducedPair* GetV0Pair(Int_t i)         const 
450     {return (i>=0 && i<fNV0candidates[1] ? (AliReducedPair*)fCandidates->At(i) : 0x0);}
451   AliReducedPair* GetDielectronPair(Int_t i) const 
452     {return (i>=0 && i<fNDielectronCandidates ? (AliReducedPair*)fCandidates->At(i+fNV0candidates[1]) : 0x0);}
453   TClonesArray* GetPairs()                              const {return fCandidates;}
454   TClonesArray* GetTracks()                             const {return fTracks;}
455   TClonesArray* GetFMD(Int_t det)                       const {return (det==0 ? fFMD1 :
456                                                                  (det==1 ? fFMD2I :
457                                                                   (det==2 ? fFMD2O :
458                                                                    (det==3 ? fFMD3I :
459                                                                     (det==4 ? fFMD3O :
460                                                                      0)))));}
461   TClonesArray* GetFMD1()                                const {return fFMD1;}
462   TClonesArray* GetFMD2I()                                const {return fFMD2I;}
463   TClonesArray* GetFMD2O()                                const {return fFMD2O;}
464   TClonesArray* GetFMD3I()                                const {return fFMD3I;}
465   TClonesArray* GetFMD3O()                                const {return fFMD3O;}
466   AliReducedFMD* GetFMD1Channel(UShort_t ch) const {return (ch<fNFMDchannels[0] ? (AliReducedFMD*)fFMD1->At(ch) : 0x0);}
467   AliReducedFMD* GetFMD2IChannel(UShort_t ch) const {return (ch<fNFMDchannels[1] ? (AliReducedFMD*)fFMD2I->At(ch) : 0x0);}
468   AliReducedFMD* GetFMD2OChannel(UShort_t ch) const {return (ch<fNFMDchannels[2] ? (AliReducedFMD*)fFMD2O->At(ch) : 0x0);}
469   AliReducedFMD* GetFMD3IChannel(UShort_t ch) const {return (ch<fNFMDchannels[3] ? (AliReducedFMD*)fFMD3I->At(ch) : 0x0);}
470   AliReducedFMD* GetFMD3OChannel(UShort_t ch) const {return (ch<fNFMDchannels[4] ? (AliReducedFMD*)fFMD3O->At(ch) : 0x0);}
471   
472   Int_t GetNCaloClusters() const {return fNCaloClusters;}
473   AliReducedCaloCluster* GetCaloCluster(Int_t i) const 
474     {return (i>=0 && i<fNCaloClusters ? (AliReducedCaloCluster*)fCaloClusters->At(i) : 0x0);}
475   
476   void  GetQvector(Double_t Qvec[][2], Int_t det, Float_t etaMin=-0.8, Float_t etaMax=+0.8, Bool_t (*IsTrackSelected)(AliReducedTrack*)=NULL);
477   Int_t GetTPCQvector(Double_t Qvec[][2], Int_t det, Float_t etaMin=-0.8, Float_t etaMax=+0.8, Bool_t (*IsTrackSelected)(AliReducedTrack*)=NULL);
478   void  GetVZEROQvector(Double_t Qvec[][2], Int_t det) ;
479   void  GetVZEROQvector(Double_t Qvec[][2], Int_t det, Float_t* vzeroMult);
480   void  GetZDCQvector(Double_t Qvec[][2], Int_t det) const;
481   void  GetZDCQvector(Double_t Qvec[][2], Int_t det, const Float_t* zdcEnergy) const;
482   void  SubtractParticleFromQvector(AliReducedTrack* particle, Double_t Qvec[][2], Int_t det,
483                                     Float_t etaMin=-0.8, Float_t etaMax=+0.8,
484                                     Bool_t (*IsTrackSelected)(AliReducedTrack*)=NULL);
485
486  private:
487   ULong64_t fEventTag;              // Event tags to be used either during analysis or to filter events
488   Int_t     fEventNumberInFile;     // Event number in ESD file
489   UInt_t    fL0TriggerInputs;       // L0 trigger inputs
490   UInt_t    fL1TriggerInputs;       // L1 trigger inputs
491   UShort_t  fL2TriggerInputs;       // L2 trigger inputs
492   Int_t     fRunNo;                 // run number
493   UShort_t  fBC;                    // bunch crossing
494   UInt_t    fTimeStamp;             // time stamp of the event                (NEW)
495   UInt_t    fEventType;             // event type                             (NEW)
496   ULong64_t fTriggerMask;           // trigger mask
497   Bool_t    fIsPhysicsSelection;    // PhysicsSelection passed event
498   Bool_t    fIsSPDPileup;           // identified as pileup event by SPD
499   Bool_t    fIsSPDPileupMultBins;   // identified as pileup event by SPD in multiplicity bins
500   Int_t     fIRIntClosestIntMap[2]; // out of bunch interactions, [0]-Int1, [1]-Int2 
501   Bool_t    fIsFMDReduced;          // FMD info, if present, is reduced       (NEW)
502   Float_t   fVtx[3];                // global event vertex vector in cm
503   Int_t     fNVtxContributors;      // global event vertex contributors
504   Float_t   fVtxTPC[3];             // TPC only event vertex           
505   Int_t     fNVtxTPCContributors;   // TPC only event vertex contributors
506   Int_t     fNpileupSPD;            // number of pileup vertices from SPD     (NEW)
507   Int_t     fNpileupTracks;         // number of pileup vertices from tracks  (NEW)
508   Int_t     fNPMDtracks;            // number of PMD tracks                   (NEW)
509   Int_t     fNTRDtracks;            // number of TRD tracks                   (NEW)
510   Int_t     fNTRDtracklets;         // number of TRD tracklets                (NEW)
511   Float_t   fCentrality[4];         // centrality; 0-VZERO, 1-SPD, 2-TPC, 3-ZEMvsZDC 
512   Int_t     fCentQuality;           // quality flag for the centrality 
513   Int_t     fNV0candidates[2];      // number of V0 candidates, [0]-total, [1]-selected for the tree
514   Int_t     fNDielectronCandidates; // number of pairs selected as dielectrons
515   Int_t     fNtracks[2];            // number of tracks, [0]-total, [1]-selected for the tree
516   Int_t     fSPDntracklets;         // number of SPD tracklets in |eta|<1.0 
517   Int_t     fSPDntrackletsEta[32];  // number of SPD tracklets in equal eta bins between -1.6 --> +1.6    
518   Int_t     fNtracksPerTrackingFlag[32];  // number of tracks for each tracking status bit                
519   
520   Float_t   fVZEROMult[64];         // VZERO multiplicity in all 64 channels
521   Float_t   fZDCnEnergy[10];         // neutron ZDC energy in all 8 channels
522   Float_t   fZDCpEnergy[10];         // neutron ZDC energy in all 8 channels
523   Float_t   fT0amplitude[26];         // T0 amplitude in all 24 channels
524   Float_t   fT0TOF[3];               // T0 timing for A&C, A, and C (first time)
525   Float_t   fT0TOFbest[3];           // T0 timing for A&C, A, and C (best time)
526   Float_t   fT0zVertex;                // T0 z vertex estimation
527   Float_t   fT0start;                // T0 timing
528   Bool_t    fT0pileup;               // TZERO pileup flag
529   Bool_t    fT0sattelite;            // TZERO flag for collisions from sattelite bunches
530     
531   TClonesArray* fTracks;            //->   array containing global tracks
532   static TClonesArray* fgTracks;    //       global tracks
533   
534   TClonesArray* fCandidates;        //->   array containing pair candidates
535   static TClonesArray* fgCandidates;  // pair candidates
536   
537   UShort_t   fNFMDchannels[5];           // number of FMD channels read out           (NEW)
538   UShort_t   fFMDtotalMult[5];           // number of FMD channels read out           (NEW)
539   TClonesArray* fFMD1;            //->   array containing fmd readout          (NEW)
540   static TClonesArray* fgFMD1;    //       fmd readout                        (NEW)
541   TClonesArray* fFMD2I;            //->   array containing fmd readout          (NEW)
542   static TClonesArray* fgFMD2I;    //       fmd readout                       (NEW)
543   TClonesArray* fFMD2O;            //->   array containing fmd readout          (NEW)
544   static TClonesArray* fgFMD2O;    //       fmd readout                       (NEW)
545   TClonesArray* fFMD3I;            //->   array containing fmd readout          (NEW)
546   static TClonesArray* fgFMD3I;    //       fmd readout                       (NEW)
547   TClonesArray* fFMD3O;            //->   array containing fmd readout          (NEW)
548   static TClonesArray* fgFMD3O;    //       fmd readout                       (NEW)
549   
550   Int_t     fNCaloClusters;         // number of calorimeter clusters  
551   TClonesArray* fCaloClusters;        //->   array containing calorimeter clusters
552   static TClonesArray* fgCaloClusters;     // calorimeter clusters
553   
554   void ClearEvent();
555   AliReducedEvent(const AliReducedEvent &c);
556   AliReducedEvent& operator= (const AliReducedEvent &c);
557
558   ClassDef(AliReducedEvent, 5);
559 };
560
561 //_______________________________________________________________________________
562 inline UShort_t AliReducedTrack::ITSncls() const
563 {
564   //
565   // ITS number of clusters from the cluster map
566   //
567   UShort_t ncls=0;
568   for(Int_t i=0; i<6; ++i) ncls += (ITSLayerHit(i) ? 1 : 0);
569   return ncls;
570 }
571
572
573 //_______________________________________________________________________________
574 inline Int_t AliReducedTrack::TPCClusterMapBitsFired()  const
575 {
576   //
577   // Count the number of bits fired in the TPC cluster map
578   //
579   Int_t nbits=0;
580   for(Int_t i=0; i<8; ++i) nbits += (TPCClusterMapBitFired(i) ? 1 : 0);
581   return nbits;
582 }
583
584
585 //_______________________________________________________________________________
586 inline Float_t AliReducedPair::Energy() const 
587 {
588   //
589   // Return the energy
590   //
591   Float_t mass=fMass[0];
592   switch (fCandidateId) {
593     case kK0sToPiPi:
594       mass = fMass[0];
595       break;
596     case kLambda0ToPPi:
597       mass = fMass[1];
598       break;
599     case kALambda0ToPPi:
600       mass = fMass[2];
601       break;
602     case kGammaConv:
603       mass = fMass[3];
604       break;
605     default:
606       mass = fMass[0];
607       break;    
608   }
609   Float_t p = P();
610   return TMath::Sqrt(mass*mass+p*p);
611 }
612
613
614 //_______________________________________________________________________________
615 inline Float_t AliReducedPair::Rapidity() const
616 {
617   //
618   // return rapidity
619   //
620   Float_t e = Energy();
621   Float_t pz = Pz();
622   if(e-TMath::Abs(pz)>1.0e-10)
623     return 0.5*TMath::Log((e+pz)/(e-pz));
624   else 
625     return -999.;
626 }
627
628
629 //_______________________________________________________________________________
630 inline Double_t AliReducedEventFriend::EventPlane(Int_t det, Int_t harmonic) const
631 {
632   //
633   // Event plane from detector "det" and harmonic "harmonic"
634   //
635   if(det<0 || det>=kNdetectors || harmonic<1 || harmonic>fgkNMaxHarmonics) return -999.;
636   return TMath::ATan2(fQvector[det][harmonic-1][1], fQvector[det][harmonic-1][0])/Double_t(harmonic);
637 }
638
639 //_______________________________________________________________________________
640 inline Bool_t AliReducedEventFriend::CheckEventPlaneStatus(Int_t det, Int_t h, EventPlaneStatus flag) const {
641   //
642   // Check the status of the event plane for a given detector and harmonic
643   //
644   if(det<0 || det>=kNdetectors || h<1 || h>fgkNMaxHarmonics) return kFALSE;
645   return (flag<kNMaxFlowFlags ? (fEventPlaneStatus[det][h-1]&(1<<flag)) : kFALSE);
646 }
647
648 #endif