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