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