]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/CaloTrackCorrelations/AliAnaParticleJetFinderCorrelation.h
move the decay bit definition from AliAnaPi0EbE to AliNeutralSelection, define there...
[u/mrichter/AliRoot.git] / PWGGA / CaloTrackCorrelations / AliAnaParticleJetFinderCorrelation.h
1 #ifndef ALIANAPARTICLEJETFINDERCORRELATION_H
2 #define ALIANAPARTICLEJETFINDERCORRELATION_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice     */
5
6 //_________________________________________________________________________
7 // Class that contains the algorithm for the analysis of particle (direct gamma) - jet 
8 // (standard jet found with JETAN) correlation 
9 // Particle and jet for correlation found by independent algorithms.
10 // For Example direct isolated photon found in AliAnaGammaDirect and the jet with JETAN
11 //
12 //-- Author: Gustavo Conesa (INFN-LNF)
13 //-- Modified by Adam Matyja (INP PAN, Krakow)
14
15 // --- ROOT system ---
16 class TH2F;
17 class TTree;
18 class TRandom2;
19 //---- Analysis system ----
20 #include "AliAnaCaloTrackCorrBaseClass.h"
21
22 class AliAnaParticleJetFinderCorrelation : public AliAnaCaloTrackCorrBaseClass {
23        
24  public:   
25   AliAnaParticleJetFinderCorrelation() ;              // default ctor
26   virtual ~AliAnaParticleJetFinderCorrelation() ; // virtual dtor
27
28   // General methods
29   
30   void     InitParameters();
31   
32   TList *  GetCreateOutputObjects();
33
34   void     MakeAnalysisFillAOD() ;
35   
36   void     MakeAnalysisFillHistograms() ;
37   
38   Int_t    SelectJet(AliAODPWG4Particle * particle, TClonesArray * aodRecJets) ;//to access non standard branch
39
40   void     Print(const Option_t * opt)           const;
41   
42   // Settings
43   
44   Bool_t   OnlyIsolated()                        const { return fSelectIsolated              ; }
45   void     SelectIsolated(Bool_t select)               { fSelectIsolated = select            ; }
46   
47   Float_t  GetConeSize()                         const { return fConeSize                    ; }
48   Float_t  GetPtThresholdInCone()                const { return fPtThresholdInCone           ; }           
49   Double_t GetDeltaPhiMaxCut()                   const { return fDeltaPhiMaxCut              ; }
50   Double_t GetDeltaPhiMinCut()                   const { return fDeltaPhiMinCut              ; }
51   Double_t GetRatioMaxCut()                      const { return fRatioMaxCut                 ; }
52   Double_t GetRatioMinCut()                      const { return fRatioMinCut                 ; }           
53   Bool_t   AreJetRefTracks()                     const { return fUseJetRefTracks             ; }
54   Bool_t   IsCorrelationMadeInHistoMaker()       const { return fMakeCorrelationInHistoMaker ; } 
55   Double_t GetJetConeSize()                      const { return fJetConeSize                 ; } 
56   Double_t GetJetMinPt()                         const { return fJetMinPt                    ; }
57   Double_t GetJetAreaFraction()                  const { return fJetAreaFraction             ; }
58   Double_t GetGammaConeSize()                    const { return fGammaConeSize               ; }
59
60   void     SetConeSize(Float_t cone)                   { fConeSize = cone                    ; }
61   void     SetPtThresholdInCone(Float_t pt)            { fPtThresholdInCone = pt             ; }           
62   void     SetDeltaPhiCutRange(Double_t phimin, Double_t phimax)
63             { fDeltaPhiMaxCut =phimax;  fDeltaPhiMinCut =phimin                              ; }
64   void     SetRatioCutRange(Double_t ratiomin, Double_t ratiomax)
65             { fRatioMaxCut =ratiomax;  fRatioMinCut = ratiomin                               ; }
66   void     UseJetRefTracks(Bool_t use)                 { fUseJetRefTracks = use              ; }        
67   void     SetMakeCorrelationInHistoMaker(Bool_t make) { fMakeCorrelationInHistoMaker = make ; }
68   void     SetJetConeSize(Double_t cone)               { fJetConeSize = cone ;                 }
69   void     SetJetMinPt(Double_t minpt)                 { fJetMinPt = minpt ;                   }
70   void     SetJetAreaFraction(Double_t areafr)         { fJetAreaFraction = areafr ;           }
71   void     SetGammaConeSize(Float_t cone)              { fGammaConeSize = cone               ; }
72
73   // Settings for non standard jet branch
74   TString  GetJetBranchName()            const { return fJetBranchName             ; }
75   void     SetJetBranchName(const char *name)  { fJetBranchName = name             ; }
76   //void     SwitchOnNonStandardJetFromReader()  { fNonStandardJetFromReader = kTRUE ; }
77   //void     SwitchOffNonStandardJetFromReader() { fNonStandardJetFromReader = kFALSE; }
78   //Bool_t   IsNonStandardJetFromReader()        { return fNonStandardJetFromReader  ; }
79
80   TString  GetBkgJetBranchName()           const { return fBkgJetBranchName          ; }
81   void     SetBkgJetBranchName(const char *name) { fBkgJetBranchName = name          ; }
82   void     SwitchOnBackgroundJetFromReader()     { fBackgroundJetFromReader = kTRUE  ; }
83   void     SwitchOffBackgroundJetFromReader()    { fBackgroundJetFromReader = kFALSE ; }
84   Bool_t   IsBackgroundJetFromReader()           { return fBackgroundJetFromReader   ; }
85
86   //switches for photons
87   void     SwitchOnBackgroundSubtractionGamma()  { fUseBackgroundSubtractionGamma = kTRUE ; }
88   void     SwitchOffBackgroundSubtractionGamma() { fUseBackgroundSubtractionGamma = kFALSE; }
89   Bool_t   IsBackgroundSubtractionGamma()        { return fUseBackgroundSubtractionGamma  ; }
90
91   void     CalculateBkg(TVector3 gamma, TVector3 jet,Double_t *vector,Int_t type);
92
93   void     FindMCgenInfo();//gives information on generated level
94
95   void     SwitchOnSaveGJTree()     { fSaveGJTree = kTRUE ; }
96   void     SwitchOffSaveGJTree()    { fSaveGJTree = kFALSE; }
97   Bool_t   IsSaveGJTree()           { return fSaveGJTree  ; }
98
99   void     SwitchOnMostEnergetic()  { fMostEnergetic = kTRUE ; fMostOpposite = kFALSE; }
100   void     SwitchOffMostEnergetic() { fMostEnergetic = kFALSE; fMostOpposite = kTRUE ; }
101   void     SwitchOffMostOpposite()  { fMostEnergetic = kTRUE ; fMostOpposite = kFALSE; }
102   void     SwitchOnMostOpposite()   { fMostEnergetic = kFALSE; fMostOpposite = kTRUE ; }
103   Bool_t   IsMostEnergetic()        { return fMostEnergetic ; }
104   Bool_t   IsMostOpposite()         { return fMostOpposite; }
105
106   //switches for histograms
107   void     SwitchOnHistogramJetBkg()     { fUseHistogramJetBkg = kTRUE ; }
108   void     SwitchOffHistogramJetBkg()    { fUseHistogramJetBkg = kFALSE; }
109   Bool_t   IsHistogramJetBkg()           { return fUseHistogramJetBkg  ; }
110
111   void     SwitchOnHistogramTracks()     { fUseHistogramTracks = kTRUE ; }
112   void     SwitchOffHistogramTracks()    { fUseHistogramTracks = kFALSE; }
113   Bool_t   IsHistogramTracks()           { return fUseHistogramTracks  ; }
114
115   void     SwitchOnHistogramJetTracks()  { fUseHistogramJetTracks = kTRUE ; }
116   void     SwitchOffHistogramJetTracks() { fUseHistogramJetTracks = kFALSE; }
117   Bool_t   IsHistogramJetTracks()        { return fUseHistogramJetTracks  ; }
118
119   void     SwitchOnMCStudies()     { fMCStudies = kTRUE ; }
120   void     SwitchOffMCStudies()    { fMCStudies = kFALSE; }
121   Bool_t   IsMCStudies()           { return fMCStudies  ; }
122
123 private:
124
125   //selection parameters  
126   Double_t   fDeltaPhiMaxCut ;    //! Minimum Delta Phi Gamma-Leading
127   Double_t   fDeltaPhiMinCut ;    //! Maximum Delta Phi Gamma-Leading
128   Double_t   fRatioMaxCut ;       //! Jet/particle Ratio cut maximum
129   Double_t   fRatioMinCut ;       //! Jet/particle Ratio cut minimum
130   
131   Double_t   fConeSize  ;         //! Jet cone size to calculate fragmentation function
132   Double_t   fPtThresholdInCone ; //! Jet pT threshold in jet cone
133   Bool_t     fUseJetRefTracks ;   //! Use track references from JETAN not the AOD tracks to calculate fragmentation function
134   Bool_t     fMakeCorrelationInHistoMaker ; //!Make particle-jet correlation in histogram maker
135   Bool_t     fSelectIsolated ;    //! Select only trigger particles isolated
136   
137   Double_t   fJetConeSize ;       //! Reconstructed jet cone size 
138   Double_t   fJetMinPt ;          //! Minumum jet pt, default 5GeV/c
139   Double_t   fJetAreaFraction ;   //! Jet area fraction X in X*pi*R^2, default 0.6
140   //Bool_t     fNonStandardJetFromReader; //! use non standard jet from reader //new
141   TString    fJetBranchName ;//! name of jet branch not set in reader part //new
142   Bool_t     fBackgroundJetFromReader; //! use background jet from reader //new
143   TString    fBkgJetBranchName ;//! name of background jet branch not set in reader part //new
144
145   Double_t   fGammaConeSize ;       //! Isolation cone radius
146   Bool_t     fUseBackgroundSubtractionGamma;//! flag to use backgrouind subtraction for photons or not
147   Bool_t     fSaveGJTree;//! flag to save gamma-jet tree
148   Bool_t     fMostEnergetic;//! flag to choose gamma-jet pairs most energetic
149   Bool_t     fMostOpposite;//! flag to choose gamma-jet pairs most opposite
150
151   Bool_t     fUseHistogramJetBkg;//! flag to save bkg jet histograms
152   Bool_t     fUseHistogramTracks;//! flag to save CTS tracks features
153   Bool_t     fUseHistogramJetTracks;//! flag to save jet tracks features
154
155   Bool_t     fMCStudies; //! flag to use MC methods
156
157   TRandom2 * fGenerator;//! pointer to random generator object
158
159   // Histograms
160   TH2F *     fhDeltaEta;           //! Difference of jet eta and trigger particle eta as function of trigger particle pT
161   //TH2F *     fhDeltaPhi;         //! Difference of jet phi and trigger particle phi as function of trigger particle pT
162   TH2F *     fhDeltaPhiCorrect;    //! Difference of jet phi and trigger particle phi as function of trigger particle pT
163   TH2F *     fhDeltaPhi0PiCorrect; //! Difference of jet phi and trigger particle phi as function of trigger particle pT
164
165   TH2F *     fhDeltaPt;           //! Difference of jet pT and trigger particle pT as function of trigger particle pT
166   TH2F *     fhPtRatio;           //! Ratio of jet pT and trigger particle pT as function of trigger particle pT
167   TH2F *     fhPt;                //! jet pT vs trigger particle pT 
168   
169   TH2F *     fhFFz ;              //! Accepted reconstructed jet fragmentation function, z=pt^particle,jet/pttrig
170   TH2F *     fhFFxi;              //! Accepted reconstructed jet fragmentation function, xsi = ln(pttrig/pt^particle,jet)
171   TH2F *     fhFFpt;              //! Jet particle pt distribution in cone
172   TH2F *     fhNTracksInCone;     //! jet multiplicity in cone
173   
174   TH2F *     fhJetFFz ;           //! Accepted reconstructed jet fragmentation function, z=pt^particle,jet/ptjet
175   TH2F *     fhJetFFxi;           //! Accepted reconstructed jet fragmentation function, xsi = ln(ptjet/pt^particle,jet)
176   TH2F *     fhJetFFpt;           //! Jet particle pt distribution in jet cone
177   TH2F *     fhJetFFzCor ;        //! Accepted reconstructed jet fragmentation function, z=pt^particle,jet*-cos(jet,trig)/ptjet
178   TH2F *     fhJetFFxiCor;        //! Accepted reconstructed jet fragmentation function, xsi = ln(ptjet/pt^particle*-cos(jet,trig),jet)
179
180   //background from RC
181   TH2F *     fhBkgFFz[5] ;              //! Background fragmentation function, z=ptjet/pttrig
182   TH2F *     fhBkgFFxi[5];              //! Background fragmentation function, xsi = ln(pttrig/ptjet)
183   TH2F *     fhBkgFFpt[5];              //! Background particle pt distribution in cone
184   TH2F *     fhBkgNTracksInCone[5];     //! Background multiplicity in cone
185   TH2F *     fhBkgSumPtInCone[5];       //! Background sum pt in cone
186   TH2F *     fhBkgSumPtnTracksInCone[5];//! Background sum pt over ntracks in cone
187
188   //temporary histograms
189   TH2F * fhNjetsNgammas; //! Number of jets vs number of photons in the event
190   TH1F * fhCuts;         //! Number of events after cuts
191
192   TH2F * fhDeltaEtaBefore;  //! Difference of jet eta and trigger particle eta as function of trigger particle pT
193   TH2F * fhDeltaPhiBefore;  //! Difference of jet phi and trigger particle phi as function of trigger particle pT
194   TH2F * fhDeltaPtBefore;   //! Difference of jet pT and trigger particle pT as function of trigger particle pT
195   TH2F * fhPtRatioBefore;   //! Ratio of jet pT and trigger particle pT as function of trigger particle pT
196   TH2F * fhPtBefore;        //! jet pT vs trigger particle pT
197   TH2F * fhDeltaPhi0PiCorrectBefore; //! Difference of jet phi and trigger particle phi (0,pi) as function of trigger particle pT
198
199   //temporary jet histograms
200   TH1F * fhJetPtBefore;           //! Pt of all jets
201   TH1F * fhJetPt;                 //! Pt of all jets after bkg correction
202   TH1F * fhJetPtMostEne;          //! Pt of the most energetic jet
203   TH1F * fhJetPhi;                    //! Phi of all jets
204   TH1F * fhJetEta;                    //! Eta of all jets
205   TH2F * fhJetEtaVsPt;            //! Eta of all jets vs pt
206   TH2F * fhJetPhiVsEta;           //! Phi vs eta of all jets
207   TH2F * fhJetEtaVsNpartInJet;    //! Eta vs number of particles in jet for all jets
208   TH2F * fhJetEtaVsNpartInJetBkg; //! Eta vs number of particles in jet for background subtracted jets
209   TH2F * fhJetChBkgEnergyVsPt;    //! background energy of each charged jet vs jet pt
210   TH2F * fhJetChAreaVsPt;         //! area of each charged jet vs jet pt
211   TH2F * fhTrackPhiVsEta;         //! Phi vs eta of all chosen tracks in all events
212   TH1F * fhTrackAveTrackPt;       //! average track pt in event
213   TH1F * fhJetNjetOverPtCut[10];  //! number of reconstructed jets in event over pT threshold
214   TH2F * fhJetChBkgEnergyVsArea;  //! area of each charged jet vs jet background
215   TH2F * fhJetRhoVsPt;            //! jet energy density vs jet pt
216   TH2F * fhJetRhoVsCentrality;    //! jet energy density vs centrality
217   TH1F * fhJetNparticlesInJet;    //! number of particles in jets
218   TH2F * fhJetDeltaEtaDeltaPhi;   //! delta eta vs delta phi for (jet-track) <-0.8,0.8>
219   TH2F * fhJetDeltaEtaDeltaPhiAllTracks;//! delta eta vs delta phi for (jet-track) <-pi,pi>
220
221
222
223   TH1F * fhJetAveTrackPt;              //! average track from jets pt in event
224   TH2F * fhJetNtracksInJetAboveThr[6]; //! number of tracks in jet with pt above 0,1,2,3,4,5GeV
225   TH2F * fhJetRatioNTrkAboveToNTrk[5]; //! ratio tracks in jet with pt above 1,2,3,4,5GeV to ntracks
226   TH2F * fhJetNtrackRatioMostEne[5];   //! the same for most energetic jet
227   TH2F * fhJetNtrackRatioJet5GeV[5];   //! the same for pt jet above 5 GeV
228   TH2F * fhJetNtrackRatioLead5GeV[5];  //! the same for jet with leading particle pt>5GeV
229
230   //temporary background jet histograms
231   TH1F * fhBkgJetBackground[4]; //! background from jet bkg branch
232   TH1F * fhBkgJetSigma[4];      //! sigma of jet in backgroud branch
233   TH1F * fhBkgJetArea[4];       //! area of jet in bkg branch
234
235   //temporary photon histograms
236   TH1F * fhPhotonPtMostEne;                       //! most pt photon
237   TH1F * fhPhotonAverageEnergy;                   //! average energy of photon
238   TH1F * fhPhotonRatioAveEneToMostEne;            //! ratio average energy to most energetic photon
239   TH1F * fhPhotonAverageEnergyMinus1;             //! average energy of photon w/o most ene photon
240   TH1F * fhPhotonRatioAveEneMinus1ToMostEne;      //! ratio average energy of photon w/o most ene photon to most energetic photon
241   TH1F * fhPhotonNgammaMoreAverageToNgamma;       //! number of gammas with ene. more than average ene divided by no. of gammas
242   TH1F * fhPhotonNgammaMoreAverageMinus1ToNgamma; //! number of gammas with ene. more than average ene (w/o most ene gamma) divided by no. of gammas
243   TH1F * fhPhotonNgammaOverPtCut[10];             //! number of photons in event over pT threshold
244   TH2F * fhPhotonBkgRhoVsNtracks;                 //! average energy in one cell vs n tracks
245   TH2F * fhPhotonBkgRhoVsNclusters;               //! average energy in one cell vs n clusters
246   TH2F * fhPhotonBkgRhoVsCentrality;              //! average energy in one cell vs centrality
247   TH2F * fhPhotonBkgRhoVsNcells;                  //! average energy in one cell vs n cells
248   TH1F * fhPhotonPt;                              //! pt of gamma before bkg correction
249   TH1F * fhPhotonPtCorrected;                     //! pt of gamma after background correction
250   TH1F * fhPhotonPtCorrectedZoom;                 //! pt of gamma after background correction in +-5 GeV/c
251   TH1F * fhPhotonPtDiff;                          //! bkg correction = n_cells * median_rho
252   TH2F * fhPhotonPtDiffVsCentrality;              //! correction vs centrality
253   TH2F * fhPhotonPtDiffVsNcells;                  //! correction vs Ncells
254   TH2F * fhPhotonPtDiffVsNtracks;                 //! correction vs Ntracks
255   TH2F * fhPhotonPtDiffVsNclusters;               //! correction vs Nclustres
256
257   TH1F * fhPhotonSumPtInCone;                     //! sum pt in cone before correction
258   TH1F * fhPhotonSumPtCorrectInCone;              //! sum pt in cone afrer correction
259   TH1F * fhPhotonSumPtChargedInCone;              //! sum pt of charged tracks in the cone before correction
260
261
262   //temporary jet histograms after selection
263   TH2F * fhSelectedJetPhiVsEta;            //! phi vs eta of selected jet
264   TH2F * fhSelectedJetChBkgEnergyVsPtJet;  //! background energy of selected charged jet vs jet pt
265   TH2F * fhSelectedJetChAreaVsPtJet;       //! area of selected charged jet vs jet pt
266   TH1F * fhSelectedJetNjet;                //! number of jets in selected event
267   TH1F * fhSelectedNtracks;                //! number of tracks in selected event
268   TH2F * fhSelectedTrackPhiVsEta;          //! Phi vs eta of all chosen tracks in selected events
269
270   TH1F * fhCuts2;                          //! efficienct cuts
271
272   //temporary photon histogram after selection
273   TH2F * fhSelectedPhotonNLMVsPt;          //! nlm vs pt for selected photons
274   TH2F * fhSelectedPhotonLambda0VsPt;      //! lambda0 vs pt for selected photons
275   TH2F * fhRandomPhiEta[5];                //! eta and phi from random generator
276
277   //MC generated histograms
278   TH1F * fhMCPhotonCuts;  //! generated photon cuts
279   TH1F * fhMCPhotonPt;    //! generated direct photon pt
280   TH2F * fhMCPhotonEtaPhi;//! generated direct photon eta vs phi
281   TH1F * fhMCJetOrigin;   //! generated origin of jet
282   TH2F * fhMCJetNPartVsPt;     //! generated N parts vs pt full jet 
283   TH2F * fhMCJetChNPartVsPt;   //! generated N parts vs pt charged jet 
284   TH2F * fhMCJetNPart150VsPt;  //! generated N parts (pt>150 MeV/c) vs pt full jet 
285   TH2F * fhMCJetChNPart150VsPt;//! generated N parts (pt>150 MeV/c) vs pt charged jet 
286   TH2F * fhMCJetChNPart150ConeVsPt;//! generated N parts (pt>150 MeV/c) vs pt charged jet R=0.4 
287   TH1F * fhMCJetRatioChFull;   //! generated ratio pt charged/full jet
288   TH1F * fhMCJetRatioCh150Ch;  //! generated ratio pt charged(pt>150MeV/c)/charged jet
289   TH2F * fhMCJetEtaPhi;     //! generated jet eta vs phi for full jet
290   TH2F * fhMCJetChEtaPhi;   //! generated jet eta vs phi for charged jet
291   TH2F * fhMCJet150EtaPhi;  //! generated jet eta vs phi full jet (pt>150 MeV/c)
292   TH2F * fhMCJetCh150EtaPhi;//! generated jet eta vs phi charged jet (pt>150 MeV/c)
293   TH2F * fhMCJetCh150ConeEtaPhi;//! generated jet eta vs phi charged jet (pt>150 MeV/c) R=0.4 
294
295   //tree with data gamma and jet
296   TTree *  fTreeGJ;       //! gamma-jet tree
297   Double_t fGamPt;        //! pt
298   Double_t fGamLambda0;   //! lambda 0
299   Int_t    fGamNLM;       //! NLM
300   Double_t fGamSumPtCh;   //! energy in isolation cone charged
301   Double_t fGamTime;      //! time
302   Int_t    fGamNcells;    //! ncells
303   Double_t fGamEta;       //! eta photon
304   Double_t fGamPhi;       //! phi photon
305   Double_t fGamSumPtNeu;  //! energy in isolation cone neutral
306   Int_t    fGamNtracks;   //! number of tracks in iso cone
307   Int_t    fGamNclusters; //! number of clusters in iso cone
308   Double_t fGamAvEne;     //! average energy of photons (without most ene)
309   Double_t fJetPhi;       //! jet phi
310   Double_t fJetEta;       //! eta phi
311   Double_t fJetPt;        //! jet pt
312   Double_t fJetBkgChEne;  //! bkg energy of jet
313   Double_t fJetArea;      //! jet area
314   Int_t    fJetNtracks;   //! number of jet tracks
315   Int_t    fJetNtracks1;  //! number of jet tracks with pt>1 GeV/c
316   Int_t    fJetNtracks2;  //! number of jet tracks with pt>2 GeV/c
317   Double_t fJetRho;       //! jet rho in event
318   Int_t    fEventNumber;  //! event number
319   Int_t    fNtracks;      //! n tracks in event
320   Double_t fZvertex;      //! z vertex
321   Double_t fCentrality;   //! centrality
322   Bool_t   fIso;          //! flag isolated or not
323   Double_t fGamRho;       //! background energy for photons per cell in EMCal
324
325   Double_t fMCGamPt;        //! MC gen pt photon
326   Double_t fMCGamEta;       //! MC gen eta photon
327   Double_t fMCGamPhi;       //! MC gen phi photon
328   Int_t    fMCPartonType;   //! MC gen parton type origin of jet
329   Double_t fMCJetPt;        //! MC gen full jet pt
330   Double_t fMCJetChPt;      //! MC gen charged jet pt
331   Double_t fMCJet150Pt;     //! MC gen full jet (pt^particles>150MeV/c) pt
332   Double_t fMCJetCh150Pt;   //! MC gen charged jet (pt^particles>150MeV/c) pt
333   Int_t    fMCJetNPart;     //! MC gen number of full jet particles
334   Int_t    fMCJetChNPart;   //! MC gen number of charged jet particles
335   Int_t    fMCJet150NPart;  //! MC gen number of full jet particles (pt>150MeV/c)
336   Int_t    fMCJetCh150NPart;//! MC gen number of charged jet particles (pt>150MeV/c)
337   Double_t fMCJetEta;       //! MC gen full jet eta
338   Double_t fMCJetPhi;       //! MC gen full jet phi
339   Double_t fMCJetChEta;     //! MC gen charged jet eta
340   Double_t fMCJetChPhi;     //! MC gen charged jet phi
341   Double_t fMCJet150Eta;    //! MC gen full jet eta (pt>150MeV/c)
342   Double_t fMCJet150Phi;    //! MC gen full jet phi (pt>150MeV/c)
343   Double_t fMCJetCh150Eta;  //! MC gen charged jet eta (pt>150MeV/c)
344   Double_t fMCJetCh150Phi;  //! MC gen charged jet phi (pt>150MeV/c)
345
346   Double_t fMCJetCh150ConePt;//! MC gen charged jet (pt^particles>150MeV/c),R=0.4 pt
347   Int_t    fMCJetCh150ConeNPart;//! MC gen number of charged jet particles (pt>150MeV/c),R=0.4
348   Double_t fMCJetCh150ConeEta;  //! MC gen charged jet eta (pt>150MeV/c),R=0.4
349   Double_t fMCJetCh150ConePhi;  //! MC gen charged jet phi (pt>150MeV/c),R=0.4
350
351   AliAnaParticleJetFinderCorrelation(              const AliAnaParticleJetFinderCorrelation & g) ; // cpy ctor
352   AliAnaParticleJetFinderCorrelation & operator = (const AliAnaParticleJetFinderCorrelation & g) ; // cpy assignment
353   
354   ClassDef(AliAnaParticleJetFinderCorrelation,3)
355   
356  } ;
357
358 #endif //ALIANAPARTICLEJETFINDERCORRELATION_H
359
360
361