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