]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGDQ/dielectron/AliDielectronVarManager.h
o update package
[u/mrichter/AliRoot.git] / PWGDQ / dielectron / AliDielectronVarManager.h
1 #ifndef ALIDIELECTRONVARMANAGER_H
2 #define ALIDIELECTRONVARMANAGER_H
3 /* Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice                               */
5
6 //#############################################################
7 //#                                                           # 
8 //#         Class AliDielectronVarManager                     #
9 //#         Class for management of available variables       #
10 //#                                                           #
11 //#  Authors:                                                 #
12 //#   Anton     Andronic, GSI / A.Andronic@gsi.de             #
13 //#   Ionut C.  Arsene,   GSI / I.C.Arsene@gsi.de             #
14 //#   Julian    Book,     Uni Ffm / Julian.Book@cern.ch       #
15 //#   Markus    Köhler,   GSI / M.Koehler@gsi.de              #
16 //#   Frederick Kramer,   Uni Ffm / Frederick.Kramer@cern.ch  #
17 //#   Magnus    Mager,    CERN / Magnus.Mager@cern.ch         #
18 //#   WooJin J. Park,     GSI / W.J.Park@gsi.de               #
19 //#   Jens      Wiechula, Uni HD / Jens.Wiechula@cern.ch      #
20 //#                                                           #
21 //#############################################################
22
23 #include <TNamed.h>
24 #include <TProfile.h>
25 #include <TProfile2D.h>
26 #include <TH3D.h>
27 #include <TFile.h>
28 #include <TDatabasePDG.h>
29 #include <TKey.h>
30
31 #include <AliVEvent.h>
32 #include <AliESDEvent.h>
33 #include <AliAODEvent.h>
34 #include <AliMCEvent.h>
35 #include <AliVVertex.h>
36 #include <AliESDVertex.h>
37 #include <AliEventplane.h>
38
39 #include <AliESDVZERO.h>
40 #include <AliAODVZERO.h>
41
42 #include <AliVParticle.h>
43 #include <AliESDtrack.h>
44 #include <AliAODTrack.h>
45 #include <AliAODPid.h>
46 #include <AliKFParticle.h>
47 #include <AliKFVertex.h>
48 #include <AliMCParticle.h>
49 #include <AliAODMCParticle.h>
50 #include <AliVTrack.h>  // ?
51
52 #include <AliExternalTrackParam.h>
53 #include <AliESDpid.h>
54 #include <AliCentrality.h>
55 #include <AliAODpidUtil.h>
56 #include <AliPID.h>
57 #include <AliPIDResponse.h>
58
59 #include "AliDielectronPair.h"
60 #include "AliDielectronMC.h"
61 #include "AliDielectronPID.h"
62 #include "AliDielectronHelper.h"
63
64 class AliVEvent;
65
66 //________________________________________________________________
67 class AliDielectronVarManager : public TNamed {
68   
69 public:
70
71   // Particle specific variables
72   enum ValueTypes {
73     kPx = 0,                 // px
74     kPy,                     // py
75     kPz,                     // pz
76     kPt,                     // transverse momentum
77     kP,                      // momentum
78     kXv,                     // vertex position in x
79     kYv,                     // vertex position in y
80     kZv,                     // vertex position in z
81     kOneOverPt,              // 1/pt
82     kPhi,                    // phi angle
83     kTheta,                  // theta angle
84     kEta,                    // pseudo-rapidity
85     kY,                      // rapidity
86     kE,                      // energy
87     kM,                      // mass
88     kCharge,                 // charge
89     kNclsITS,                // number of clusters assigned in the ITS
90     kITSchi2Cl,              // chi2/cl in the ITS
91     kNclsTPC,                // number of clusters assigned in the TPC
92     kNclsSTPC,                // number of shared clusters assigned in the TPC
93     kNclsSFracTPC,           // fraction of shared clusters assigned in the TPC
94     kNclsTPCiter1,           // number of clusters assigned in the TPC after first iteration
95     kNFclsTPC,               // number of findable clusters in the TPC
96     kNFclsTPCr,              // number of findable clusters in the TPC with more robust definition
97     kNFclsTPCrFrac,          // number of found/findable clusters in the TPC with more robust definition
98     kTPCsignalN,             // number of points used for dEdx
99     kTPCsignalNfrac,         // fraction of points used for dEdx / cluster used for tracking
100     kTPCchi2Cl,              // chi2/cl in TPC
101     kTPCclsDiff,             // TPC cluster difference
102     kTrackStatus,            // track status bits
103     
104     kNclsTRD,                // number of clusters assigned in the TRD
105     kTRDntracklets,          // number of TRD tracklets used for tracking/PID TODO: correct getter
106     kTRDpidQuality,          // number of TRD tracklets used for PID
107     kTRDprobEle,             // TRD electron pid probability
108     kTRDprobPio,             // TRD electron pid probability
109     kTRDphi,                 // Phi angle of the track at the entrance of the TRD
110     kTRDpidEffLeg,           // TRD pid efficiency from conversion electrons
111       
112     kImpactParXY,            // Impact parameter in XY plane
113     kImpactParZ,             // Impact parameter in Z
114     kTrackLength,            // Track length
115
116     kPdgCode,                // PDG code
117     kPdgCodeMother,          // PDG code of the mother
118     kPdgCodeGrandMother,     // PDG code of the grandmother
119     kNumberOfDaughters,      // number of daughters
120     kHaveSameMother,         // check that particles have the same mother (MC)
121     kIsJpsiPrimary,          // check if the particle is primary (MC)
122     kNumberOfJPsis,          // number of generated inclusive jpsis per event (MC)
123     kNumberOfJPsisPrompt,    // number of generated prompt jpsis per event (MC)
124     kNumberOfJPsisNPrompt,   // number of generated non-prompt jpsis per event (MC)
125
126     kITSsignal,              // ITS dE/dx signal
127     kITSsignalSSD1,          // SSD1 dE/dx signal
128     kITSsignalSSD2,          // SSD2 dE/dx signal
129     kITSsignalSDD1,          // SDD1 dE/dx signal
130     kITSsignalSDD2,          // SDD2 dE/dx signal
131     kITSclusterMap,          // ITS cluster map
132     kITSLayerFirstCls,       // No of innermost ITS layer with a cluster of a track
133     kITSnSigmaEle,           // number of sigmas to the dE/dx electron line in the ITS
134     kITSnSigmaPio,           // number of sigmas to the dE/dx pion line in the ITS
135     kITSnSigmaMuo,           // number of sigmas to the dE/dx muon line in the ITS
136     kITSnSigmaKao,           // number of sigmas to the dE/dx kaon line in the ITS
137     kITSnSigmaPro,           // number of sigmas to the dE/dx proton line in the ITS
138
139     kPIn,                    // momentum at inner wall of TPC (if available), used for PID
140     kPOut,                   // momentum at outer wall of TPC, used for TRD studies
141     kYsignedIn,              // signed local y at inner wall of TPC
142     kTPCsignal,              // TPC dE/dx signal
143     
144     kTOFsignal,              // TOF signal
145     kTOFbeta,                // TOF beta
146     kTOFPIDBit,              // TOF PID bit (1:set, 0:TOF not available)
147         
148     kTPCnSigmaEle,           // number of sigmas to the dE/dx electron line in the TPC
149     kTPCnSigmaPio,           // number of sigmas to the dE/dx pion line in the TPC
150     kTPCnSigmaMuo,           // number of sigmas to the dE/dx muon line in the TPC
151     kTPCnSigmaKao,           // number of sigmas to the dE/dx kaon line in the TPC
152     kTPCnSigmaPro,           // number of sigmas to the dE/dx proton line in the TPC
153       
154     kTOFnSigmaEle,           // number of sigmas to the electron line in the TOF
155     kTOFnSigmaPio,           // number of sigmas to the pion line in the TOF
156     kTOFnSigmaMuo,           // number of sigmas to the muon line in the TOF
157     kTOFnSigmaKao,           // number of sigmas to the kaon line in the TOF
158     kTOFnSigmaPro,           // number of sigmas to the proton line in the TOF
159
160     kEMCALnSigmaEle,         // number of sigmas to the proton line in the TOF
161     kEMCALEoverP,            // E over P from EMCAL
162     kEMCALNCells,            // NCells from EMCAL
163     kEMCALM02,               // M02 showershape parameter
164     kEMCALM20,               // M20 showershape parameter
165     kEMCALDispersion,        // Dispersion paramter
166     
167     kKinkIndex0,             // kink index 0
168       
169     kParticleMax,             //
170     // TODO: kRNClusters ??
171     // AliDielectronPair specific variables
172     kChi2NDF = kParticleMax, // Chi^2/NDF
173     kDecayLength,            // decay length
174     kR,                      // distance to the origin
175     kOpeningAngle,           // opening angle
176     // helicity picture: Z-axis is considered the direction of the mother's 3-momentum vector
177     kThetaHE,                // theta in mother's rest frame in the helicity picture 
178     kPhiHE,                  // phi in mother's rest frame in the helicity picture
179     kCos2PhiHE,              // Cosine of 2*phi in mother's rest frame in the helicity picture
180     kCosTilPhiHE,            // Sin(phi) +/- Cos(phi) in mother's rest frame in the helicity picture; Sign is sign(Cos(phi))
181     // Collins-Soper picture: Z-axis is considered the direction of the vectorial difference between 
182     // the 3-mom vectors of target and projectile beams
183     kThetaCS,                // theta in mother's rest frame in Collins-Soper picture
184     kPhiCS,                  // phi in mother's rest frame in Collins-Soper picture
185     kCos2PhiCS,              // Cosine of 2*phi in mother's rest frame in the Collins-Soper picture
186     kCosTilPhiCS,            // Sin(phi) +/- Cos(phi) in mother's rest frame in Collins-Soper picture; Sign is sign(Cos(phi))
187     kDeltaPhiV0ArpH2,        // Delta phi of the pair with respect to the 2nd order harmonic reaction plane from V0-A
188     kDeltaPhiV0CrpH2,        // Delta phi of the pair with respect to the 2nd order harmonic reaction plane from V0-C
189     kDeltaPhiV0ACrpH2,       // Delta phi of the pair with respect to the 2nd order harmonic reaction plane from V0-A + V0-C
190     kV0ArpH2FlowV2,          // v2 coefficient with respect to the 2nd order reaction plane from V0-A
191     kV0CrpH2FlowV2,          // v2 coefficient with respect to the 2nd order reaction plane from V0-C
192     kV0ACrpH2FlowV2,         // v2 coefficient with respect to the 2nd order reaction plane from V0-A + V0-C
193     kLegDist,                // distance of the legs
194     kLegDistXY,              // distance of the legs in XY
195     kDeltaEta,         // Absolute value of Delta Eta for the legs
196     kDeltaPhi,           // Absolute value of Delta Phi for the legs
197     kMerr,                   // error of mass calculation
198     kDCA,                    // distance of closest approach TODO: not implemented yet
199     kPairType,               // type of the pair, like like sign ++ unlikesign ...
200     kPseudoProperTime,       // pseudo proper time
201     kPseudoProperTimeResolution,     // resolution for pseudo proper decay time (reconstructed - MC truth)
202     kPseudoProperTimePull,   // normalizd resolution for pseudo proper time = (reco - MC truth)/dReco
203     kTRDpidEffPair,          // TRD pid efficieny from conversion electrons
204     kPairMax,                 //
205   // Event specific variables
206     kXvPrim=kPairMax,        // prim vertex
207     kYvPrim,                 // prim vertex
208     kZvPrim,                 // prim vertex
209     kXRes,                   // primary vertex x-resolution
210     kYRes,                   // primary vertex y-resolution
211     kZRes,                   // primary vertex z-resolution
212
213     // v0 reaction plane quantities from AliEPSelectionTaks
214     kv0ArpH2,                  // VZERO-A reaction plane of the Q vector for 2nd harmonic
215     kv0CrpH2,                  //         reaction plane
216     kv0ACrpH2,                 // VZERO-AC reaction plane of the Q vector for 2nd harmonic
217     kv0ATPCDiffH2,             // V0A-TPC reaction plane difference for 2nd harmonic
218     kv0CTPCDiffH2,             // V0C-TPC reaction plane difference for 2nd harmonic
219     kv0Av0CDiffH2,             // V0A-V0C reaction plane difference for 2nd harmonic
220
221     kMultV0A,                // VZERO multiplicity and ADC amplitudes
222     kMultV0C,
223     kMultV0,
224     kAdcV0A,
225     kAdcV0C,
226     kAdcV0,
227     kVZEROchMult,
228     // VZERO reaction plane quantities
229     kV0AxH2=kVZEROchMult+64,   // VZERO-A x-component of the Q vector for 2nd harmonic
230     kV0AyH2,                   // VZERO-A y-component of the Q vector for 2nd harmonic
231     kV0ArpH2,                  // VZERO-A reaction plane of the Q vector for 2nd harmonic
232     kV0CxH2,                   // VZERO-C x-component of the Q vector for 2nd harmonic
233     kV0CyH2,                   //         y-component
234     kV0CrpH2,                  //         reaction plane
235     kV0ACxH2,                  // VZERO-AC x-component of the Q vector for 2nd harmonic
236     kV0ACyH2,                  // VZERO-AC y-component of the Q vector for 2nd harmonic
237     kV0ACrpH2,                 // VZERO-AC reaction plane of the Q vector for 2nd harmonic
238     kV0ArpResH2,               // 2nd harmonic reaction plane resolution for V0A
239     kV0CrpResH2,               //                               V0C
240     kV0ACrpResH2,              //                             V0A+V0C
241     kV0XaXcH2,                 // Correlation quantities to check V0 reaction plane quality
242     kV0XaYaH2,
243     kV0XaYcH2,
244     kV0YaXcH2,
245     kV0YaYcH2,
246     kV0XcYcH2,
247     kV0ATPCDiffH2,             // V0A-TPC reaction plane difference for 2nd harmonic
248     kV0CTPCDiffH2,             // V0C-TPC reaction plane difference for 2nd harmonic
249     kV0AV0CDiffH2,             // V0A-V0C reaction plane difference for 2nd harmonic
250     // TPC reaction plane quantities
251     kTPCxH2,                  // TPC x-component of the Q vector for 2nd harmonic (corrected)
252     kTPCyH2,                  // TPC y-component of the Q vector for 2nd harmonic (corrected)
253     kTPCrpH2,                 // TPC reaction plane of the Q vector for 2nd harmonic (corrected)
254     kTPCsub1xH2,              // TPC x-component of the Q vector for 2nd harmonic (corrected, sub event 1) 
255     kTPCsub1yH2,              // TPC y-component of the Q vector for 2nd harmonic (corrected, sub event 1)
256     kTPCsub1rpH2,             // TPC reaction plane of the Q vector for 2nd harmonic (corrected, sub event 1)
257     kTPCsub2xH2,              // TPC x-component of the Q vector for 2nd harmonic (corrected, sub event 2)
258     kTPCsub2yH2,              // TPC y-component of the Q vector for 2nd harmonic (corrected, sub event 2)
259     kTPCsub2rpH2,             // TPC reaction plane of the Q vector for 2nd harmonic (corrected, sub event 2)
260     kTPCsub12DiffH2,          // TPC reaction plane difference of sub event 1,2 for 2nd harmonic
261     kTPCsub12DiffH2Sin,       // TPC reaction plane difference of sub event 1,2 for 2nd harmonic, sinus term
262
263     kTPCxH2uc,                  // TPC x-component of the Q vector for 2nd harmonic (uncorrected)
264     kTPCyH2uc,                  // TPC y-component of the Q vector for 2nd harmonic (uncorrected)
265     kTPCrpH2uc,                 // TPC reaction plane of the Q vector for 2nd harmonic (uncorrected)
266     kTPCsub1xH2uc,              // TPC x-component of the Q vector for 2nd harmonic (uncorrected, sub event 1) 
267     kTPCsub1yH2uc,              // TPC y-component of the Q vector for 2nd harmonic (uncorrected, sub event 1)
268     kTPCsub1rpH2uc,             // TPC reaction plane of the Q vector for 2nd harmonic (uncorrected, sub event 1)
269     kTPCsub2xH2uc,              // TPC x-component of the Q vector for 2nd harmonic (uncorrected, sub event 2)
270     kTPCsub2yH2uc,              // TPC y-component of the Q vector for 2nd harmonic (uncorrected, sub event 2)
271     kTPCsub2rpH2uc,             // TPC reaction plane of the Q vector for 2nd harmonic (uncorrected, sub event 2)
272     kTPCsub12DiffH2uc,          // TPC reaction plane difference of sub event 1,2 for 2nd harmonic (uncorrected)
273
274     kNTrk,                   // number of tracks (or tracklets) TODO: ambiguous
275     kTracks,                 // ESD tracks TODO: ambiguous
276     kNVtxContrib,             // number of primary vertex contibutors
277     kNacc,                   // Number of accepted tracks
278     kNaccTrcklts,            // number of accepted SPD tracklets in |eta|<1.6        
279     kNaccTrcklts0916,        // number of accepted SPD tracklets in 0.9<|eta|<1.6
280     
281     kNaccTrckltsEsd05,       // number of accepted SPD tracklets in |eta|<0.5 (AliESDEvent::EstimateMultiplicity())
282     kNaccTrckltsEsd10,       // number of accepted SPD tracklets in |eta|<1.0 (AliESDEvent::EstimateMultiplicity())
283     kNaccTrckltsEsd16,       // number of accepted SPD tracklets in |eta|<1.6 (AliESDEvent::EstimateMultiplicity())
284     kNaccTrckltsEsd05Corr,   //
285     kNaccTrckltsEsd10Corr,   //
286     kNaccTrckltsEsd16Corr,   //
287     kNaccItsTpcEsd05,        // ITS-TPC tracks + ITS SA complementary tracks + tracklets from unassigned tracklets in |eta|<0.5 (AliESDEvent::EstimateMultiplicity())
288     kNaccItsTpcEsd10,        // ITS-TPC tracks + ITS SA complementary tracks + tracklets from unassigned tracklets in |eta|<1.0 (AliESDEvent::EstimateMultiplicity())
289     kNaccItsTpcEsd16,        // ITS-TPC tracks + ITS SA complementary tracks + tracklets from unassigned tracklets in |eta|<1.6 (AliESDEvent::EstimateMultiplicity())
290     kNaccItsTpcEsd05Corr,        // 
291     kNaccItsTpcEsd10Corr,        // 
292     kNaccItsTpcEsd16Corr,        // 
293     
294     kNaccItsPureEsd05,       // ITS SA tracks + tracklets from unassigned tracklets in |eta|<0.5 (AliESDEvent::EstimateMultiplicity())
295     kNaccItsPureEsd10,       // ITS SA tracks + tracklets from unassigned tracklets in |eta|<1.0 (AliESDEvent::EstimateMultiplicity())
296     kNaccItsPureEsd16,       // ITS SA tracks + tracklets from unassigned tracklets in |eta|<1.6 (AliESDEvent::EstimateMultiplicity())
297     kNaccItsPureEsd05Corr,   // 
298     kNaccItsPureEsd10Corr,   // 
299     kNaccItsPureEsd16Corr,   // 
300     
301     kNch,                    // MC true number of charged particles in |eta|<1.6
302     kNch05,                  // MC true number of charged particles in |eta|<0.5
303     kNch10,                  // MC true number of charged particles in |eta|<1.0
304
305     kCentrality,             // event centrality fraction
306     kCentralitySPD,          // centrality using SPD
307     kNevents,                // event counter
308     kRunNumber,              // run number
309     kMixingBin,
310     kNMaxValues              //
311     // TODO: (for A+A) ZDCEnergy, impact parameter, Iflag??
312   };
313   
314
315   AliDielectronVarManager();
316   AliDielectronVarManager(const char* name, const char* title);
317   virtual ~AliDielectronVarManager();
318   static void Fill(const TObject* particle, Double_t * const values);
319   static void FillVarMCParticle2(const AliVParticle *p1, const AliVParticle *p2, Double_t * const values);
320
321   static void InitESDpid(Int_t type=0);
322   static void InitAODpidUtil(Int_t type=0);
323   static void InitEstimatorAvg(const Char_t* filename);
324   static void InitTRDpidEffHistograms(const Char_t* filename);
325   static void SetVZEROCalibrationFile(const Char_t* filename) {fgVZEROCalibrationFile = filename;}
326   
327   static void SetVZERORecenteringFile(const Char_t* filename) {fgVZERORecenteringFile = filename;}
328   static void SetPIDResponse(AliPIDResponse *pidResponse) {fgPIDResponse=pidResponse;}
329   static AliPIDResponse* GetPIDResponse() { return fgPIDResponse; }
330   static void SetEvent(AliVEvent * const ev);
331   static void SetEventData(const Double_t data[AliDielectronVarManager::kNMaxValues]);
332   static Bool_t GetDCA(const AliAODTrack *track, Double_t d0z0[2]);
333   static void SetTPCEventPlane(AliEventplane *const evplane);
334   static void GetVzeroRP(const AliVEvent* event, Double_t* qvec, Int_t sideOption);      // 0- V0A; 1- V0C; 2- V0A+V0C
335
336   static TProfile* GetEstimatorHistogram(Int_t period, Int_t type) {return fgMultEstimatorAvg[period][type];}
337   static Double_t GetTRDpidEfficiency(Int_t runNo, Double_t centrality, Double_t eta, Double_t trdPhi, Double_t pout, Double_t& effErr);
338
339   static const AliKFVertex* GetKFVertex() {return fgKFVertex;}
340   
341   static const char* GetValueName(Int_t i) { return (i>=0&&i<kNMaxValues)?fgkParticleNames[i]:""; }
342
343   static const Double_t* GetData() {return fgData;}
344
345   static Double_t GetValue(ValueTypes var) {return fgData[var];}
346   static void SetValue(ValueTypes var, Double_t val) { fgData[var]=val; }
347   
348 private:
349
350   static const char* fgkParticleNames[kNMaxValues];  //variable names
351
352   static void FillVarVParticle(const AliVParticle *particle,         Double_t * const values);
353   static void FillVarESDtrack(const AliESDtrack *particle,           Double_t * const values);
354   static void FillVarAODTrack(const AliAODTrack *particle,           Double_t * const values);
355   static void FillVarMCParticle(const AliMCParticle *particle,       Double_t * const values);
356   static void FillVarAODMCParticle(const AliAODMCParticle *particle, Double_t * const values);
357   static void FillVarDielectronPair(const AliDielectronPair *pair,   Double_t * const values);
358   static void FillVarKFParticle(const AliKFParticle *pair,   Double_t * const values);
359   
360   static void FillVarVEvent(const AliVEvent *event,                  Double_t * const values);
361   static void FillVarESDEvent(const AliESDEvent *event,              Double_t * const values);
362   static void FillVarAODEvent(const AliAODEvent *event,              Double_t * const values);
363   static void FillVarMCEvent(const AliMCEvent *event,                Double_t * const values);
364   static void FillVarTPCEventPlane(const AliEventplane *evplane,     Double_t * const values);
365
366   static void InitVZEROCalibrationHistograms(Int_t runNo);
367   static void InitVZERORecenteringHistograms(Int_t runNo);
368   
369   static AliPIDResponse  *fgPIDResponse;        // PID response object
370   static AliVEvent       *fgEvent;              // current event pointer
371   static AliEventplane   *fgTPCEventPlane;      // current event plane pointer
372   static AliKFVertex     *fgKFVertex;           // kf vertex
373   static TProfile        *fgMultEstimatorAvg[4][9];  // multiplicity estimator averages (4 periods x 9 estimators)
374   static Double_t         fgTRDpidEffCentRanges[10][4];   // centrality ranges for the TRD pid efficiency histograms
375   static TH3D            *fgTRDpidEff[10][4];   // TRD pid efficiencies from conversion electrons
376   static TString          fgVZEROCalibrationFile;  // file with VZERO channel-by-channel calibrations
377   static TString          fgVZERORecenteringFile;  // file with VZERO Q-vector averages needed for event plane recentering
378   static TProfile2D      *fgVZEROCalib[64];           // 1 histogram per VZERO channel
379   static TProfile2D      *fgVZERORecentering[2][2];   // 2 VZERO sides x 2 Q-vector components
380   static Int_t            fgCurrentRun;
381   
382   static Double_t fgData[kNMaxValues];        //! data
383   
384   AliDielectronVarManager(const AliDielectronVarManager &c);
385   AliDielectronVarManager &operator=(const AliDielectronVarManager &c);
386   
387   ClassDef(AliDielectronVarManager,1);
388 };
389
390
391 //Inline functions
392 inline void AliDielectronVarManager::Fill(const TObject* object, Double_t * const values)
393 {
394   //
395   // Main function to fill all available variables according to the type of particle
396   //
397 //   if (!object) return;
398   if      (object->IsA() == AliESDtrack::Class())       FillVarESDtrack(static_cast<const AliESDtrack*>(object), values);
399   else if (object->IsA() == AliAODTrack::Class())       FillVarAODTrack(static_cast<const AliAODTrack*>(object), values);
400   else if (object->IsA() == AliMCParticle::Class())     FillVarMCParticle(static_cast<const AliMCParticle*>(object), values);
401   else if (object->IsA() == AliAODMCParticle::Class())  FillVarAODMCParticle(static_cast<const AliAODMCParticle*>(object), values);
402   else if (object->IsA() == AliDielectronPair::Class()) FillVarDielectronPair(static_cast<const AliDielectronPair*>(object), values);
403   else if (object->IsA() == AliKFParticle::Class())     FillVarKFParticle(static_cast<const AliKFParticle*>(object),values);
404   // Main function to fill all available variables according to the type of event
405   
406   else if (object->IsA() == AliVEvent::Class())         FillVarVEvent(static_cast<const AliVEvent*>(object), values);
407   else if (object->IsA() == AliESDEvent::Class())       FillVarESDEvent(static_cast<const AliESDEvent*>(object), values);
408   else if (object->IsA() == AliAODEvent::Class())       FillVarAODEvent(static_cast<const AliAODEvent*>(object), values);
409   else if (object->IsA() == AliMCEvent::Class())        FillVarMCEvent(static_cast<const AliMCEvent*>(object), values);
410   else if (object->IsA() == AliEventplane::Class())     FillVarTPCEventPlane(static_cast<const AliEventplane*>(object), values);
411 //   else printf(Form("AliDielectronVarManager::Fill: Type %s is not supported by AliDielectronVarManager!", object->ClassName())); //TODO: implement without object needed
412 }
413
414 inline void AliDielectronVarManager::FillVarVParticle(const AliVParticle *particle, Double_t * const values)
415 {
416   //
417   // Fill track information available in AliVParticle into an array
418   //
419   values[AliDielectronVarManager::kPx]        = particle->Px();
420   values[AliDielectronVarManager::kPy]        = particle->Py();
421   values[AliDielectronVarManager::kPz]        = particle->Pz();
422   values[AliDielectronVarManager::kPt]        = particle->Pt();
423   values[AliDielectronVarManager::kP]         = particle->P();
424
425   values[AliDielectronVarManager::kXv]        = particle->Xv();
426   values[AliDielectronVarManager::kYv]        = particle->Yv();
427   values[AliDielectronVarManager::kZv]        = particle->Zv();
428
429   values[AliDielectronVarManager::kOneOverPt] = (particle->Pt()>1.0e-3 ? particle->OneOverPt() : 0.0);
430   values[AliDielectronVarManager::kPhi]       = particle->Phi();
431   values[AliDielectronVarManager::kTheta]     = particle->Theta();
432   values[AliDielectronVarManager::kEta]       = particle->Eta();
433   values[AliDielectronVarManager::kY]         = particle->Y();
434   
435   values[AliDielectronVarManager::kE]         = particle->E();
436   values[AliDielectronVarManager::kM]         = particle->M();
437   values[AliDielectronVarManager::kCharge]    = particle->Charge();
438   
439   values[AliDielectronVarManager::kPdgCode]   = particle->PdgCode();
440     
441 //   if ( fgEvent ) AliDielectronVarManager::Fill(fgEvent, values);
442   for (Int_t i=AliDielectronVarManager::kPairMax; i<AliDielectronVarManager::kNMaxValues; ++i)
443     values[i]=fgData[i];
444 }
445
446 inline void AliDielectronVarManager::FillVarESDtrack(const AliESDtrack *particle, Double_t * const values)
447 {
448   //
449   // Fill track information available for histogramming into an array
450   //
451
452   // Fill common AliVParticle interface information
453   FillVarVParticle(particle, values);
454
455   AliESDtrack *esdTrack=0x0;
456   Double_t origdEdx=particle->GetTPCsignal();
457   
458   // apply ETa correction, remove once this is in the tender
459   esdTrack=const_cast<AliESDtrack*>(particle);
460   if (!esdTrack) return;
461   esdTrack->SetTPCsignal(origdEdx/AliDielectronPID::GetEtaCorr(esdTrack)/AliDielectronPID::GetCorrValdEdx(),esdTrack->GetTPCsignalSigma(),esdTrack->GetTPCsignalN());
462   
463   Double_t pidProbs[AliPID::kSPECIES];
464   // Fill AliESDtrack interface specific information
465   Double_t tpcNcls=particle->GetTPCNcls();
466   Double_t tpcNclsS = particle->GetTPCnclsS(); 
467   Double_t itsNcls=particle->GetNcls(0);
468   Double_t tpcSignalN=particle->GetTPCsignalN();
469   values[AliDielectronVarManager::kNclsITS]       = itsNcls; // TODO: get rid of the plain numbers
470   values[AliDielectronVarManager::kNclsTPC]       = tpcNcls; // TODO: get rid of the plain numbers
471   values[AliDielectronVarManager::kNclsSTPC]      = tpcNclsS;
472   values[AliDielectronVarManager::kNclsSFracTPC]  = tpcNcls>0?tpcNclsS/tpcNcls:0;
473   values[AliDielectronVarManager::kNclsTPCiter1]  = particle->GetTPCNclsIter1(); // TODO: get rid of the plain numbers
474   values[AliDielectronVarManager::kNFclsTPC]      = particle->GetTPCNclsF();
475   values[AliDielectronVarManager::kNFclsTPCr]     = particle->GetTPCClusterInfo(2,1);
476   values[AliDielectronVarManager::kNFclsTPCrFrac] = particle->GetTPCClusterInfo(2);
477   values[AliDielectronVarManager::kTPCsignalN]    = tpcSignalN;
478   values[AliDielectronVarManager::kTPCsignalNfrac]= tpcNcls>0?tpcSignalN/tpcNcls:0;
479   values[AliDielectronVarManager::kNclsTRD]       = particle->GetNcls(2); // TODO: get rid of the plain numbers
480   values[AliDielectronVarManager::kTRDntracklets] = particle->GetTRDntracklets(); // TODO: GetTRDtracklets/GetTRDntracklets?
481   values[AliDielectronVarManager::kTRDpidQuality] = particle->GetTRDpidQuality();
482   values[AliDielectronVarManager::kTPCclsDiff]    = tpcSignalN-tpcNcls;
483   values[AliDielectronVarManager::kTrackStatus]   = (Double_t)particle->GetStatus();
484   
485   
486   values[AliDielectronVarManager::kTPCchi2Cl] = -1;
487   if (tpcNcls>0) values[AliDielectronVarManager::kTPCchi2Cl] = particle->GetTPCchi2() / tpcNcls;
488   values[AliDielectronVarManager::kITSchi2Cl] = -1;
489   if (itsNcls>0) values[AliDielectronVarManager::kITSchi2Cl] = particle->GetITSchi2() / itsNcls;
490   //TRD pidProbs
491   particle->GetTRDpid(pidProbs);
492   values[AliDielectronVarManager::kTRDprobEle]    = pidProbs[AliPID::kElectron];
493   values[AliDielectronVarManager::kTRDprobPio]    = pidProbs[AliPID::kPion];
494
495   values[AliDielectronVarManager::kKinkIndex0]    = particle->GetKinkIndex(0);
496   
497   Float_t impactParXY, impactParZ;
498   particle->GetImpactParameters(impactParXY, impactParZ);
499   values[AliDielectronVarManager::kImpactParXY]   = impactParXY;
500   values[AliDielectronVarManager::kImpactParZ]    = impactParZ;
501
502
503   values[AliDielectronVarManager::kPdgCode]=-1;
504   values[AliDielectronVarManager::kPdgCodeMother]=-1;
505   values[AliDielectronVarManager::kPdgCodeGrandMother]=-1;
506   
507   values[AliDielectronVarManager::kNumberOfDaughters]=-999;
508   
509   AliDielectronMC *mc=AliDielectronMC::Instance();
510   
511   if (mc->HasMC()){
512     if (mc->GetMCTrack(particle))
513       values[AliDielectronVarManager::kPdgCode]=mc->GetMCTrack(particle)->PdgCode();
514     
515     AliMCParticle *motherMC=mc->GetMCTrackMother(particle); //mother
516     if (motherMC){
517       values[AliDielectronVarManager::kPdgCodeMother]=motherMC->PdgCode();
518       
519       motherMC=mc->GetMCTrackMother(motherMC);  //grand motherMC
520       if (motherMC) values[AliDielectronVarManager::kPdgCodeGrandMother]=motherMC->PdgCode();;
521     }
522       
523     values[AliDielectronVarManager::kNumberOfDaughters]=mc->NumberOfDaughters(particle);
524   } //if(mc->HasMC())
525   
526
527
528   values[AliDielectronVarManager::kITSsignal]   =   particle->GetITSsignal();
529   
530   Double_t itsdEdx[4];
531   particle->GetITSdEdxSamples(itsdEdx);
532
533   values[AliDielectronVarManager::kITSsignalSSD1]   =   itsdEdx[0];
534   values[AliDielectronVarManager::kITSsignalSSD2]   =   itsdEdx[1];
535   values[AliDielectronVarManager::kITSsignalSDD1]   =   itsdEdx[2];
536   values[AliDielectronVarManager::kITSsignalSDD2]   =   itsdEdx[3];
537   values[AliDielectronVarManager::kITSclusterMap]   =   particle->GetITSClusterMap();
538   values[AliDielectronVarManager::kITSLayerFirstCls] = -1.;
539
540   for (Int_t iC=0; iC<6; iC++) {
541     if (((particle->GetITSClusterMap()) & (1<<(iC))) > 0) {
542       values[AliDielectronVarManager::kITSLayerFirstCls] = iC;
543       break;
544     }
545   }
546
547   
548   values[AliDielectronVarManager::kTrackLength]   = particle->GetIntegratedLength();
549   //dEdx information
550   Double_t mom = particle->GetP();
551   const AliExternalTrackParam *in=particle->GetInnerParam();
552   Double_t ysignedIn=-100;
553   if (in) {
554     mom = in->GetP();
555     ysignedIn=particle->Charge()*in->GetY();
556   }
557   values[AliDielectronVarManager::kPIn]=mom;
558   values[AliDielectronVarManager::kYsignedIn]=ysignedIn;
559   const AliExternalTrackParam *out=particle->GetOuterParam();
560   if(out) values[AliDielectronVarManager::kPOut] = out->GetP();
561   else values[AliDielectronVarManager::kPOut] = mom;
562   if(out && fgEvent) {
563     Double_t localCoord[3]={0.0};
564     Bool_t localCoordGood = out->GetXYZAt(298.0, ((AliESDEvent*)fgEvent)->GetMagneticField(), localCoord);
565     values[AliDielectronVarManager::kTRDphi] = (localCoordGood && TMath::Abs(localCoord[0])>1.0e-6 && TMath::Abs(localCoord[1])>1.0e-6 ? TMath::ATan2(localCoord[1], localCoord[0]) : -999.);
566   }
567   if(mc->HasMC() && fgTRDpidEff[0][0]) {
568     Int_t runNo = (fgEvent ? fgEvent->GetRunNumber() : -1);
569     Float_t centrality=-1.0;
570     AliCentrality *esdCentrality = (fgEvent ? fgEvent->GetCentrality() : 0x0);
571     if(esdCentrality) centrality = esdCentrality->GetCentralityPercentile("V0M");
572     Double_t effErr=0.0;
573     values[kTRDpidEffLeg] = GetTRDpidEfficiency(runNo, centrality, values[AliDielectronVarManager::kEta], 
574                                                 values[AliDielectronVarManager::kTRDphi], 
575                                                 values[AliDielectronVarManager::kPOut], effErr);
576   }
577   values[AliDielectronVarManager::kTPCsignal]=particle->GetTPCsignal();
578   
579   values[AliDielectronVarManager::kTOFsignal]=particle->GetTOFsignal();
580   
581   Double_t l = particle->GetIntegratedLength();  // cm
582   Double_t t = particle->GetTOFsignal();
583   Double_t t0 = fgPIDResponse->GetTOFResponse().GetTimeZero(); // ps
584
585   if( (l < 360. || l > 800.) || (t <= 0.) || (t0 >999990.0) ) {
586         values[AliDielectronVarManager::kTOFbeta]=0.0;
587   }
588   else {
589         t -= t0; // subtract the T0
590         l *= 0.01;  // cm ->m
591         t *= 1e-12; //ps -> s
592     
593         Double_t v = l / t;
594         Float_t beta = v / TMath::C();
595         values[AliDielectronVarManager::kTOFbeta]=beta;
596   }
597   values[AliDielectronVarManager::kTOFPIDBit]=(particle->GetStatus()&AliESDtrack::kTOFpid? 1: 0);
598   
599   // nsigma to Electron band
600   // TODO: for the moment we set the bethe bloch parameters manually
601   //       this should be changed in future!
602   
603   values[AliDielectronVarManager::kTPCnSigmaEle]=fgPIDResponse->NumberOfSigmasTPC(particle,AliPID::kElectron)-AliDielectronPID::GetCorrVal();
604   values[AliDielectronVarManager::kTPCnSigmaPio]=fgPIDResponse->NumberOfSigmasTPC(particle,AliPID::kPion);
605   values[AliDielectronVarManager::kTPCnSigmaMuo]=fgPIDResponse->NumberOfSigmasTPC(particle,AliPID::kMuon);
606   values[AliDielectronVarManager::kTPCnSigmaKao]=fgPIDResponse->NumberOfSigmasTPC(particle,AliPID::kKaon);
607   values[AliDielectronVarManager::kTPCnSigmaPro]=fgPIDResponse->NumberOfSigmasTPC(particle,AliPID::kProton);
608
609   values[AliDielectronVarManager::kITSnSigmaEle]=fgPIDResponse->NumberOfSigmasITS(particle,AliPID::kElectron);
610   values[AliDielectronVarManager::kITSnSigmaPio]=fgPIDResponse->NumberOfSigmasITS(particle,AliPID::kPion);
611   values[AliDielectronVarManager::kITSnSigmaMuo]=fgPIDResponse->NumberOfSigmasITS(particle,AliPID::kMuon);
612   values[AliDielectronVarManager::kITSnSigmaKao]=fgPIDResponse->NumberOfSigmasITS(particle,AliPID::kKaon);
613   values[AliDielectronVarManager::kITSnSigmaPro]=fgPIDResponse->NumberOfSigmasITS(particle,AliPID::kProton);
614
615   values[AliDielectronVarManager::kTOFnSigmaEle]=fgPIDResponse->NumberOfSigmasTOF(particle,AliPID::kElectron);
616   values[AliDielectronVarManager::kTOFnSigmaPio]=fgPIDResponse->NumberOfSigmasTOF(particle,AliPID::kPion);
617   values[AliDielectronVarManager::kTOFnSigmaMuo]=fgPIDResponse->NumberOfSigmasTOF(particle,AliPID::kMuon);
618   values[AliDielectronVarManager::kTOFnSigmaKao]=fgPIDResponse->NumberOfSigmasTOF(particle,AliPID::kKaon);
619   values[AliDielectronVarManager::kTOFnSigmaPro]=fgPIDResponse->NumberOfSigmasTOF(particle,AliPID::kProton);
620
621   //EMCAL PID information
622   Double_t eop=0;
623   Double_t showershape[4]={0.,0.,0.,0.};
624 //   values[AliDielectronVarManager::kEMCALnSigmaEle]  = fgPIDResponse->NumberOfSigmasEMCAL(particle,AliPID::kElectron);
625   values[AliDielectronVarManager::kEMCALnSigmaEle]  = fgPIDResponse->NumberOfSigmasEMCAL(particle,AliPID::kElectron,eop,showershape);
626   values[AliDielectronVarManager::kEMCALEoverP]     = eop;
627   values[AliDielectronVarManager::kEMCALNCells]     = showershape[0];
628   values[AliDielectronVarManager::kEMCALM02]        = showershape[1];
629   values[AliDielectronVarManager::kEMCALM20]        = showershape[2];
630   values[AliDielectronVarManager::kEMCALDispersion] = showershape[3];
631   
632   
633   //restore TPC signal if it was changed
634   if (esdTrack) esdTrack->SetTPCsignal(origdEdx,esdTrack->GetTPCsignalSigma(),esdTrack->GetTPCsignalN());
635 }
636
637 inline void AliDielectronVarManager::FillVarAODTrack(const AliAODTrack *particle, Double_t * const values)
638 {
639   //
640   // Fill track information available for histogramming into an array
641   //
642
643   // Fill common AliVParticle interface information
644   FillVarVParticle(particle, values);
645   Double_t tpcNcls=particle->GetTPCNcls();
646
647   //GetNclsS not present in AODtrack
648   //Replace with method as soon as available
649   TBits tpcSharedMap = particle->GetTPCSharedMap();   
650   Double_t tpcNclsS=  tpcSharedMap.CountBits(0)-tpcSharedMap.CountBits(159);
651   // Reset AliESDtrack interface specific information
652   values[AliDielectronVarManager::kNclsITS]       = 0;
653   values[AliDielectronVarManager::kITSchi2Cl]     = -1;
654   values[AliDielectronVarManager::kNclsTPC]       = tpcNcls;
655   values[AliDielectronVarManager::kNclsSTPC]       = tpcNclsS;
656   values[AliDielectronVarManager::kNclsSFracTPC]       = tpcNcls>0?tpcNclsS/tpcNcls:0;
657   values[AliDielectronVarManager::kNclsTPCiter1]  = tpcNcls; // not really available in AOD
658   values[AliDielectronVarManager::kNFclsTPC]      = 0;
659   values[AliDielectronVarManager::kNFclsTPCr]     = 0;
660   values[AliDielectronVarManager::kNFclsTPCrFrac] = 0;
661   values[AliDielectronVarManager::kNclsTRD]       = 0;
662   values[AliDielectronVarManager::kTRDntracklets] = 0;
663   values[AliDielectronVarManager::kTRDpidQuality] = 0;
664   
665   values[AliDielectronVarManager::kTPCchi2Cl] = -1;
666   values[AliDielectronVarManager::kTrackStatus]   = (Double_t)particle->GetStatus();
667   
668   //TRD pidProbs
669   //TODO: set correctly
670   values[AliDielectronVarManager::kTRDprobEle]    = 0;
671   values[AliDielectronVarManager::kTRDprobPio]    = 0;
672   
673   values[AliDielectronVarManager::kTPCsignalN]    = 0;
674   values[AliDielectronVarManager::kTPCsignalNfrac]= 0;
675   
676   // Fill AliAODTrack interface information
677   //
678
679   values[AliDielectronVarManager::kKinkIndex0]    = particle->GetProdVertex()->GetType()==AliAODVertex::kKink ? 1 : 0;
680
681   Double_t d0z0[2];
682   GetDCA(particle, d0z0);
683   values[AliDielectronVarManager::kImpactParXY]   = d0z0[0];
684   values[AliDielectronVarManager::kImpactParZ]    = d0z0[1];
685
686   values[AliDielectronVarManager::kPIn]=0;
687   values[AliDielectronVarManager::kTPCsignal]=0;
688   
689   values[AliDielectronVarManager::kTOFsignal]=0;
690   values[AliDielectronVarManager::kTOFbeta]=0;
691
692   values[AliDielectronVarManager::kTPCnSigmaEle]=0;
693   values[AliDielectronVarManager::kTPCnSigmaPio]=0;
694   values[AliDielectronVarManager::kTPCnSigmaMuo]=0;
695   values[AliDielectronVarManager::kTPCnSigmaKao]=0;
696   values[AliDielectronVarManager::kTPCnSigmaPro]=0;
697
698   values[AliDielectronVarManager::kITSclusterMap]   =   particle->GetITSClusterMap();
699   
700   AliAODPid *pid=particle->GetDetPid();
701   if (pid){
702     Double_t mom =pid->GetTPCmomentum();
703     Double_t tpcSignalN=pid->GetTPCsignalN();
704     values[AliDielectronVarManager::kTPCsignalN] = tpcSignalN;
705     values[AliDielectronVarManager::kTPCsignalN] = tpcNcls>0?tpcSignalN/tpcNcls:0;
706     values[AliDielectronVarManager::kTPCclsDiff] = tpcSignalN-tpcNcls;
707     Double_t tpcNsigmaEle=fgPIDResponse->NumberOfSigmasTPC(particle,AliPID::kElectron);
708     Double_t tpcNsigmaPio=fgPIDResponse->NumberOfSigmasTPC(particle,AliPID::kPion);
709     Double_t tpcNsigmaMuo=fgPIDResponse->NumberOfSigmasTPC(particle,AliPID::kMuon);
710     Double_t tpcNsigmaKao=fgPIDResponse->NumberOfSigmasTPC(particle,AliPID::kKaon);
711     Double_t tpcNsigmaPro=fgPIDResponse->NumberOfSigmasTPC(particle,AliPID::kProton);
712     
713     values[AliDielectronVarManager::kPIn]=mom;
714     values[AliDielectronVarManager::kTPCsignal]=pid->GetTPCsignal();
715
716     values[AliDielectronVarManager::kTPCnSigmaEle]=tpcNsigmaEle;
717     values[AliDielectronVarManager::kTPCnSigmaPio]=tpcNsigmaPio;
718     values[AliDielectronVarManager::kTPCnSigmaMuo]=tpcNsigmaMuo;
719     values[AliDielectronVarManager::kTPCnSigmaKao]=tpcNsigmaKao;
720     values[AliDielectronVarManager::kTPCnSigmaPro]=tpcNsigmaPro;
721     
722     values[AliDielectronVarManager::kTRDntracklets] = 0;
723     values[AliDielectronVarManager::kTRDpidQuality] = 0;
724
725     Double_t tofNsigmaEle=fgPIDResponse->NumberOfSigmasTOF(particle,AliPID::kElectron);
726     Double_t tofNsigmaPio=fgPIDResponse->NumberOfSigmasTOF(particle,AliPID::kPion);
727     Double_t tofNsigmaMuo=fgPIDResponse->NumberOfSigmasTOF(particle,AliPID::kMuon);
728     Double_t tofNsigmaKao=fgPIDResponse->NumberOfSigmasTOF(particle,AliPID::kKaon);
729     Double_t tofNsigmaPro=fgPIDResponse->NumberOfSigmasTOF(particle,AliPID::kProton);
730     
731     values[AliDielectronVarManager::kTOFnSigmaEle]=tofNsigmaEle;
732     values[AliDielectronVarManager::kTOFnSigmaPio]=tofNsigmaPio;
733     values[AliDielectronVarManager::kTOFnSigmaMuo]=tofNsigmaMuo;
734     values[AliDielectronVarManager::kTOFnSigmaKao]=tofNsigmaKao;
735     values[AliDielectronVarManager::kTOFnSigmaPro]=tofNsigmaPro;
736     
737   }
738
739   //EMCAL PID information
740   Double_t eop=0;
741   Double_t showershape[4]={0.,0.,0.,0.};
742 //   values[AliDielectronVarManager::kEMCALnSigmaEle]  = fgPIDResponse->NumberOfSigmasEMCAL(particle,AliPID::kElectron);
743   values[AliDielectronVarManager::kEMCALnSigmaEle]  = fgPIDResponse->NumberOfSigmasEMCAL(particle,AliPID::kElectron,eop,showershape);
744   values[AliDielectronVarManager::kEMCALEoverP]     = eop;
745   values[AliDielectronVarManager::kEMCALNCells]     = showershape[0];
746   values[AliDielectronVarManager::kEMCALM02]        = showershape[1];
747   values[AliDielectronVarManager::kEMCALM20]        = showershape[2];
748   values[AliDielectronVarManager::kEMCALDispersion] = showershape[3];
749
750   values[AliDielectronVarManager::kPdgCode]=-1;
751   values[AliDielectronVarManager::kPdgCodeMother]=-1;
752   values[AliDielectronVarManager::kPdgCodeGrandMother]=-1;
753   
754   values[AliDielectronVarManager::kNumberOfDaughters]=-999;
755   
756   AliDielectronMC *mc=AliDielectronMC::Instance();
757   
758   if (mc->HasMC()){
759     if (mc->GetMCTrack(particle))
760       values[AliDielectronVarManager::kPdgCode]=mc->GetMCTrack(particle)->PdgCode();
761     
762     AliAODMCParticle *motherMC=mc->GetMCTrackMother(particle); //mother
763     if (motherMC){
764       values[AliDielectronVarManager::kPdgCodeMother]=motherMC->PdgCode();
765       
766       motherMC=mc->GetMCTrackMother(motherMC);  //grand motherMC
767       if (motherMC) values[AliDielectronVarManager::kPdgCodeGrandMother]=motherMC->PdgCode();;
768     }
769     
770     values[AliDielectronVarManager::kNumberOfDaughters]=mc->NumberOfDaughters(particle);
771   } //if(mc->HasMC())
772   
773   values[AliDielectronVarManager::kTOFPIDBit]=(particle->GetStatus()&AliESDtrack::kTOFpid? 1: 0);
774 }
775
776 inline void AliDielectronVarManager::FillVarMCParticle(const AliMCParticle *particle, Double_t * const values)
777 {
778   //
779   // Fill track information available for histogramming into an array
780   //
781
782   values[AliDielectronVarManager::kNclsITS]       = 0;
783   values[AliDielectronVarManager::kITSchi2Cl]     = 0;
784   values[AliDielectronVarManager::kNclsTPC]       = 0;
785   values[AliDielectronVarManager::kNclsSTPC]       = 0;
786   values[AliDielectronVarManager::kNclsSFracTPC]       = 0;
787   values[AliDielectronVarManager::kNclsTPCiter1]  = 0; 
788   values[AliDielectronVarManager::kNFclsTPC]      = 0;
789   values[AliDielectronVarManager::kNFclsTPCr]     = 0;
790   values[AliDielectronVarManager::kNFclsTPCrFrac] = 0;
791   values[AliDielectronVarManager::kNclsTRD]       = 0;
792   values[AliDielectronVarManager::kTRDntracklets] = 0;
793   values[AliDielectronVarManager::kTRDpidQuality] = 0;
794   values[AliDielectronVarManager::kTPCchi2Cl]     = 0;
795   values[AliDielectronVarManager::kTrackStatus]   = 0;
796   values[AliDielectronVarManager::kTRDprobEle]    = 0;
797   values[AliDielectronVarManager::kTRDprobPio]    = 0;
798   values[AliDielectronVarManager::kTPCsignalN]    = 0;
799   values[AliDielectronVarManager::kTPCclsDiff]    = 0;
800   values[AliDielectronVarManager::kTPCsignalNfrac]    = 0;
801   values[AliDielectronVarManager::kImpactParXY]   = 0;
802   values[AliDielectronVarManager::kImpactParZ]    = 0;
803   values[AliDielectronVarManager::kPIn]           = 0;
804   values[AliDielectronVarManager::kYsignedIn]     = 0;
805   values[AliDielectronVarManager::kTPCsignal]     = 0;
806   values[AliDielectronVarManager::kTOFsignal]     = 0;
807   values[AliDielectronVarManager::kTOFbeta]       = 0;
808   values[AliDielectronVarManager::kTPCnSigmaEle]  = 0;
809   values[AliDielectronVarManager::kTPCnSigmaPio]  = 0;
810   values[AliDielectronVarManager::kTPCnSigmaMuo]  = 0;
811   values[AliDielectronVarManager::kTPCnSigmaKao]  = 0;
812   values[AliDielectronVarManager::kTPCnSigmaPro]  = 0;
813   values[AliDielectronVarManager::kITSclusterMap] = 0;
814   
815   values[AliDielectronVarManager::kPdgCode]       = -1;
816   values[AliDielectronVarManager::kPdgCodeMother] = -1;
817   values[AliDielectronVarManager::kPdgCodeGrandMother] = -1;
818   
819   // Fill common AliVParticle interface information
820   FillVarVParticle(particle, values);
821   
822   AliDielectronMC *mc=AliDielectronMC::Instance();
823
824   // Fill AliMCParticle interface specific information
825   values[AliDielectronVarManager::kPdgCode] = particle->PdgCode();
826
827   AliMCParticle *motherMC = mc->GetMCTrackMother(particle);
828   if (motherMC){
829     values[AliDielectronVarManager::kPdgCodeMother]=motherMC->PdgCode();
830     motherMC=mc->GetMCTrackMother(motherMC);  //grand mother
831     if (motherMC) values[AliDielectronVarManager::kPdgCodeGrandMother]=motherMC->PdgCode();;
832   }
833   
834   values[AliDielectronVarManager::kIsJpsiPrimary] = mc->IsJpsiPrimary(particle);
835   values[AliDielectronVarManager::kNumberOfDaughters]=mc->NumberOfDaughters(particle);
836 }
837
838
839 inline void AliDielectronVarManager::FillVarMCParticle2(const AliVParticle *p1, const AliVParticle *p2, Double_t * const values) {
840   //
841   // fill 2 track information starting from MC legs
842   //
843
844   values[AliDielectronVarManager::kNclsITS]       = 0;
845   values[AliDielectronVarManager::kITSchi2Cl]     = -1;
846   values[AliDielectronVarManager::kNclsTPC]       = 0;
847   values[AliDielectronVarManager::kNclsSTPC]       = 0;
848   values[AliDielectronVarManager::kNclsSFracTPC]  = 0;
849   values[AliDielectronVarManager::kNclsTPCiter1]  = 0; 
850   values[AliDielectronVarManager::kNFclsTPC]      = 0;
851   values[AliDielectronVarManager::kNFclsTPCr]     = 0;
852   values[AliDielectronVarManager::kNFclsTPCrFrac] = 0;
853   values[AliDielectronVarManager::kNclsTRD]       = 0;
854   values[AliDielectronVarManager::kTRDntracklets] = 0;
855   values[AliDielectronVarManager::kTRDpidQuality] = 0;
856   values[AliDielectronVarManager::kTPCchi2Cl]     = 0;
857   values[AliDielectronVarManager::kTrackStatus]   = 0;
858   values[AliDielectronVarManager::kTRDprobEle]    = 0;
859   values[AliDielectronVarManager::kTRDprobPio]    = 0;
860   values[AliDielectronVarManager::kTPCsignalN]    = 0;
861   values[AliDielectronVarManager::kTPCclsDiff]    = 0;
862   values[AliDielectronVarManager::kTPCsignalNfrac]    = 0;
863   values[AliDielectronVarManager::kImpactParXY]   = 0;
864   values[AliDielectronVarManager::kImpactParZ]    = 0;
865   values[AliDielectronVarManager::kPIn]           = 0;
866   values[AliDielectronVarManager::kYsignedIn]     = 0;
867   values[AliDielectronVarManager::kTPCsignal]     = 0;
868   values[AliDielectronVarManager::kTPCnSigmaEle]  = 0;
869   values[AliDielectronVarManager::kTPCnSigmaPio]  = 0;
870   values[AliDielectronVarManager::kTPCnSigmaMuo]  = 0;
871   values[AliDielectronVarManager::kTPCnSigmaKao]  = 0;
872   values[AliDielectronVarManager::kTPCnSigmaPro]  = 0;
873   values[AliDielectronVarManager::kITSclusterMap] = 0;
874   
875   values[AliDielectronVarManager::kPdgCode]       = 0;
876   values[AliDielectronVarManager::kPdgCodeMother] = 0;
877
878   AliDielectronMC *mc=AliDielectronMC::Instance();
879   AliVParticle* mother=0x0;
880   Int_t mLabel1 = mc->GetMothersLabel(p1->GetLabel());
881   Int_t mLabel2 = mc->GetMothersLabel(p2->GetLabel());
882   if(mLabel1==mLabel2)
883     mother = mc->GetMCTrackFromMCEvent(mLabel1);
884
885   values[AliDielectronVarManager::kPseudoProperTime] = -2e10;
886   if(mother) {    // same mother
887     FillVarVParticle(mother, values);
888     Double_t lxy = ((mother->Xv()- values[AliDielectronVarManager::kXvPrim]) * mother->Px() + 
889                     (mother->Yv()- values[AliDielectronVarManager::kYvPrim]) * mother->Py() )/mother->Pt();
890     values[AliDielectronVarManager::kPseudoProperTime] = lxy*(TDatabasePDG::Instance()->GetParticle(443)->Mass())/mother->Pt();
891   }
892   else {         // the 2 particles come from different mothers so 2-particles quantities are calculated here
893     // AliVParticle part
894     values[AliDielectronVarManager::kPx]        = p1->Px()+p2->Px();
895     values[AliDielectronVarManager::kPy]        = p1->Py()+p2->Py();
896     values[AliDielectronVarManager::kPz]        = p1->Pz()+p2->Pz();
897     values[AliDielectronVarManager::kPt]        = TMath::Sqrt(values[AliDielectronVarManager::kPx]*
898                                                               values[AliDielectronVarManager::kPx]+
899                                                               values[AliDielectronVarManager::kPy]*
900                                                               values[AliDielectronVarManager::kPy]);
901     values[AliDielectronVarManager::kP]         = TMath::Sqrt(values[AliDielectronVarManager::kPt]*
902                                                               values[AliDielectronVarManager::kPt]+
903                                                               values[AliDielectronVarManager::kPz]*
904                                                               values[AliDielectronVarManager::kPz]);
905     
906     values[AliDielectronVarManager::kXv]        = 0;
907     values[AliDielectronVarManager::kYv]        = 0;
908     values[AliDielectronVarManager::kZv]        = 0;
909     
910     values[AliDielectronVarManager::kOneOverPt] = (values[AliDielectronVarManager::kPt]>1.0e-6 ? 1.0/values[AliDielectronVarManager::kPt] : 0.0);
911     values[AliDielectronVarManager::kPhi]       = TMath::ATan2(values[AliDielectronVarManager::kPy],values[AliDielectronVarManager::kPx]);
912     values[AliDielectronVarManager::kTheta]     = TMath::ATan2(values[AliDielectronVarManager::kPt],values[AliDielectronVarManager::kPz]);
913     values[AliDielectronVarManager::kEta]       = ((values[AliDielectronVarManager::kP]-values[AliDielectronVarManager::kPz])>1.0e-6 && (values[AliDielectronVarManager::kP]+values[AliDielectronVarManager::kPz])>1.0e-6 ? 0.5*TMath::Log((values[AliDielectronVarManager::kP]+values[AliDielectronVarManager::kPz])/(values[AliDielectronVarManager::kP]-values[AliDielectronVarManager::kPz])) : -9999.);
914     values[AliDielectronVarManager::kE]         = p1->E()+p2->E();
915     values[AliDielectronVarManager::kY]         = ((values[AliDielectronVarManager::kE]-values[AliDielectronVarManager::kPz])>1.0e-6 && (values[AliDielectronVarManager::kE]+values[AliDielectronVarManager::kPz])>1.0e-6 ? 0.5*TMath::Log((values[AliDielectronVarManager::kE]+values[AliDielectronVarManager::kPz])/(values[AliDielectronVarManager::kE]-values[AliDielectronVarManager::kPz])) : -9999.);
916     values[AliDielectronVarManager::kCharge]    = p1->Charge()+p2->Charge();
917
918     values[AliDielectronVarManager::kM]         = TMath::Sqrt(p1->M()*p1->M()+p2->M()*p2->M()+
919       2.0*(p1->E()*p2->E()-p1->Px()*p2->Px()-p1->Py()*p2->Py()-p1->Pz()*p2->Pz()));
920
921     if ( fgEvent ) AliDielectronVarManager::Fill(fgEvent, values);
922   }
923
924   values[AliDielectronVarManager::kThetaHE]   = AliDielectronPair::ThetaPhiCM(p1,p2,kTRUE,  kTRUE);
925   values[AliDielectronVarManager::kPhiHE]     = AliDielectronPair::ThetaPhiCM(p1,p2,kTRUE,  kFALSE);
926   values[AliDielectronVarManager::kCos2PhiHE] = TMath::Cos(2*values[AliDielectronVarManager::kPhiHE]);
927   values[AliDielectronVarManager::kThetaCS]   = AliDielectronPair::ThetaPhiCM(p1,p2,kFALSE, kTRUE);
928   values[AliDielectronVarManager::kPhiCS]     = AliDielectronPair::ThetaPhiCM(p1,p2,kFALSE, kFALSE);
929   values[AliDielectronVarManager::kCos2PhiCS] = TMath::Cos(2*values[AliDielectronVarManager::kPhiCS]);
930   values[AliDielectronVarManager::kCosTilPhiHE]  = (TMath::Cos(values[AliDielectronVarManager::kPhiHE])>0)?(TMath::Cos(values[AliDielectronVarManager::kPhiHE]-TMath::Pi()/4.)):(TMath::Cos(values[AliDielectronVarManager::kPhiHE]-3*TMath::Pi()/4.));
931   values[AliDielectronVarManager::kCosTilPhiCS]  = (TMath::Cos(values[AliDielectronVarManager::kPhiCS])>0)?(TMath::Cos(values[AliDielectronVarManager::kPhiCS]-TMath::Pi()/4.)):(TMath::Cos(values[AliDielectronVarManager::kPhiCS]-3*TMath::Pi()/4.));
932 }
933
934
935 inline void AliDielectronVarManager::FillVarAODMCParticle(const AliAODMCParticle *particle, Double_t * const values)
936 {
937   //
938   // Fill track information available for histogramming into an array
939   //
940
941   values[AliDielectronVarManager::kNclsITS]       = 0;
942   values[AliDielectronVarManager::kITSchi2Cl]     = -1;
943   values[AliDielectronVarManager::kNclsTPC]       = 0;
944   values[AliDielectronVarManager::kNclsSTPC]      = 0;
945   values[AliDielectronVarManager::kNclsSFracTPC]  = 0;
946   values[AliDielectronVarManager::kNclsTPCiter1]  = 0;
947   values[AliDielectronVarManager::kNFclsTPC]      = 0;
948   values[AliDielectronVarManager::kNclsTRD]       = 0;
949   values[AliDielectronVarManager::kTRDntracklets] = 0;
950   values[AliDielectronVarManager::kTRDpidQuality] = 0;
951   values[AliDielectronVarManager::kTPCchi2Cl]     = 0;
952   values[AliDielectronVarManager::kTrackStatus]   = 0;
953   values[AliDielectronVarManager::kTRDprobEle]    = 0;
954   values[AliDielectronVarManager::kTRDprobPio]    = 0;
955   values[AliDielectronVarManager::kTPCsignalN]    = 0;
956   values[AliDielectronVarManager::kTPCclsDiff]    = 0;
957   values[AliDielectronVarManager::kTPCsignalNfrac]= 0;
958   values[AliDielectronVarManager::kImpactParXY]   = 0;
959   values[AliDielectronVarManager::kImpactParZ]    = 0;
960   values[AliDielectronVarManager::kPIn]           = 0;
961   values[AliDielectronVarManager::kYsignedIn]     = 0;
962   values[AliDielectronVarManager::kTPCsignal]     = 0;
963   values[AliDielectronVarManager::kTPCnSigmaEle]  = 0;
964   values[AliDielectronVarManager::kTPCnSigmaPio]  = 0;
965   values[AliDielectronVarManager::kTPCnSigmaMuo]  = 0;
966   values[AliDielectronVarManager::kTPCnSigmaKao]  = 0;
967   values[AliDielectronVarManager::kTPCnSigmaPro]  = 0;
968   values[AliDielectronVarManager::kITSclusterMap] = 0;
969   
970   values[AliDielectronVarManager::kPdgCode]       = -1;
971   values[AliDielectronVarManager::kPdgCodeMother] = -1;
972   values[AliDielectronVarManager::kPdgCodeGrandMother] = -1;
973   
974   // Fill common AliVParticle interface information
975   FillVarVParticle(particle, values);
976   
977   AliDielectronMC *mc=AliDielectronMC::Instance();
978
979
980   // Fill AliAODMCParticle interface specific information
981   values[AliDielectronVarManager::kPdgCode] = particle->PdgCode();
982
983   AliAODMCParticle *motherMC = mc->GetMCTrackMother(particle);
984   if (motherMC){
985     values[AliDielectronVarManager::kPdgCodeMother]=motherMC->PdgCode();
986     motherMC=mc->GetMCTrackMother(motherMC);  //grand mother
987     if (motherMC) values[AliDielectronVarManager::kPdgCodeGrandMother]=motherMC->PdgCode();;
988   }
989   
990   
991   values[AliDielectronVarManager::kIsJpsiPrimary] = mc->IsJpsiPrimary(particle);
992
993   values[AliDielectronVarManager::kNumberOfDaughters]=mc->NumberOfDaughters(particle);
994 }
995
996 inline void AliDielectronVarManager::FillVarDielectronPair(const AliDielectronPair *pair, Double_t * const values)
997 {
998   //
999   // Fill pair information available for histogramming into an array
1000   //
1001   
1002   values[AliDielectronVarManager::kPdgCode]=-1;
1003   values[AliDielectronVarManager::kPdgCodeMother]=-1;
1004   values[AliDielectronVarManager::kPdgCodeGrandMother]=-1;
1005   
1006   // Fill common AliVParticle interface information
1007   FillVarVParticle(pair, values);
1008
1009   // Fill AliDielectronPair specific information
1010   const AliKFParticle &kfPair = pair->GetKFParticle();
1011
1012   Double_t thetaHE=0;
1013   Double_t phiHE=0;
1014   Double_t thetaCS=0;
1015   Double_t phiCS=0;
1016
1017   pair->GetThetaPhiCM(thetaHE,phiHE,thetaCS,phiCS);
1018   
1019   values[AliDielectronVarManager::kChi2NDF]      = kfPair.GetChi2()/kfPair.GetNDF();
1020   values[AliDielectronVarManager::kDecayLength]  = kfPair.GetDecayLength();
1021   values[AliDielectronVarManager::kR]            = kfPair.GetR();
1022   values[AliDielectronVarManager::kOpeningAngle] = pair->OpeningAngle();
1023   values[AliDielectronVarManager::kThetaHE]      = thetaHE;
1024   values[AliDielectronVarManager::kPhiHE]        = phiHE;
1025   values[AliDielectronVarManager::kCos2PhiHE]    = TMath::Cos(2.0*phiHE);
1026   values[AliDielectronVarManager::kCosTilPhiHE]  = (TMath::Cos(phiHE)>0)?(TMath::Cos(phiHE-TMath::Pi()/4.)):(TMath::Cos(phiHE-3*TMath::Pi()/4.));
1027   values[AliDielectronVarManager::kThetaCS]      = thetaCS;
1028   values[AliDielectronVarManager::kPhiCS]        = phiCS;
1029   values[AliDielectronVarManager::kCos2PhiCS]    = TMath::Cos(2.0*phiCS);
1030   values[AliDielectronVarManager::kCosTilPhiCS]  = (TMath::Cos(phiCS)>0)?(TMath::Cos(phiCS-TMath::Pi()/4.)):(TMath::Cos(phiCS-3*TMath::Pi()/4.));
1031   values[AliDielectronVarManager::kLegDist]      = pair->DistanceDaughters();
1032   values[AliDielectronVarManager::kLegDistXY]    = pair->DistanceDaughtersXY();
1033   values[AliDielectronVarManager::kDeltaEta]     = pair->DeltaEta();
1034   values[AliDielectronVarManager::kDeltaPhi]     = pair->DeltaPhi();
1035   values[AliDielectronVarManager::kMerr]         = kfPair.GetErrMass()>1e-30&&kfPair.GetMass()>1e-30?kfPair.GetErrMass()/kfPair.GetMass():1000000;
1036   values[AliDielectronVarManager::kPairType]     = pair->GetType();
1037
1038   Double_t errPseudoProperTime2 = -1;
1039   values[AliDielectronVarManager::kPseudoProperTime] = fgEvent ? pair->GetKFParticle().GetPseudoProperDecayTime(*(fgEvent->GetPrimaryVertex()), TDatabasePDG::Instance()->GetParticle(443)->Mass(), &errPseudoProperTime2 ) : -1e10;
1040   // values[AliDielectronVarManager::kPseudoProperTime] = fgEvent ? pair->GetPseudoProperTime(fgEvent->GetPrimaryVertex()): -1e10;
1041
1042   // Flow quantities
1043   Double_t delta=0.0;
1044   // v2 with respect to VZERO-A event plane
1045   delta = values[AliDielectronVarManager::kPhi] - fgData[AliDielectronVarManager::kV0ArpH2];
1046   if(delta>TMath::Pi()) delta -= 2.0*TMath::Pi();             // keep the [-pi,+pi] interval
1047   if(delta<-1.0*TMath::Pi()) delta += 2.0*TMath::Pi();
1048   values[AliDielectronVarManager::kV0ArpH2FlowV2] = TMath::Cos(2.0*delta);  // 2nd harmonic flow coefficient
1049   values[AliDielectronVarManager::kDeltaPhiV0ArpH2] = delta;
1050   // v2 with respect to VZERO-C event plane
1051   delta = values[AliDielectronVarManager::kPhi] - fgData[AliDielectronVarManager::kV0CrpH2];
1052   if(delta>TMath::Pi()) delta -= 2.0*TMath::Pi();             // keep the [-pi,+pi] interval
1053   if(delta<-1.0*TMath::Pi()) delta += 2.0*TMath::Pi();
1054   values[AliDielectronVarManager::kV0CrpH2FlowV2] = TMath::Cos(2.0*delta);  // 2nd harmonic flow coefficient
1055   values[AliDielectronVarManager::kDeltaPhiV0CrpH2] = delta;
1056   // v2 with respect to the combined VZERO-A and VZERO-C event plane
1057   delta = values[AliDielectronVarManager::kPhi] - fgData[AliDielectronVarManager::kV0ACrpH2];
1058   if(delta>TMath::Pi()) delta -= 2.0*TMath::Pi();             // keep the [-pi,+pi] interval
1059   if(delta<-1.0*TMath::Pi()) delta += 2.0*TMath::Pi();
1060   values[AliDielectronVarManager::kV0ACrpH2FlowV2] = TMath::Cos(2.0*delta);  // 2nd harmonic flow coefficient
1061   values[AliDielectronVarManager::kDeltaPhiV0ACrpH2] = delta;
1062
1063              
1064   AliDielectronMC *mc=AliDielectronMC::Instance();
1065   
1066   if (mc->HasMC()){
1067     values[AliDielectronVarManager::kPseudoProperTimeResolution] = -10.0e+10;
1068     Bool_t samemother =  mc->HaveSameMother(pair);
1069     values[AliDielectronVarManager::kIsJpsiPrimary] = mc->IsJpsiPrimary(pair);
1070     values[AliDielectronVarManager::kHaveSameMother] = samemother ;
1071
1072       // fill kPseudoProperTimeResolution
1073     values[AliDielectronVarManager::kPseudoProperTimeResolution] = -1e10;
1074     // values[AliDielectronVarManager::kPseudoProperTimePull] = -1e10;
1075     if(samemother && fgEvent) {
1076       if(pair->GetFirstDaughter()->GetLabel() > 0) {
1077         const AliVParticle* d1 = mc->GetMCTrackFromMCEvent(pair->GetFirstDaughter()->GetLabel());
1078         const AliVParticle* motherMC = mc->GetMCTrackFromMCEvent(((AliMCParticle*)d1)->GetMother());
1079         const AliVVertex* mcVtx = mc->GetMCEvent()->GetPrimaryVertex();
1080         if(motherMC && mcVtx) {
1081           const Double_t lxyMC = ( (motherMC->Xv() - mcVtx->GetX()) * motherMC->Px() +
1082                                    (motherMC->Yv() - mcVtx->GetY()) * motherMC->Py()   ) / motherMC->Pt();
1083           const Double_t pseudoMC = lxyMC * (TDatabasePDG::Instance()->GetParticle(443)->Mass())/motherMC->Pt();
1084           values[AliDielectronVarManager::kPseudoProperTimeResolution] = values[AliDielectronVarManager::kPseudoProperTime] - pseudoMC;
1085           if (errPseudoProperTime2 > 0)
1086             values[AliDielectronVarManager::kPseudoProperTimePull] = values[AliDielectronVarManager::kPseudoProperTimeResolution]/sqrt(errPseudoProperTime2);
1087         }
1088       }
1089     }
1090     Double_t valuesLeg1[AliDielectronVarManager::kNMaxValues];
1091     Double_t valuesLeg2[AliDielectronVarManager::kNMaxValues];
1092     AliVParticle* leg1 = pair->GetFirstDaughter();
1093     AliVParticle* leg2 = pair->GetSecondDaughter();
1094     Fill(leg1, valuesLeg1);
1095     Fill(leg2, valuesLeg2);
1096     values[AliDielectronVarManager::kTRDpidEffPair] = valuesLeg1[AliDielectronVarManager::kTRDpidEffLeg]*valuesLeg2[AliDielectronVarManager::kTRDpidEffLeg];
1097   }//if (mc->HasMC())
1098
1099
1100 }
1101
1102 inline void AliDielectronVarManager::FillVarKFParticle(const AliKFParticle *particle, Double_t * const values)
1103 {
1104   //
1105   // Fill track information available in AliVParticle into an array
1106   //
1107   values[AliDielectronVarManager::kPx]        = particle->GetPx();
1108   values[AliDielectronVarManager::kPy]        = particle->GetPy();
1109   values[AliDielectronVarManager::kPz]        = particle->GetPz();
1110   values[AliDielectronVarManager::kPt]        = particle->GetPt();
1111   values[AliDielectronVarManager::kP]         = particle->GetP();
1112   
1113   values[AliDielectronVarManager::kXv]        = particle->GetX();
1114   values[AliDielectronVarManager::kYv]        = particle->GetY();
1115   values[AliDielectronVarManager::kZv]        = particle->GetZ();
1116   
1117   values[AliDielectronVarManager::kOneOverPt] = 0;
1118   values[AliDielectronVarManager::kPhi]       = particle->GetPhi();
1119   values[AliDielectronVarManager::kTheta]     = 0.;
1120   values[AliDielectronVarManager::kEta]       = particle->GetEta();
1121   values[AliDielectronVarManager::kY]         = ((particle->GetE()*particle->GetE()-particle->GetPx()*particle->GetPx()-particle->GetPy()*particle->GetPy()-particle->GetPz()*particle->GetPz())>0.) ? TLorentzVector(particle->GetPx(),particle->GetPy(),particle->GetPz(),particle->GetE()).Rapidity() : -1111.;
1122   
1123   values[AliDielectronVarManager::kE]         = particle->GetE();
1124   values[AliDielectronVarManager::kM]         = particle->GetMass();
1125   values[AliDielectronVarManager::kCharge]    = particle->GetQ();
1126   
1127   values[AliDielectronVarManager::kNclsITS]       = 0;
1128   values[AliDielectronVarManager::kITSchi2Cl]     = -1;
1129   values[AliDielectronVarManager::kNclsTPC]       = 0;
1130   values[AliDielectronVarManager::kNclsSTPC]      = 0;
1131   values[AliDielectronVarManager::kNclsSFracTPC]  = 0;
1132   values[AliDielectronVarManager::kNclsTPCiter1]  = 0;
1133   values[AliDielectronVarManager::kNFclsTPC]      = 0;
1134   values[AliDielectronVarManager::kNclsTRD]       = 0;
1135   values[AliDielectronVarManager::kTRDntracklets] = 0;
1136   values[AliDielectronVarManager::kTRDpidQuality] = 0;
1137   values[AliDielectronVarManager::kTPCchi2Cl]     = 0;
1138   values[AliDielectronVarManager::kTrackStatus]   = 0;
1139   values[AliDielectronVarManager::kTRDprobEle]    = 0;
1140   values[AliDielectronVarManager::kTRDprobPio]    = 0;
1141   values[AliDielectronVarManager::kTPCsignalN]    = 0;
1142   values[AliDielectronVarManager::kTPCclsDiff]    = 0;
1143   values[AliDielectronVarManager::kTPCsignalNfrac]= 0;
1144   values[AliDielectronVarManager::kImpactParXY]   = 0;
1145   values[AliDielectronVarManager::kImpactParZ]    = 0;
1146   values[AliDielectronVarManager::kPIn]           = 0;
1147   values[AliDielectronVarManager::kYsignedIn]     = 0;
1148   values[AliDielectronVarManager::kTPCsignal]     = 0;
1149   values[AliDielectronVarManager::kTOFsignal]     = 0;
1150   values[AliDielectronVarManager::kTOFbeta]       = 0;
1151   values[AliDielectronVarManager::kTPCnSigmaEle]  = 0;
1152   values[AliDielectronVarManager::kTPCnSigmaPio]  = 0;
1153   values[AliDielectronVarManager::kTPCnSigmaMuo]  = 0;
1154   values[AliDielectronVarManager::kTPCnSigmaKao]  = 0;
1155   values[AliDielectronVarManager::kTPCnSigmaPro]  = 0;
1156   values[AliDielectronVarManager::kITSclusterMap] = 0;
1157   
1158   values[AliDielectronVarManager::kPdgCode]       = -1;
1159   values[AliDielectronVarManager::kPdgCodeMother] = -1;
1160   values[AliDielectronVarManager::kPdgCodeGrandMother] = -1;
1161   
1162   
1163 //   if ( fgEvent ) AliDielectronVarManager::Fill(fgEvent, values);
1164   for (Int_t i=AliDielectronVarManager::kPairMax; i<AliDielectronVarManager::kNMaxValues; ++i)
1165     values[i]=fgData[i];
1166   
1167 }
1168
1169 inline void AliDielectronVarManager::FillVarVEvent(const AliVEvent *event, Double_t * const values)
1170 {
1171   //
1172   // Fill event information available for histogramming into an array
1173   //
1174   values[AliDielectronVarManager::kRunNumber]    = event->GetRunNumber();
1175   if(fgCurrentRun!=event->GetRunNumber()) {
1176     if(fgVZEROCalibrationFile.Contains(".root")) InitVZEROCalibrationHistograms(event->GetRunNumber());
1177     if(fgVZERORecenteringFile.Contains(".root")) InitVZERORecenteringHistograms(event->GetRunNumber());
1178     fgCurrentRun=event->GetRunNumber();
1179   }
1180   const AliVVertex *primVtx = event->GetPrimaryVertex();
1181   
1182   values[AliDielectronVarManager::kXvPrim]       = 0;
1183   values[AliDielectronVarManager::kYvPrim]       = 0;
1184   values[AliDielectronVarManager::kZvPrim]       = 0;
1185 //   values[AliDielectronVarManager::kChi2NDF]      = 0; //This is the pair value!!!
1186
1187   values[AliDielectronVarManager::kNTrk]         = 0;
1188   values[AliDielectronVarManager::kNVtxContrib]   = 0;
1189   values[AliDielectronVarManager::kNacc]         = 0;
1190   values[AliDielectronVarManager::kNaccTrcklts]  = 0;
1191   values[AliDielectronVarManager::kNaccTrcklts0916] = 0;
1192   values[AliDielectronVarManager::kNevents]      = 0; //always fill bin 0;
1193   
1194   if (!primVtx) return;
1195   
1196   values[AliDielectronVarManager::kXvPrim]       = primVtx->GetX();
1197   values[AliDielectronVarManager::kYvPrim]       = primVtx->GetY();
1198   values[AliDielectronVarManager::kZvPrim]       = primVtx->GetZ();
1199 //   values[AliDielectronVarManager::kChi2NDF]      = primVtx->GetChi2perNDF(); //this is the pair value
1200
1201   values[AliDielectronVarManager::kNTrk]            = event->GetNumberOfTracks();
1202   values[AliDielectronVarManager::kNVtxContrib]     = primVtx->GetNContributors();
1203   values[AliDielectronVarManager::kNacc]            = AliDielectronHelper::GetNacc(event);
1204   values[AliDielectronVarManager::kNaccTrcklts]     = AliDielectronHelper::GetNaccTrcklts(event);      // etaRange = 1.6 (default)
1205   values[AliDielectronVarManager::kNaccTrcklts0916] = AliDielectronHelper::GetNaccTrcklts(event,1.6)-AliDielectronHelper::GetNaccTrcklts(event,.9);
1206   //  values[AliDielectronVarManager::kNaccTrcklts05]   = AliDielectronHelper::GetNaccTrcklts(event, 0.5);
1207   //  values[AliDielectronVarManager::kNaccTrcklts10]   = AliDielectronHelper::GetNaccTrcklts(event, 1.0);
1208   //  values[AliDielectronVarManager::kNaccTrckltsCorr] = AliDielectronHelper::GetNaccTrckltsCorrected(event, values[AliDielectronVarManager::kNaccTrcklts], values[AliDielectronVarManager::kZvPrim]);
1209
1210   // TPC event plane (corrected)
1211
1212   // VZERO event plane quantities from the AliEPSelectionTask
1213   AliEventplane *ep = const_cast<AliVEvent*>(event)->GetEventplane();
1214   values[AliDielectronVarManager::kv0ACrpH2]  = TVector2::Phi_mpi_pi(ep->GetEventplane("V0", event, 2));
1215   values[AliDielectronVarManager::kv0ArpH2]   = TVector2::Phi_mpi_pi(ep->GetEventplane("V0A",event, 2));
1216   values[AliDielectronVarManager::kv0CrpH2]   = TVector2::Phi_mpi_pi(ep->GetEventplane("V0C",event, 2));
1217
1218   // ESD VZERO information
1219   AliVVZERO* vzeroData = event->GetVZEROData();
1220   values[AliDielectronVarManager::kMultV0A] = 0.0;
1221   values[AliDielectronVarManager::kMultV0C] = 0.0;
1222   values[AliDielectronVarManager::kAdcV0A]  = 0.0;
1223   values[AliDielectronVarManager::kAdcV0C]  = 0.0;
1224   for(Int_t i=0; i<32; ++i) {
1225     values[AliDielectronVarManager::kVZEROchMult+i] = vzeroData->GetMultiplicity(i);
1226     values[AliDielectronVarManager::kVZEROchMult+32+i] = vzeroData->GetMultiplicity(i+32);
1227     //values[AliDielectronVarManager::kVZEROchMult+i] = event->GetVZEROEqMultiplicity(i);
1228     //values[AliDielectronVarManager::kVZEROchMult+32+i] = event->GetVZEROEqMultiplicity(i+32);
1229     values[AliDielectronVarManager::kMultV0A] += vzeroData->GetMultiplicityV0A(i);
1230     values[AliDielectronVarManager::kMultV0C] += vzeroData->GetMultiplicityV0C(i);
1231     //values[AliDielectronVarManager::kAdcV0A] += vzeroData->GetAdcV0A(i);
1232     //values[AliDielectronVarManager::kAdcV0C] += vzeroData->GetAdcV0C(i);
1233   }
1234   values[AliDielectronVarManager::kMultV0] = values[AliDielectronVarManager::kMultV0A] + values[AliDielectronVarManager::kMultV0C];
1235   values[AliDielectronVarManager::kAdcV0] = values[AliDielectronVarManager::kAdcV0A] + values[AliDielectronVarManager::kAdcV0C];
1236   // VZERO event plane quantities
1237   Double_t qvec[3]={0.0};
1238   GetVzeroRP(event, qvec,0);      // V0-A
1239   values[AliDielectronVarManager::kV0AxH2] = qvec[0]; values[AliDielectronVarManager::kV0AyH2] = qvec[1]; 
1240   values[AliDielectronVarManager::kV0ArpH2] = qvec[2];
1241   qvec[0]=0.0; qvec[1]=0.0; qvec[2]=0.0;
1242   GetVzeroRP(event, qvec,1);      // V0-C
1243   values[AliDielectronVarManager::kV0CxH2] = qvec[0]; values[AliDielectronVarManager::kV0CyH2] = qvec[1]; 
1244   values[AliDielectronVarManager::kV0CrpH2] = qvec[2];
1245   qvec[0]=0.0; qvec[1]=0.0; qvec[2]=0.0;
1246   GetVzeroRP(event, qvec,2);      // V0-A and V0-C combined
1247   values[AliDielectronVarManager::kV0ACxH2] = qvec[0]; values[AliDielectronVarManager::kV0ACyH2] = qvec[1]; 
1248   values[AliDielectronVarManager::kV0ACrpH2] = qvec[2];
1249   // VZERO event plane resolution
1250   values[AliDielectronVarManager::kV0ArpResH2] = 1.0;
1251   values[AliDielectronVarManager::kV0CrpResH2] = 1.0;
1252   values[AliDielectronVarManager::kV0ACrpResH2] = 1.0;
1253   // Q vector components correlations  
1254   values[AliDielectronVarManager::kV0XaXcH2] = values[AliDielectronVarManager::kV0AxH2]*values[AliDielectronVarManager::kV0CxH2];
1255   values[AliDielectronVarManager::kV0XaYaH2] = values[AliDielectronVarManager::kV0AxH2]*values[AliDielectronVarManager::kV0AyH2];
1256   values[AliDielectronVarManager::kV0XaYcH2] = values[AliDielectronVarManager::kV0AxH2]*values[AliDielectronVarManager::kV0CyH2];
1257   values[AliDielectronVarManager::kV0YaXcH2] = values[AliDielectronVarManager::kV0AyH2]*values[AliDielectronVarManager::kV0CxH2];
1258   values[AliDielectronVarManager::kV0YaYcH2] = values[AliDielectronVarManager::kV0AyH2]*values[AliDielectronVarManager::kV0CyH2];
1259   values[AliDielectronVarManager::kV0XcYcH2] = values[AliDielectronVarManager::kV0CxH2]*values[AliDielectronVarManager::kV0CyH2];
1260
1261   // TPC event plane quantities (uncorrected)
1262   AliEventplane *evplane = const_cast<AliVEvent*>(event)->GetEventplane();
1263   TVector2 *qstd  = evplane->GetQVector();  // This is the "standard" Q-Vector
1264   TVector2 *qsub1 = evplane->GetQsub1();
1265   TVector2 *qsub2 = evplane->GetQsub2();
1266   if(!qstd || !qsub1 || !qsub2)  return;
1267      
1268   values[AliDielectronVarManager::kTPCxH2uc]   = qstd->X();
1269   values[AliDielectronVarManager::kTPCyH2uc]   = qstd->Y();
1270   values[AliDielectronVarManager::kTPCrpH2uc]  = ((TMath::Abs(qstd->X())>1.0e-10) ? TMath::ATan2(qstd->Y(),qstd->X())/2.0 : 0.0);
1271
1272   values[AliDielectronVarManager::kTPCsub1xH2uc]   = qsub1->X();
1273   values[AliDielectronVarManager::kTPCsub1yH2uc]   = qsub1->Y();
1274   values[AliDielectronVarManager::kTPCsub1rpH2uc]  = ((TMath::Abs(qsub1->X())>1.0e-10) ? TMath::ATan2(qsub1->Y(),qsub1->X())/2.0 : 0.0);
1275
1276   values[AliDielectronVarManager::kTPCsub2xH2uc]   = qsub2->X();
1277   values[AliDielectronVarManager::kTPCsub2yH2uc]   = qsub2->Y();
1278   values[AliDielectronVarManager::kTPCsub2rpH2uc]  = ((TMath::Abs(qsub2->X())>1.0e-10) ? TMath::ATan2(qsub2->Y(),qsub2->X())/2.0 : 0.0);
1279
1280   values[AliDielectronVarManager::kTPCsub12DiffH2uc] = TMath::Cos( 2.*(values[AliDielectronVarManager::kTPCsub1rpH2uc] - 
1281                                                                        values[AliDielectronVarManager::kTPCsub2rpH2uc]) );
1282   
1283   // using corrected tpc quantities
1284   values[AliDielectronVarManager::kV0ATPCDiffH2]   = TMath::Cos( 2.*(values[AliDielectronVarManager::kV0ArpH2] - 
1285                                                                      values[AliDielectronVarManager::kTPCrpH2]) ); 
1286   
1287   values[AliDielectronVarManager::kV0CTPCDiffH2]   = TMath::Cos( 2.*(values[AliDielectronVarManager::kV0CrpH2] - 
1288                                                                      values[AliDielectronVarManager::kTPCrpH2]) ); 
1289   
1290   values[AliDielectronVarManager::kV0AV0CDiffH2]   = TMath::Cos( 2.*(values[AliDielectronVarManager::kV0ArpH2] - 
1291                                                                      values[AliDielectronVarManager::kV0CrpH2]) ); 
1292
1293   values[AliDielectronVarManager::kv0ATPCDiffH2]   = TMath::Cos( 2.*(values[AliDielectronVarManager::kv0ArpH2] - 
1294                                                                      values[AliDielectronVarManager::kTPCrpH2]) ); 
1295   
1296   values[AliDielectronVarManager::kv0CTPCDiffH2]   = TMath::Cos( 2.*(values[AliDielectronVarManager::kv0CrpH2] - 
1297                                                                      values[AliDielectronVarManager::kTPCrpH2]) ); 
1298   
1299   values[AliDielectronVarManager::kv0Av0CDiffH2]   = TMath::Cos( 2.*(values[AliDielectronVarManager::kv0ArpH2] - 
1300                                                                      values[AliDielectronVarManager::kv0CrpH2]) ); 
1301
1302   values[AliDielectronVarManager::kMixingBin]=0;
1303 }
1304
1305 inline void AliDielectronVarManager::FillVarESDEvent(const AliESDEvent *event, Double_t * const values)
1306 {
1307   //
1308   // Fill event information available for histogramming into an array
1309   // 
1310   
1311   // Fill common AliVEvent interface information
1312   FillVarVEvent(event, values);
1313
1314   Double_t centralityF=-1; Double_t centralitySPD=-1;
1315   AliCentrality *esdCentrality = const_cast<AliESDEvent*>(event)->GetCentrality();
1316   if (esdCentrality) centralityF = esdCentrality->GetCentralityPercentile("V0M");
1317   if (esdCentrality) centralitySPD = esdCentrality->GetCentralityPercentile("CL1");
1318   
1319   // Fill AliESDEvent interface specific information
1320   const AliESDVertex *primVtx = event->GetPrimaryVertex();
1321   values[AliDielectronVarManager::kXRes]       = primVtx->GetXRes();
1322   values[AliDielectronVarManager::kYRes]       = primVtx->GetYRes();
1323   values[AliDielectronVarManager::kZRes]       = primVtx->GetZRes();
1324   values[AliDielectronVarManager::kCentrality] = centralityF;
1325   values[AliDielectronVarManager::kCentralitySPD] = centralitySPD;
1326   
1327   // Event multiplicity estimators
1328   Int_t nTrSPD05=0; Int_t nTrITSTPC05=0; Int_t nTrITSSA05=0;
1329   event->EstimateMultiplicity(nTrSPD05, nTrITSTPC05, nTrITSSA05, 0.5);
1330   values[AliDielectronVarManager::kNaccTrckltsEsd05] = nTrSPD05;
1331   values[AliDielectronVarManager::kNaccItsTpcEsd05] = nTrITSTPC05;
1332   values[AliDielectronVarManager::kNaccItsPureEsd05] = nTrITSSA05;
1333   values[AliDielectronVarManager::kNaccTrckltsEsd05Corr] = 
1334     AliDielectronHelper::GetNaccTrckltsCorrected(event,Double_t(nTrSPD05),values[AliDielectronVarManager::kZvPrim],0);
1335   values[AliDielectronVarManager::kNaccItsTpcEsd05Corr] = 
1336     AliDielectronHelper::GetNaccTrckltsCorrected(event,Double_t(nTrITSTPC05),values[AliDielectronVarManager::kZvPrim],3);
1337   values[AliDielectronVarManager::kNaccItsPureEsd05Corr] = 
1338     AliDielectronHelper::GetNaccTrckltsCorrected(event,Double_t(nTrITSSA05),values[AliDielectronVarManager::kZvPrim],6);
1339   
1340   Int_t nTrSPD10=0; Int_t nTrITSTPC10=0; Int_t nTrITSSA10=0;
1341   event->EstimateMultiplicity(nTrSPD10, nTrITSTPC10, nTrITSSA10, 1.0);
1342   values[AliDielectronVarManager::kNaccTrckltsEsd10] = nTrSPD10;
1343   values[AliDielectronVarManager::kNaccItsTpcEsd10] = nTrITSTPC10;
1344   values[AliDielectronVarManager::kNaccItsPureEsd10] = nTrITSSA10;
1345   values[AliDielectronVarManager::kNaccTrckltsEsd10Corr] =
1346     AliDielectronHelper::GetNaccTrckltsCorrected(event,Double_t(nTrSPD10),values[AliDielectronVarManager::kZvPrim],1);
1347   values[AliDielectronVarManager::kNaccItsTpcEsd10Corr] =
1348     AliDielectronHelper::GetNaccTrckltsCorrected(event,Double_t(nTrITSTPC10),values[AliDielectronVarManager::kZvPrim],4);
1349   values[AliDielectronVarManager::kNaccItsPureEsd10Corr] =
1350     AliDielectronHelper::GetNaccTrckltsCorrected(event,Double_t(nTrITSSA10),values[AliDielectronVarManager::kZvPrim],7); 
1351
1352   Int_t nTrSPD16=0; Int_t nTrITSTPC16=0; Int_t nTrITSSA16=0;
1353   event->EstimateMultiplicity(nTrSPD16, nTrITSTPC16, nTrITSSA16, 1.6);
1354   values[AliDielectronVarManager::kNaccTrckltsEsd16] = nTrSPD16;
1355   values[AliDielectronVarManager::kNaccItsTpcEsd16] = nTrITSTPC16;
1356   values[AliDielectronVarManager::kNaccItsPureEsd16] = nTrITSSA16;
1357   values[AliDielectronVarManager::kNaccTrckltsEsd16Corr] =
1358     AliDielectronHelper::GetNaccTrckltsCorrected(event,Double_t(nTrSPD16),values[AliDielectronVarManager::kZvPrim],2);
1359   values[AliDielectronVarManager::kNaccItsTpcEsd16Corr] =
1360     AliDielectronHelper::GetNaccTrckltsCorrected(event,Double_t(nTrITSTPC16),values[AliDielectronVarManager::kZvPrim],5);
1361   values[AliDielectronVarManager::kNaccItsPureEsd16Corr] =
1362     AliDielectronHelper::GetNaccTrckltsCorrected(event,Double_t(nTrITSSA16),values[AliDielectronVarManager::kZvPrim],8);
1363  
1364 }
1365   
1366 inline void AliDielectronVarManager::FillVarAODEvent(const AliAODEvent *event, Double_t * const values)
1367 {
1368   //
1369   // Fill event information available for histogramming into an array
1370   //   
1371
1372   // Fill common AliVEvent interface information
1373   FillVarVEvent(event, values);
1374
1375   // Fill AliAODEvent interface specific information
1376   AliAODHeader *header = event->GetHeader();
1377   
1378   Double_t centralityF=-1; Double_t centralitySPD=-1;
1379   AliCentrality *aodCentrality = header->GetCentralityP();
1380   if (aodCentrality) centralityF = aodCentrality->GetCentralityPercentile("V0M");
1381   if (aodCentrality) centralitySPD = aodCentrality->GetCentralityPercentile("CL1");
1382   values[AliDielectronVarManager::kCentrality] = centralityF;
1383   values[AliDielectronVarManager::kCentralitySPD] = centralitySPD;
1384   
1385   //const AliAODVertex *primVtx = event->GetPrimaryVertex();
1386 }
1387   
1388 inline void AliDielectronVarManager::FillVarMCEvent(const AliMCEvent *event, Double_t * const values)
1389
1390   //
1391   // Fill event information available for histogramming into an array
1392   //   
1393         
1394   // Fill common AliVEvent interface information
1395   //  FillVarVEvent(event, values);
1396   const AliVVertex* vtx = event->GetPrimaryVertex();
1397   values[AliDielectronVarManager::kXvPrim]       = (vtx ? vtx->GetX() : 0.0);
1398   values[AliDielectronVarManager::kYvPrim]       = (vtx ? vtx->GetY() : 0.0);
1399   values[AliDielectronVarManager::kZvPrim]       = (vtx ? vtx->GetZ() : 0.0);
1400   // Fill AliMCEvent interface specific information
1401   values[AliDielectronVarManager::kNch]   = AliDielectronHelper::GetNch(event, 1.6);
1402   values[AliDielectronVarManager::kNch05] = AliDielectronHelper::GetNch(event, 0.5);
1403   values[AliDielectronVarManager::kNch10] = AliDielectronHelper::GetNch(event, 1.0);
1404   
1405   values[AliDielectronVarManager::kNumberOfJPsis] = AliDielectronHelper::GetNMothers(event, 0.9, 443, 11);
1406   values[AliDielectronVarManager::kNumberOfJPsisPrompt]  = AliDielectronHelper::GetNMothers(event, 0.9, 443, 11, 1);
1407   values[AliDielectronVarManager::kNumberOfJPsisNPrompt] = AliDielectronHelper::GetNMothers(event, 0.9, 443, 11, 0);
1408 }
1409
1410 inline void AliDielectronVarManager::FillVarTPCEventPlane(const AliEventplane *evplane, Double_t * const values)
1411 {
1412   //
1413   // Fill TPC event plane information after correction
1414   //
1415   if(!evplane) return;
1416   TVector2 *qcorr  = const_cast<AliEventplane *>(evplane)->GetQVector();  // This is the "corrected" Q-Vector
1417   TVector2 *qcsub1 = const_cast<AliEventplane *>(evplane)->GetQsub1();
1418   TVector2 *qcsub2 = const_cast<AliEventplane *>(evplane)->GetQsub2();
1419   if(!qcorr || !qcsub1 || !qcsub2) return;
1420
1421   values[AliDielectronVarManager::kTPCxH2]   = qcorr->X();
1422   values[AliDielectronVarManager::kTPCyH2]   = qcorr->Y();
1423   values[AliDielectronVarManager::kTPCrpH2]  = ((TMath::Abs(qcorr->X())>1.0e-10) ? TMath::ATan2(qcorr->Y(),qcorr->X())/2.0 : 0.0);
1424  
1425   values[AliDielectronVarManager::kTPCsub1xH2]   = qcsub1->X();
1426   values[AliDielectronVarManager::kTPCsub1yH2]   = qcsub1->Y();
1427   values[AliDielectronVarManager::kTPCsub1rpH2]  = ((TMath::Abs(qcsub1->X())>1.0e-10) ? TMath::ATan2(qcsub1->Y(),qcsub1->X())/2.0 : 0.0);
1428
1429   values[AliDielectronVarManager::kTPCsub2xH2]   = qcsub2->X();
1430   values[AliDielectronVarManager::kTPCsub2yH2]   = qcsub2->Y();
1431   values[AliDielectronVarManager::kTPCsub2rpH2]  = ((TMath::Abs(qcsub2->X())>1.0e-10) ? TMath::ATan2(qcsub2->Y(),qcsub2->X())/2.0 : 0.0);
1432
1433   values[AliDielectronVarManager::kTPCsub12DiffH2] = TMath::Cos( 2.*(values[AliDielectronVarManager::kTPCsub1rpH2] - 
1434                                                                      values[AliDielectronVarManager::kTPCsub2rpH2]) );
1435   values[AliDielectronVarManager::kTPCsub12DiffH2Sin] = TMath::Sin( 2.*(values[AliDielectronVarManager::kTPCsub1rpH2] - 
1436                                                                         values[AliDielectronVarManager::kTPCsub2rpH2]) );
1437 }
1438   
1439
1440 inline void AliDielectronVarManager::InitESDpid(Int_t type)
1441 {
1442   //
1443   // initialize PID parameters
1444   // type=0 is simulation
1445   // type=1 is data
1446
1447   if (!fgPIDResponse) fgPIDResponse=new AliESDpid((Bool_t)(type==0));
1448   Double_t alephParameters[5];
1449   // simulation
1450   alephParameters[0] = 2.15898e+00/50.;
1451   alephParameters[1] = 1.75295e+01;
1452   alephParameters[2] = 3.40030e-09;
1453   alephParameters[3] = 1.96178e+00;
1454   alephParameters[4] = 3.91720e+00;
1455   fgPIDResponse->GetTOFResponse().SetTimeResolution(80.);
1456   
1457   // data
1458   if (type==1){    
1459     alephParameters[0] = 0.0283086/0.97;
1460     alephParameters[1] = 2.63394e+01;
1461     alephParameters[2] = 5.04114e-11;
1462     alephParameters[3] = 2.12543e+00;
1463     alephParameters[4] = 4.88663e+00;
1464     fgPIDResponse->GetTOFResponse().SetTimeResolution(130.);
1465     fgPIDResponse->GetTPCResponse().SetMip(50.);
1466   }
1467
1468   fgPIDResponse->GetTPCResponse().SetBetheBlochParameters(
1469     alephParameters[0],alephParameters[1],alephParameters[2],
1470     alephParameters[3],alephParameters[4]);
1471   
1472   fgPIDResponse->GetTPCResponse().SetSigma(3.79301e-03, 2.21280e+04);
1473 }
1474
1475 inline void AliDielectronVarManager::InitAODpidUtil(Int_t type)
1476 {
1477   if (!fgPIDResponse) fgPIDResponse=new AliAODpidUtil;
1478   Double_t alephParameters[5];
1479   // simulation
1480   alephParameters[0] = 2.15898e+00/50.;
1481   alephParameters[1] = 1.75295e+01;
1482   alephParameters[2] = 3.40030e-09;
1483   alephParameters[3] = 1.96178e+00;
1484   alephParameters[4] = 3.91720e+00;
1485   fgPIDResponse->GetTOFResponse().SetTimeResolution(80.);
1486   
1487   // data
1488   if (type==1){
1489     alephParameters[0] = 0.0283086/0.97;
1490     alephParameters[1] = 2.63394e+01;
1491     alephParameters[2] = 5.04114e-11;
1492     alephParameters[3] = 2.12543e+00;
1493     alephParameters[4] = 4.88663e+00;
1494     fgPIDResponse->GetTOFResponse().SetTimeResolution(130.);
1495     fgPIDResponse->GetTPCResponse().SetMip(50.);
1496   }
1497   
1498   fgPIDResponse->GetTPCResponse().SetBetheBlochParameters(
1499     alephParameters[0],alephParameters[1],alephParameters[2],
1500     alephParameters[3],alephParameters[4]);
1501   
1502   fgPIDResponse->GetTPCResponse().SetSigma(3.79301e-03, 2.21280e+04);
1503 }
1504
1505
1506 inline void AliDielectronVarManager::InitEstimatorAvg(const Char_t* filename)
1507 {
1508   //
1509   // initialize the profile histograms neccessary for the correction of the multiplicity estimators in pp collisions
1510   //
1511   
1512   const Char_t* estimatorNames[9] = {"SPDmult05","SPDmult10","SPDmult16",
1513                                      "ITSTPC05", "ITSTPC10", "ITSTPC16", 
1514                                      "ITSSA05",  "ITSSA10",  "ITSSA16"};
1515   const Char_t* periodNames[4] = {"LHC10b", "LHC10c", "LHC10d", "LHC10e"};
1516   TFile* file=TFile::Open(filename);
1517   if(!file) return;
1518   
1519   for(Int_t ip=0; ip<4; ++ip) {
1520     for(Int_t ie=0; ie<9; ++ie) {
1521       fgMultEstimatorAvg[ip][ie] = (TProfile*)(file->Get(Form("%s_%s",estimatorNames[ie],periodNames[ip]))->Clone(Form("%s_%s_clone",estimatorNames[ie],periodNames[ip])));
1522     }
1523   }
1524 }
1525
1526
1527 inline void AliDielectronVarManager::InitTRDpidEffHistograms(const Char_t* filename)
1528 {
1529   //
1530   // initialize the 3D histograms with the TRD pid efficiency histograms
1531   //
1532   
1533   // reset the centrality ranges and the efficiency histograms
1534   for(Int_t i=0; i<10; ++i) {         // centrality ranges
1535     for(Int_t j=0; j<4; ++j) fgTRDpidEffCentRanges[i][j] = -1.;
1536     if(fgTRDpidEff[i][0]) {
1537       delete fgTRDpidEff[i][0];
1538       fgTRDpidEff[i][0] = 0x0;
1539     }
1540     if(fgTRDpidEff[i][1]) {
1541       delete fgTRDpidEff[i][1];
1542       fgTRDpidEff[i][1] = 0x0;
1543     }
1544   }
1545   
1546   TFile* file=TFile::Open(filename);
1547   TList* keys=file->GetListOfKeys();
1548   Int_t idxp=0; Int_t idxn=0;
1549   for(Int_t i=0; i<keys->GetEntries(); ++i) {
1550     if(idxp>=10) continue;
1551     if(idxn>=10) continue;
1552     TString name=((TKey*)keys->At(i))->ReadObj()->GetName();    
1553     // Name of histograms should be in the format:
1554     // TRDeff<field>_cent_<centLow>_<centHigh>
1555     // <field> is either "BPLUS" or "BMINUS"
1556     if(!(name.Contains("BPLUS") || name.Contains("BMINUS"))) continue;
1557     TObjArray* arr = name.Tokenize("_");
1558     Bool_t isBplus = kTRUE;
1559     if(name.Contains("BMINUS")) isBplus = kFALSE;
1560     TString centMinStr = arr->At(2)->GetName(); TString centMaxStr = arr->At(3)->GetName();
1561     if(isBplus) {
1562       fgTRDpidEffCentRanges[idxp][2] = centMinStr.Atof();
1563       fgTRDpidEffCentRanges[idxp][3] = centMaxStr.Atof();
1564       fgTRDpidEff[idxp][1] = (TH3D*)(file->Get(name.Data())->Clone(Form("%s_clone",name.Data())));
1565       ++idxp;
1566     }
1567     else {
1568       fgTRDpidEffCentRanges[idxn][0] = centMinStr.Atof();
1569       fgTRDpidEffCentRanges[idxn][1] = centMaxStr.Atof();
1570       fgTRDpidEff[idxn][0] = (TH3D*)(file->Get(name.Data())->Clone(Form("%s_clone",name.Data())));
1571       ++idxn;
1572     }
1573   }
1574 }
1575
1576
1577 inline void AliDielectronVarManager::InitVZEROCalibrationHistograms(Int_t runNo) {
1578   //
1579   // Initialize the VZERO channel-by-channel calibration histograms
1580   //
1581
1582   //initialize only once
1583   if(fgVZEROCalib[0]) return;
1584   
1585   for(Int_t i=0; i<64; ++i) 
1586     if(fgVZEROCalib[i]) {
1587       delete fgVZEROCalib[i];
1588       fgVZEROCalib[i] = 0x0;
1589     }
1590   
1591   TFile file(fgVZEROCalibrationFile.Data());
1592   
1593   for(Int_t i=0; i<64; ++i){
1594     fgVZEROCalib[i] = (TProfile2D*)(file.Get(Form("RUN%d_ch%d_VtxCent", runNo, i)));
1595     if (fgVZEROCalib[i]) fgVZEROCalib[i]->SetDirectory(0x0);
1596   }
1597 }
1598
1599
1600 inline void AliDielectronVarManager::InitVZERORecenteringHistograms(Int_t runNo) {
1601   //
1602   // Initialize the VZERO event plane recentering histograms
1603   //
1604
1605   //initialize only once
1606   if(fgVZERORecentering[0][0]) return;
1607   
1608   for(Int_t i=0; i<2; ++i)
1609     for(Int_t j=0; j<2; ++j)
1610       if(fgVZERORecentering[i][j]) {
1611         delete fgVZERORecentering[i][j];
1612         fgVZERORecentering[i][j] = 0x0;
1613       }
1614   
1615   TFile file(fgVZERORecenteringFile.Data());
1616   if (!file.IsOpen()) return;
1617   
1618   fgVZERORecentering[0][0] = (TProfile2D*)(file.Get(Form("RUN%d_QxA_CentVtx", runNo)));
1619   fgVZERORecentering[0][1] = (TProfile2D*)(file.Get(Form("RUN%d_QyA_CentVtx", runNo)));
1620   fgVZERORecentering[1][0] = (TProfile2D*)(file.Get(Form("RUN%d_QxC_CentVtx", runNo)));
1621   fgVZERORecentering[1][1] = (TProfile2D*)(file.Get(Form("RUN%d_QyC_CentVtx", runNo)));
1622
1623   if (fgVZERORecentering[0][0]) fgVZERORecentering[0][0]->SetDirectory(0x0);
1624   if (fgVZERORecentering[0][1]) fgVZERORecentering[0][1]->SetDirectory(0x0);
1625   if (fgVZERORecentering[1][0]) fgVZERORecentering[1][0]->SetDirectory(0x0);
1626   if (fgVZERORecentering[1][1]) fgVZERORecentering[1][1]->SetDirectory(0x0);
1627   
1628 }
1629
1630
1631 inline Double_t AliDielectronVarManager::GetTRDpidEfficiency(Int_t runNo, Double_t centrality, 
1632                                                              Double_t eta, Double_t trdPhi, Double_t pout,
1633                                                              Double_t& effErr) {
1634   //
1635   // return the efficiency in the given phase space cell
1636   //
1637   // LHC10h data----------------------------------------------
1638   Bool_t isBplus = kTRUE;
1639   if(runNo<=138275) isBplus = kFALSE;
1640   // TODO: check magnetic polarity for runs in 2011 data
1641   // ---------------------------------------------------------
1642   Int_t centIdx = -1;
1643   for(Int_t icent=0; icent<10; ++icent) {
1644     if(isBplus) {
1645       if(centrality>=fgTRDpidEffCentRanges[icent][2] && centrality<fgTRDpidEffCentRanges[icent][3]) {
1646         centIdx = icent;
1647         break;
1648       }
1649     }
1650     else {
1651       if(centrality>=fgTRDpidEffCentRanges[icent][0] && centrality<fgTRDpidEffCentRanges[icent][1]) {
1652         centIdx = icent;
1653         break;
1654       }
1655     }
1656   }
1657   //TODO: chek logick
1658   if (centIdx<0) return 1;
1659   
1660   TH3D* effH = fgTRDpidEff[centIdx][(isBplus ? 1 : 0)];
1661   if(!effH) {effErr=0x0; return 1.0;}
1662   Int_t etaBin = effH->GetXaxis()->FindBin(eta);
1663   if(eta<effH->GetXaxis()->GetXmin()) etaBin=1;
1664   if(eta>effH->GetXaxis()->GetXmax()) etaBin=effH->GetXaxis()->GetNbins();
1665   Int_t phiBin = effH->GetYaxis()->FindBin(trdPhi);
1666   if(trdPhi<effH->GetYaxis()->GetXmin()) phiBin=1;
1667   if(trdPhi>effH->GetYaxis()->GetXmax()) phiBin=effH->GetYaxis()->GetNbins();
1668   Int_t poutBin = effH->GetZaxis()->FindBin(pout);
1669   if(pout<effH->GetZaxis()->GetXmin()) poutBin=1;
1670   if(pout>effH->GetZaxis()->GetXmax()) poutBin=effH->GetZaxis()->GetNbins();
1671   Double_t eff = effH->GetBinContent(etaBin, phiBin, poutBin);
1672   effErr = effH->GetBinError(etaBin, phiBin, poutBin);
1673   if(eff<-0.0001) {
1674     effErr = 0.0;
1675     eff = 1.0;
1676   }
1677   return eff;
1678 }
1679
1680
1681 inline void AliDielectronVarManager::SetEvent(AliVEvent * const ev)
1682 {
1683   
1684   fgEvent = ev;
1685   if (fgKFVertex) delete fgKFVertex;
1686   fgKFVertex=0x0;
1687   if (ev && ev->GetPrimaryVertex()) fgKFVertex=new AliKFVertex(*ev->GetPrimaryVertex());
1688
1689   for (Int_t i=0; i<AliDielectronVarManager::kNMaxValues;++i) fgData[i]=0.;
1690   AliDielectronVarManager::Fill(fgEvent, fgData);
1691 }
1692
1693 inline void AliDielectronVarManager::SetEventData(const Double_t data[AliDielectronVarManager::kNMaxValues])
1694 {
1695   for (Int_t i=0; i<kNMaxValues;++i) fgData[i]=0.;
1696   for (Int_t i=kPairMax; i<kNMaxValues;++i) fgData[i]=data[i];
1697 }
1698
1699
1700 inline Bool_t AliDielectronVarManager::GetDCA(const AliAODTrack *track, Double_t d0z0[2])
1701 {
1702   if(track->TestBit(AliAODTrack::kIsDCA)){
1703     d0z0[0]=track->DCA();
1704     d0z0[1]=track->ZAtDCA();
1705     return kTRUE;
1706   }
1707   
1708   Double_t covd0z0[3];
1709   AliAODTrack copy(*track);
1710   AliAODVertex *vtx =(AliAODVertex*)(fgEvent->GetPrimaryVertex());
1711   Double_t fBzkG = fgEvent->GetMagneticField(); // z componenent of field in kG
1712   Bool_t ok = copy.PropagateToDCA(vtx,fBzkG,kVeryBig,d0z0,covd0z0);
1713   if(!ok){
1714     d0z0[0]=-999.;
1715     d0z0[1]=-999.;
1716   }
1717   return ok;
1718 }
1719
1720 inline void AliDielectronVarManager::SetTPCEventPlane(AliEventplane *const evplane)
1721 {
1722   
1723   fgTPCEventPlane = evplane;
1724   //  for (Int_t i=0; i<AliDielectronVarManager::kNMaxValues;++i) fgData[i]=0.;
1725   //  AliDielectronVarManager::Fill(fgEvent, fgData);
1726 }
1727
1728
1729 //_________________________________________________________________
1730 inline void AliDielectronVarManager::GetVzeroRP(const AliVEvent* event, Double_t* qvec, Int_t sideOption) {
1731   //
1732   // Get the reaction plane from the VZERO detector for a given harmonic
1733   //
1734   // sideOption = 0- V0A, 1- V0C, 2-both
1735   //  Q{x,y} = SUM_i mult(i) * {cos(n*phi_i), sin(n*phi_i)} 
1736   //  phi_i - phi angle of the VZERO sector i
1737   //          Each sector covers 45 degrees(8 sectors per ring). Middle of sector 0 is at 45/2
1738   //        channel 0: 22.5
1739   //                1: 22.5+45
1740   //                2: 22.5+45*2
1741   //               ...
1742   //        at the next ring continues the same
1743   //        channel 8: 22.5
1744   //        channel 9: 22.5 + 45
1745   //               ... 
1746   const Double_t kX[8] = {0.92388, 0.38268, -0.38268, -0.92388, -0.92388, -0.38268, 0.38268, 0.92388};    // cosines of the angles of the VZERO sectors (8 per ring)
1747   const Double_t kY[8] = {0.38268, 0.92388, 0.92388, 0.38268, -0.38268, -0.92388, -0.92388, -0.38268};    // sines     -- " --
1748   Int_t phi;
1749   Float_t mult;
1750   
1751   // get centrality and vertex for this event
1752   Double_t centralitySPD = -1; Double_t vtxZ = -999.;
1753   if(event->IsA() == AliESDEvent::Class()) {
1754     const AliESDEvent* esdEv = static_cast<const AliESDEvent*>(event);
1755     AliCentrality *esdCentrality = const_cast<AliESDEvent*>(esdEv)->GetCentrality();
1756     if(esdCentrality) centralitySPD = esdCentrality->GetCentralityPercentile("CL1");
1757   }
1758   if(event->IsA() == AliAODEvent::Class()) {
1759     const AliAODEvent* aodEv = static_cast<const AliAODEvent*>(event);
1760     AliAODHeader *header = aodEv->GetHeader();
1761     AliCentrality *aodCentrality = header->GetCentralityP();
1762     if(aodCentrality) centralitySPD = aodCentrality->GetCentralityPercentile("CL1");
1763   }
1764   const AliVVertex *primVtx = event->GetPrimaryVertex();
1765   if(!primVtx) return;
1766   vtxZ = primVtx->GetZ();
1767   if(TMath::Abs(vtxZ)>10.) return;
1768   if(centralitySPD<0. || centralitySPD>80.) return;
1769   
1770   Int_t binCent = -1; Int_t binVtx = -1;
1771   if(fgVZEROCalib[0]) {
1772     binVtx = fgVZEROCalib[0]->GetXaxis()->FindBin(vtxZ);
1773     binCent = fgVZEROCalib[0]->GetYaxis()->FindBin(centralitySPD);
1774   }
1775   
1776   AliVVZERO* vzero = event->GetVZEROData();
1777   Double_t average = 0.0;
1778   for(Int_t iChannel=0; iChannel<64; ++iChannel) {
1779     if(iChannel<32 && sideOption==0) continue;
1780     if(iChannel>=32 && sideOption==1) continue;
1781     phi=iChannel%8;
1782     mult = vzero->GetMultiplicity(iChannel);
1783     if(fgVZEROCalib[iChannel])
1784       average = fgVZEROCalib[iChannel]->GetBinContent(binVtx, binCent);
1785     if(average>1.0e-10 && mult>0.5) 
1786       mult /= average;
1787     else
1788       mult = 0.0;
1789     //  2nd harmonic
1790     qvec[0] += mult*(2.0*TMath::Power(kX[phi],2.0)-1);
1791     qvec[1] += mult*(2.0*kX[phi]*kY[phi]);
1792   }    // end loop over channels 
1793   
1794   // do recentering
1795   if(fgVZERORecentering[0][0]) {
1796 //     printf("vzero: %p\n",fgVZERORecentering[0][0]);
1797     Int_t binCentRecenter = -1; Int_t binVtxRecenter = -1;
1798     binCentRecenter = fgVZERORecentering[0][0]->GetXaxis()->FindBin(centralitySPD);
1799     binVtxRecenter = fgVZERORecentering[0][0]->GetYaxis()->FindBin(vtxZ);
1800     if(sideOption==0) {  // side A
1801       qvec[0] -= fgVZERORecentering[0][0]->GetBinContent(binCentRecenter, binVtxRecenter);
1802       qvec[1] -= fgVZERORecentering[0][1]->GetBinContent(binCentRecenter, binVtxRecenter);
1803     }
1804     if(sideOption==1) {  // side C
1805       qvec[0] -= fgVZERORecentering[1][0]->GetBinContent(binCentRecenter, binVtxRecenter);
1806       qvec[1] -= fgVZERORecentering[1][1]->GetBinContent(binCentRecenter, binVtxRecenter);
1807     }
1808     if(sideOption==2) {  // side A and C together
1809       qvec[0] -= fgVZERORecentering[0][0]->GetBinContent(binCentRecenter, binVtxRecenter);
1810       qvec[0] -= fgVZERORecentering[1][0]->GetBinContent(binCentRecenter, binVtxRecenter);
1811       qvec[1] -= fgVZERORecentering[0][1]->GetBinContent(binCentRecenter, binVtxRecenter);
1812       qvec[1] -= fgVZERORecentering[1][1]->GetBinContent(binCentRecenter, binVtxRecenter);
1813     }
1814   }
1815   
1816   // calculate the reaction plane
1817   if(TMath::Abs(qvec[0])>1.0e-10)
1818     qvec[2] = TMath::ATan2(qvec[1],qvec[0])/2.0;
1819 }
1820
1821
1822
1823 /*
1824 inline void AliDielectronVarManager::FillValues(const TParticle *particle, Double_t *values)
1825 {
1826   //
1827   // Fill track information available for histogramming into an array
1828   //
1829
1830   // Fill TParticle interface information
1831   values[AliDielectronVarManager::kPx]     = particle->Px();
1832   values[AliDielectronVarManager::kPy]     = particle->Py();
1833   values[AliDielectronVarManager::kPz]     = particle->Pz();
1834   values[AliDielectronVarManager::kPt]     = particle->Pt();
1835   values[AliDielectronVarManager::kP]      = particle->P();
1836
1837   values[AliDielectronVarManager::kXv]     = particle->Vx();
1838   values[AliDielectronVarManager::kYv]     = particle->Vy();
1839   values[AliDielectronVarManager::kZv]     = particle->Vz();
1840
1841   values[AliDielectronVarManager::kOneOverPt] = 1./particle->Pt();
1842   values[AliDielectronVarManager::kPhi]    = particle->Phi();
1843   values[AliDielectronVarManager::kTheta]  = 
1844   values[AliDielectronVarManager::kEta]    = particle->Eta();
1845   values[AliDielectronVarManager::kY]      = 
1846
1847   values[AliDielectronVarManager::kE]      = particle->Energy();
1848   values[AliDielectronVarManager::kM]      = particle->GetMass();
1849
1850   values[AliDielectronVarManager::kCharge] = particle->GetPDG()->Charge()/3; // uggly
1851
1852 }*/
1853
1854 #endif
1855