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