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