]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/AliAnalysisTaskV0QA.cxx
d00c05401017597b1a27913827b5619f6fe9721a
[u/mrichter/AliRoot.git] / PWG1 / AliAnalysisTaskV0QA.cxx
1 #define AliAnalysisTaskV0QA_cxx
2 #include "Riostream.h"
3 #include "TROOT.h"
4 #include "TChain.h"
5 #include "TH1.h"
6 #include "TH2.h"
7
8 #include "TCanvas.h"
9 #include "TSystem.h"
10 #include "TLorentzVector.h"
11
12 #include "AliAnalysisTask.h"
13 #include "AliTrackReference.h"
14 #include "AliKFParticle.h"
15 #include "AliKFVertex.h"
16 #include "AliVertexerTracks.h"
17 #include "AliESDEvent.h"
18 #include "AliESDInputHandler.h"
19 #include "AliMCEvent.h"
20 #include "AliMCEventHandler.h"
21 #include "AliStack.h"
22 #include "AliESD.h"
23 #include "AliLog.h"
24
25 #include "AliAnalysisTaskV0QA.h"
26
27 ClassImp(AliAnalysisTaskV0QA)
28
29 //________________________________________________________________________
30 AliAnalysisTaskV0QA::AliAnalysisTaskV0QA(const char *name) :AliAnalysisTask(name,""), 
31 fESD(0), 
32 stack(0),
33 mctruth(0),
34 fChain(0),
35 fOutputContainer(0),
36 fSparseV0(0),
37 fSparseK0(0),
38 fSparseL(0),
39 fSparseAL(0),
40 nEv(0),
41 nConvGamGeant(-1),
42 gConvGamGeantIndex(0),
43 eNegConvGamGeantIndex(0),
44 ePosConvGamGeantIndex(0),
45 eNegConvGamGeantLength(0),
46 ePosConvGamGeantLength(0),
47 eNegConvGamSingleRecIndex(0),
48 ePosConvGamSingleRecIndex(0),
49 eNegConvGamV0RecIndex(0),
50 ePosConvGamV0RecIndex(0),
51 ConvGamV0RecIndexPos(0),
52 ConvGamV0RecIndexNeg(0),
53 gDim(50),
54 nDecayLGeant(-1),
55 lDecayLGeantIndex(0),
56 piNegDecayLGeantIndex(0),
57 pPosDecayLGeantIndex(0),
58 piNegDecayLGeantLength(0),
59 pPosDecayLGeantLength(0),
60 piNegDecayLSingleRecIndex(0),
61 pPosDecayLSingleRecIndex(0),
62 piNegDecayLV0RecIndex(0),
63 pPosDecayLV0RecIndex(0),
64 DecayLV0RecIndexPos(0),
65 DecayLV0RecIndexNeg(0),
66 nDecayALGeant(-1),
67 alDecayALGeantIndex(0),
68 piPosDecayALGeantIndex(0),
69 apNegDecayALGeantIndex(0),
70 piPosDecayALGeantLength(0),
71 apNegDecayALGeantLength(0),
72 piPosDecayALSingleRecIndex(0),
73 apNegDecayALSingleRecIndex(0),
74 piPosDecayALV0RecIndex(0),
75 apNegDecayALV0RecIndex(0),
76 DecayALV0RecIndexPos(0),
77 DecayALV0RecIndexNeg(0),
78 nDecayK0Geant(-1),
79 K0DecayK0GeantIndex(0),
80 piNegDecayK0GeantIndex(0),
81 piPosDecayK0GeantIndex(0),
82 piNegDecayK0GeantLength(0),
83 piPosDecayK0GeantLength(0),
84 piNegDecayK0SingleRecIndex(0),
85 piPosDecayK0SingleRecIndex(0),
86 piNegDecayK0V0RecIndex(0),
87 piPosDecayK0V0RecIndex(0),
88 DecayK0V0RecIndexPos(0),
89 DecayK0V0RecIndexNeg(0),
90 piPosK0Index(-1),
91 piNegK0Index(-1),
92 nTracksPrim(-1),
93 tpcRefit(0),
94 itsRefit(0),
95 trdRefit(0),
96 trdOut(0),
97 fValueL(0),
98 fValueAL(0),
99 fValueK0(0),
100 fValueV0(0),
101 xminV0(0),
102 xmaxV0(0),
103 binsV0(0),
104 fDim(37),
105 fRefTPC(0),
106 clRefsN(0),
107 clRefsP(0)
108
109  {
110   // Constructor.
111   // Input slot #0 works with an Ntuple
112   DefineInput(0, TChain::Class());
113   // Output slot #0 writes into a TH1 container
114   //  DefineOutput(0,TObjArray::Class());
115   DefineOutput(0,TList::Class());
116
117   // Reconstructed arrays
118
119   nEv=0;
120   fDim=37; 
121   fValueK0 = new Double_t[fDim];
122   fValueL = new Double_t[fDim];
123   fValueAL = new Double_t[fDim];
124   fValueV0 = new Double_t[fDim];
125   xminV0 = new Double_t[fDim];
126   xmaxV0 = new Double_t[fDim];
127   binsV0 = new Int_t[fDim];
128
129
130   gDim=50; 
131   gConvGamGeantIndex = new Int_t[gDim];
132   eNegConvGamGeantIndex = new Int_t[gDim];
133   ePosConvGamGeantIndex = new Int_t[gDim];
134   eNegConvGamGeantLength = new Float_t[gDim];
135   ePosConvGamGeantLength = new Float_t[gDim];
136
137   eNegConvGamSingleRecIndex = new Int_t[gDim];
138   ePosConvGamSingleRecIndex = new Int_t[gDim];
139
140   eNegConvGamV0RecIndex = new Int_t[gDim];
141   ePosConvGamV0RecIndex = new Int_t[gDim];
142
143   ConvGamV0RecIndexPos = new Int_t[gDim];
144   ConvGamV0RecIndexNeg = new Int_t[gDim];
145
146   // Lambda to proton pi-
147   lDecayLGeantIndex = new Int_t[gDim];
148   piNegDecayLGeantIndex = new Int_t[gDim];
149   pPosDecayLGeantIndex = new Int_t[gDim];
150   piNegDecayLGeantLength = new Float_t[gDim];
151   pPosDecayLGeantLength = new Float_t[gDim];
152
153   piNegDecayLSingleRecIndex = new Int_t[gDim];
154   pPosDecayLSingleRecIndex = new Int_t[gDim];
155
156   piNegDecayLV0RecIndex = new Int_t[gDim];
157   pPosDecayLV0RecIndex = new Int_t[gDim];
158
159   DecayLV0RecIndexPos = new Int_t[gDim];
160   DecayLV0RecIndexNeg = new Int_t[gDim];
161
162   //K0S to pi+ pi-
163   K0DecayK0GeantIndex = new Int_t[gDim];
164   piNegDecayK0GeantIndex = new Int_t[gDim];
165   piPosDecayK0GeantIndex = new Int_t[gDim];
166   piNegDecayK0GeantLength = new Float_t[gDim];
167   piPosDecayK0GeantLength = new Float_t[gDim];
168
169   piNegDecayK0SingleRecIndex = new Int_t[gDim];
170   piPosDecayK0SingleRecIndex = new Int_t[gDim];
171
172   piNegDecayK0V0RecIndex = new Int_t[gDim];
173   piPosDecayK0V0RecIndex = new Int_t[gDim];
174
175   DecayK0V0RecIndexPos = new Int_t[gDim];
176   DecayK0V0RecIndexNeg = new Int_t[gDim];
177
178   //Antilambda to antiproton piplus
179   alDecayALGeantIndex = new Int_t[gDim];
180   piPosDecayALGeantIndex = new Int_t[gDim];
181   apNegDecayALGeantIndex = new Int_t[gDim];
182   piPosDecayALGeantLength = new Float_t[gDim];
183   apNegDecayALGeantLength = new Float_t[gDim];
184
185   piPosDecayALSingleRecIndex = new Int_t[gDim];
186   apNegDecayALSingleRecIndex = new Int_t[gDim];
187
188   piPosDecayALV0RecIndex = new Int_t[gDim];
189   apNegDecayALV0RecIndex = new Int_t[gDim];
190
191   DecayALV0RecIndexPos = new Int_t[gDim];
192   DecayALV0RecIndexNeg = new Int_t[gDim];
193
194
195   //////////////////////////////////////////////
196   clRefsP = new TClonesArray("AliTrackReference");
197   clRefsN = new TClonesArray("AliTrackReference");
198
199   //  SetESDtrackCuts();
200
201
202   AliLog::SetGlobalLogLevel(AliLog::kError);
203
204 }
205
206 //________________________________________________________________________
207 void AliAnalysisTaskV0QA::ConnectInputData(Option_t *) {
208   printf("   ConnectInputData %s\n", GetName());
209
210
211   
212   AliESDInputHandler* esdH = (AliESDInputHandler*) ((AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler());
213   fESD = esdH->GetEvent();
214   
215   fChain = (TChain*)GetInputData(0);
216
217 }
218
219 //_____________________________________________________
220 AliAnalysisTaskV0QA::~AliAnalysisTaskV0QA()
221 {
222   // Remove all pointers
223
224
225   delete [] clRefsP;
226   delete [] clRefsN;
227
228
229   delete [] fValueK0;
230   delete [] fValueL;
231   delete [] fValueAL;
232   delete [] fValueV0;
233   delete [] binsV0;
234   delete [] xminV0;
235   delete [] xmaxV0;
236
237   delete [] gConvGamGeantIndex;
238   delete [] eNegConvGamGeantIndex;
239   delete [] ePosConvGamGeantIndex;
240
241   delete [] eNegConvGamSingleRecIndex;
242   delete [] ePosConvGamSingleRecIndex;
243
244   delete [] eNegConvGamV0RecIndex;
245   delete [] ePosConvGamV0RecIndex;
246   delete [] ConvGamV0RecIndexPos;
247   delete [] ConvGamV0RecIndexNeg;
248
249   delete [] lDecayLGeantIndex;
250   delete [] piNegDecayLGeantIndex;
251   delete [] pPosDecayLGeantIndex;
252
253   delete [] piNegDecayLGeantLength;
254   delete [] pPosDecayLGeantLength;
255   delete [] piNegDecayLSingleRecIndex;
256   delete [] pPosDecayLSingleRecIndex;
257
258   delete [] piNegDecayLV0RecIndex;
259   delete [] pPosDecayLV0RecIndex;
260   delete [] DecayLV0RecIndexPos;
261   delete [] DecayLV0RecIndexNeg;
262
263   delete [] alDecayALGeantIndex;
264   delete [] piPosDecayALGeantIndex;
265   delete [] apNegDecayALGeantIndex;
266
267   delete [] piPosDecayALGeantLength;
268   delete [] apNegDecayALGeantLength;
269   delete [] piPosDecayALSingleRecIndex;
270   delete [] apNegDecayALSingleRecIndex;
271
272   delete [] piPosDecayALV0RecIndex;
273   delete [] apNegDecayALV0RecIndex;
274   delete [] DecayALV0RecIndexPos;
275   delete [] DecayALV0RecIndexNeg;
276
277
278   delete [] piNegDecayK0GeantIndex;
279   delete [] piPosDecayK0GeantIndex;
280
281   delete [] piNegDecayK0GeantLength;
282   delete [] piPosDecayK0GeantLength;
283   delete [] piNegDecayK0SingleRecIndex;
284   delete [] piPosDecayK0SingleRecIndex;
285
286   delete [] piNegDecayK0V0RecIndex;
287   delete [] piPosDecayK0V0RecIndex;
288
289   delete [] DecayK0V0RecIndexPos;
290   delete [] DecayK0V0RecIndexNeg;
291
292 }
293
294
295 //________________________________________________________________________
296 void AliAnalysisTaskV0QA::CreateOutputObjects() {
297
298   for(Int_t d=0;d<fDim;d++){
299     binsV0[d]=70;
300   }
301   xminV0[0]=   0;     // 1/sqrt(pt) Gamma geant
302   xmaxV0[0]=   8;
303
304
305   xminV0[1]=-2.5;   // eta Gamma Geant
306   xmaxV0[1]= 1.5;
307
308
309   xminV0[2]=-2*TMath::Pi();   // phi Gamma geant
310   xmaxV0[2]= TMath::Pi();
311
312
313   xminV0[3]=   0;     // r geant
314   xmaxV0[3]= 200;
315
316
317   xminV0[4]=-250;   // z geant
318   xmaxV0[4]= 250;
319
320
321   xminV0[5]=   0;     // 1/sqrt(pt) Geant Pos
322   xmaxV0[5]=   8;
323
324
325   xminV0[6]=-2.5;   // eta geant Pos
326   xmaxV0[6]= 1.5;
327
328
329   xminV0[7]=-2*TMath::Pi();   // phi Geant Pos
330   xmaxV0[7]= TMath::Pi();
331
332   
333   xminV0[8]=0;   // Track Length TPC Geant Pos
334   xmaxV0[8]= 200;
335
336
337   xminV0[9]=   0;     // 1/sqrt(pt) Geant Neg
338   xmaxV0[9]=   8;
339
340
341   xminV0[10]=-2.5;   // eta Geant Neg
342   xmaxV0[10]= 1.5;
343
344
345   xminV0[11]=-2*TMath::Pi();   // phi Geant Neg
346   xmaxV0[11]= TMath::Pi();
347
348
349   xminV0[12]=0;   // Track Length TPC Geant Neg
350   xmaxV0[12]= 200;
351
352
353
354   //-----------Rec single variables
355
356   xminV0[13]=   -0.5;     // (pt-ptGeant)/ptGeant rec Pos
357   xmaxV0[13]=   0.5;
358
359
360   xminV0[14]=-2.5;   // eta  rec Pos
361   xmaxV0[14]= 1.5;
362
363
364   xminV0[15]=-2*TMath::Pi();   // phi rec Pos
365   xmaxV0[15]= TMath::Pi();
366
367
368   xminV0[16]=   0;     // Impact parameter rec Pos
369   xmaxV0[16]=   100;
370
371
372   xminV0[17]=   0;     // nsigmas Impact parameter rec Pos
373   xmaxV0[17]=   100;
374
375
376
377   xminV0[18]=   -1;     // Ncls ITS rec Pos
378   xmaxV0[18]=   6;
379
380
381   xminV0[19]=   -1;     // Ncls TPC rec Pos
382   xmaxV0[19]=   180;
383
384
385   xminV0[20]=   -2;     // Status Single  TPC rec Pos
386   xmaxV0[20]=   2;
387
388
389   xminV0[21]=   -0.5;     // (pt-ptGeant)/ptGeant rec Neg
390   xmaxV0[21]=   0.5;
391
392
393   xminV0[22]=-2.5;   // eta  rec Neg
394   xmaxV0[22]= 1.5;
395
396
397   xminV0[23]=-2*TMath::Pi();   // phi rec Neg
398   xmaxV0[23]= TMath::Pi();
399
400
401   xminV0[24]=   0;     // Impact parameter rec Neg
402   xmaxV0[24]=   100;
403
404
405   xminV0[25]=   0;     // Sigmas Impact parameter rec Neg
406   xmaxV0[25]=   100;
407
408
409
410   xminV0[26]=   -1;     // Ncls ITS rec Neg
411   xmaxV0[26]=   6;
412
413
414   xminV0[27]=   -1;     // Ncls TPC rec Neg
415   xmaxV0[27]=   180;
416
417
418   xminV0[28]=   -2;     // Status Single  TPC rec Neg
419   xmaxV0[28]=   2;
420
421   // ------------------Rec V0 variables
422
423
424
425   xminV0[29]=   -0.5;     // (pt-ptGeant)/ptGeant rec V0 Pos
426   xmaxV0[29]=   0.5;
427
428
429   xminV0[30]=-2.5;   // eta  rec V0 Pos
430   xmaxV0[30]= 1.5;
431
432
433   xminV0[31]=-2*TMath::Pi();   // phi rec V0 Pos
434   xmaxV0[31]= TMath::Pi();
435
436
437   xminV0[32]=   -2;     // Status V0 TPC rec Pos
438   xmaxV0[32]=   2;
439
440
441
442
443   xminV0[33]=   -0.5;     // 1/sqrt(pt) rec V0 Neg
444   xmaxV0[33]=   0.5;
445
446
447   xminV0[34]=-2.5;   // eta  rec V0 Neg
448   xmaxV0[34]= 1.5;
449
450
451   xminV0[35]=-2*TMath::Pi();   // phi rec V0 Neg
452   xmaxV0[35]= TMath::Pi();
453
454
455
456   xminV0[36]=   -2;     // Status V0 TPC rec Neg
457   xmaxV0[36]=   2;
458
459
460
461   TString  axisName[37]={"ptGammaGeant", 
462                          "etaGammaGeant",
463                          "phiGammaGeant",
464                          "rGeant",
465                          "zGeant",
466                          "ptEPlusGeant",  
467                          "etaEPlusGeant",
468                          "phiEPlusGeant",
469                          "TPCTrackLengthEPlusGeant",
470                          "ptEMinusGeant", 
471                          "etaEMinusGeant",
472                          "phiEMinusGeant",
473                          "TPCTrackLengthEMinusGeant",
474                          "ptResEPlusRecSingle", 
475                          "etaEPlusRecSingle",
476                          "phiEPlusRecSingle",
477                          "bXYZEPlusRecSingle",
478                          "sigbXYZEPlusRecSingle",
479                          "NclsITSEPlusRecSingle",
480                          "NclsTPCEPlusRecSingle",
481                          "statusRecSinglePos",
482                          "ptResEMinusRecSingle", 
483                          "etaEMinusRecSingle",
484                          "phiEMinusRecSingle",
485                          "bXYZEMinusRecSingle",
486                          "sigbXYZEMinusRecSingle",
487                          "NclsITSEMinusRecSingle",
488                          "NclsTPCEMinusRecSingle",
489                          "statusRecSingleNeg",
490                          "ptResEPlusRecV0", 
491                          "etaEPlusRecV0",
492                          "phiEPlusRecV0",
493                          "statusV0SinglePos",
494                          "ptResEMinusRecV0", 
495                          "etaEMinusRecV0",
496                          "phiEMinusRecV0",
497                          "statusV0SingleNeg"};
498
499
500   fSparseV0= new THnSparseF("sparseV0","sparseV0",fDim,binsV0,xminV0,xmaxV0);
501
502   for (Int_t iaxis=0; iaxis<fDim; iaxis++){
503    fSparseV0->GetAxis(iaxis)->SetName(axisName[iaxis]);
504    fSparseV0->GetAxis(iaxis)->SetTitle(axisName[iaxis]);
505   }
506
507   TString  axisNameK0[37]={"ptK0Geant", 
508                          "etaK0Geant",
509                          "phiK0Geant",
510                          "rGeant",
511                          "zGeant",
512                          "ptPiPlusGeant",  
513                          "etaPiPlusGeant",
514                          "phiPiPlusGeant",
515                          "TPCTrackLengthPiPlusGeant",
516                          "ptPiMinusGeant", 
517                          "etaPiMinusGeant",
518                          "phiPiMinusGeant",
519                          "TPCTrackLengthPiMinusGeant",
520                          "ptResPiPlusRecSingle", 
521                          "etaPiPlusRecSingle",
522                          "phiPiPlusRecSingle",
523                          "bXYZPiPlusRecSingle",
524                          "sigbXYZPiPlusRecSingle",
525                          "NclsITSPiPlusRecSingle",
526                          "NclsTPCPiPlusRecSingle",
527                          "statusRecSinglePos",
528                          "ptResPiMinusRecSingle", 
529                          "etaPiMinusRecSingle",
530                          "phiPiMinusRecSingle",
531                          "bXYZPiMinusRecSingle",
532                          "sigbXYZPiMinusRecSingle",
533                          "NclsITSPiMinusRecSingle",
534                          "NclsTPCPiMinusRecSingle",
535                          "statusRecSingleNeg",
536                          "ptResPiPlusRecV0", 
537                          "etaPiPlusRecV0",
538                          "phiPiPlusRecV0",
539                          "statusRecV0Pos",
540                          "ptResPiMinusRecV0", 
541                          "etaPiMinusRecV0",
542                          "phiPiMinusRecV0",
543                          "statusRecV0Neg"};
544
545
546
547   fSparseK0= new THnSparseF("sparseK0","sparseK0",fDim,binsV0,xminV0,xmaxV0);
548   for (Int_t iaxis=0; iaxis<fDim; iaxis++){
549    fSparseK0->GetAxis(iaxis)->SetName(axisNameK0[iaxis]);
550    fSparseK0->GetAxis(iaxis)->SetTitle(axisNameK0[iaxis]);
551   }
552
553   TString  axisNameL[37]={"ptLGeant", 
554                          "etaLGeant",
555                          "phiLGeant",
556                          "rGeant",
557                          "zGeant",
558                          "ptPPlusGeant",  
559                          "etaPPlusGeant",
560                          "phiPPlusGeant",
561                          "TPCTrackLengthPPlusGeant",
562                          "ptPiMinusGeant", 
563                          "etaPiMinusGeant",
564                          "phiPiMinusGeant",
565                          "TPCTrackLengthPiMinusGeant",
566                          "ptResPPlusRecSingle", 
567                          "etaPPlusRecSingle",
568                          "phiPPlusRecSingle",
569                          "bXYZPPlusRecSingle",
570                          "sigbXYZPPlusRecSingle",
571                          "NclsITSPPlusRecSingle",
572                          "NclsTPCPPlusRecSingle",
573                          "statusRecSinglePos",
574                          "ptResPiMinusRecSingle", 
575                          "etaPiMinusRecSingle",
576                          "phiPiMinusRecSingle",
577                          "bXYZPiMinusRecSingle",
578                          "sigbXYZPiMinusRecSingle",
579                          "NclsITSPiMinusRecSingle",
580                          "NclsTPCPiMinusRecSingle",
581                          "statusRecSingleNeg",
582                          "ptResPPlusRecV0", 
583                          "etaPPlusRecV0",
584                          "phiPPlusRecV0",
585                          "statusRecV0Pos",
586                          "ptResPiMinusRecV0", 
587                          "etaPiMinusRecV0",
588                          "phiPiMinusRecV0",
589                          "statusRecV0Neg"};
590
591
592   fSparseL= new THnSparseF("sparseL","sparseL",fDim,binsV0,xminV0,xmaxV0);
593   for (Int_t iaxis=0; iaxis<fDim; iaxis++){
594    fSparseL->GetAxis(iaxis)->SetName(axisNameL[iaxis]);
595    fSparseL->GetAxis(iaxis)->SetTitle(axisNameL[iaxis]);
596   }
597
598   TString  axisNameAL[37]={"ptALGeant", 
599                          "etaALGeant",
600                          "phiALGeant",
601                          "rGeant",
602                          "zGeant",
603                          "ptPiPluusGeant",  
604                          "etaPiPlusGeant",
605                          "phiPiPlusGeant",
606                          "TPCTrackLengthPiPlusGeant",
607                          "ptAPMinusGeant", 
608                          "etaAPMinusGeant",
609                          "phiAPMinusGeant",
610                          "TPCTrackLengthAPMinusGeant",
611                          "ptResPiPlusRecSingle", 
612                          "etaPiPlusRecSingle",
613                          "phiPiPlusRecSingle",
614                          "bXYZPiPlusRecSingle",
615                          "sigbXYZPiPlusRecSingle",
616                          "NclsITSPiPlusRecSingle",
617                          "NclsTPCPiPlusRecSingle",
618                          "statusRecSinglePos",
619                          "ptResAPMinusRecSingle", 
620                          "etaAPMinusRecSingle",
621                          "phiAPMinusRecSingle",
622                          "bXYZAPMinusRecSingle",
623                          "sigbXYZAPMinusRecSingle",
624                          "NclsITSAPMinusRecSingle",
625                          "NclsTPCAPMinusRecSingle",
626                          "statusRecSingleNeg",
627                          "ptResPiPlusRecV0", 
628                          "etaPiPlusRecV0",
629                          "phiPiPlusRecV0",
630                          "statusRecV0Pos",
631                          "ptResAPMinusRecV0", 
632                          "etaAPMinusRecV0",
633                          "phiAPMinusRecV0",
634                          "statusRecV0Neg"};
635
636
637   fSparseAL= new THnSparseF("sparseAL","sparseAL",fDim,binsV0,xminV0,xmaxV0);
638   for (Int_t iaxis=0; iaxis<fDim; iaxis++){
639    fSparseAL->GetAxis(iaxis)->SetName(axisNameAL[iaxis]);
640    fSparseAL->GetAxis(iaxis)->SetTitle(axisNameAL[iaxis]);
641   }
642
643   // create output container
644  
645   fOutputContainer = new TList() ;
646   fOutputContainer->SetName(GetName()) ;
647
648
649   fOutputContainer->Add(fSparseV0);
650   fOutputContainer->Add(fSparseK0);
651   fOutputContainer->Add(fSparseL);
652   fOutputContainer->Add(fSparseAL);
653   
654 }
655
656 //________________________________________________________________________
657 void AliAnalysisTaskV0QA::Exec(Option_t *) {
658
659   if (!fESD) {
660     cout<< "not a tree"<< endl;
661     return;
662   }
663
664   nEv++;
665
666
667   //Get MC data 
668   mctruth = (AliMCEventHandler*)((AliAnalysisManager::GetAnalysisManager())->GetMCtruthEventHandler());
669
670   //  Double_t vertex[3];
671   Double_t MaxVertex=150.;
672   Double_t MaxEta=1.2;
673   Double_t LineCutZRSlope=0.662486;
674   Double_t LineCutZValue=7.;
675   Int_t elecGIndex=-1;
676   Int_t posiGIndex=-1;
677   Int_t pPosLIndex=-1;
678   Int_t piNegLIndex=-1;
679
680   Int_t apNegALIndex=-1;
681   Int_t piPosALIndex=-1;
682
683   nConvGamGeant=-1;
684   nDecayK0Geant=-1;
685   nDecayLGeant=-1;
686   nDecayALGeant=-1;
687
688   for(Int_t i=0; i<gDim;i++){
689     gConvGamGeantIndex[i] = -1;
690     eNegConvGamGeantIndex[i] = -1;
691     ePosConvGamGeantIndex[i] = -1;
692     
693     eNegConvGamSingleRecIndex[i] = -1;
694     ePosConvGamSingleRecIndex[i] = -1;
695
696     eNegConvGamV0RecIndex[i] = -1;
697     ePosConvGamV0RecIndex[i] = -1;
698     ConvGamV0RecIndexPos[i] = -1;
699     ConvGamV0RecIndexNeg[i] = -1;
700
701
702     K0DecayK0GeantIndex[i] = -1;
703     piNegDecayK0GeantIndex[i] = -1;
704     piPosDecayK0GeantIndex[i] = -1;
705     
706     piNegDecayK0SingleRecIndex[i] = -1;
707     piPosDecayK0SingleRecIndex[i] = -1;
708
709     piNegDecayK0V0RecIndex[i] = -1;
710     piPosDecayK0V0RecIndex[i] = -1;
711     DecayK0V0RecIndexPos[i] = -1;
712     DecayK0V0RecIndexNeg[i] = -1;
713
714
715     lDecayLGeantIndex[i] = -1;
716     piNegDecayLGeantIndex[i] = -1;
717     pPosDecayLGeantIndex[i] = -1;
718     
719     piNegDecayLSingleRecIndex[i] = -1;
720     pPosDecayLSingleRecIndex[i] = -1;
721
722     piNegDecayLV0RecIndex[i] = -1;
723     pPosDecayLV0RecIndex[i] = -1;
724     DecayLV0RecIndexPos[i] = -1;
725     DecayLV0RecIndexNeg[i] = -1;
726
727     // Antilambda
728     alDecayALGeantIndex[i] = -1;
729     piPosDecayALGeantIndex[i] = -1;
730     apNegDecayALGeantIndex[i] = -1;
731     
732     piPosDecayALSingleRecIndex[i] = -1;
733     apNegDecayALSingleRecIndex[i] = -1;
734
735     piPosDecayALV0RecIndex[i] = -1;
736     apNegDecayALV0RecIndex[i] = -1;
737     DecayALV0RecIndexPos[i] = -1;
738     DecayALV0RecIndexNeg[i] = -1;
739
740
741   }
742
743   Int_t doMC=1;
744
745   AliKFVertex primVtx(*(fESD->GetPrimaryVertex()));
746   nTracksPrim=primVtx.GetNContributors();
747
748
749   if(mctruth && nTracksPrim>0){
750
751    stack = mctruth->MCEvent()->Stack();
752
753
754    if ( doMC){
755
756     for (Int_t iTracks = 0; iTracks < mctruth->MCEvent()->GetNumberOfTracks(); iTracks++) {
757       
758
759      TParticle* particle = stack->Particle(iTracks);
760
761
762
763      if (!particle) {
764        Printf("ERROR: Could not receive particle %d (mc loop)", iTracks);
765        continue;
766      }
767      
768      if(particle->Pt()<0.050) continue;
769      if(TMath::Abs(particle->Eta())> 1.2) continue;
770
771      
772      if (particle->GetPdgCode()== 22){
773        
774        
775        if(particle->GetMother(0) >-1 && stack->Particle(particle->GetMother(0))->GetPdgCode() == 22){
776          continue; // no photon as mothers!
777        }
778        
779        if(particle->GetMother(0) >= stack->GetNprimary()){
780          continue; // the gamma has a mother, and it is not a primary particle
781        }
782        
783        TParticle* ePos = NULL;
784        TParticle* eNeg = NULL;
785        elecGIndex=-1;
786        posiGIndex=-1;
787   
788        if(particle->GetNDaughters() >= 2){
789          for(Int_t daughterIndex=particle->GetFirstDaughter();daughterIndex<=particle->GetLastDaughter();daughterIndex++){
790            TParticle *tmpDaughter = stack->Particle(daughterIndex);
791            if(tmpDaughter->GetUniqueID() == 5){
792              if(tmpDaughter->GetPdgCode() == 11){
793                eNeg = tmpDaughter;
794                elecGIndex=daughterIndex;
795              }
796              else if(tmpDaughter->GetPdgCode() == -11){
797                ePos = tmpDaughter;
798                posiGIndex=daughterIndex;
799              }
800            }
801          }
802        }
803        
804        
805        if(ePos == NULL || eNeg == NULL){ // means we do not have two daughters from pair production
806          continue;
807        }
808        
809        if(TMath::Abs(ePos->Eta())> MaxEta || TMath::Abs(eNeg->Eta())> MaxEta){
810          continue;
811        }        
812        
813        if(ePos->R()> MaxVertex ){
814          continue; // cuts on distance from collision point
815        }
816       
817       
818        if( (TMath::Abs(ePos->Vz()) * LineCutZRSlope - LineCutZValue)  > ePos->R() ){
819          continue;               // line cut to exclude regions where we do not reconstruct
820        }                
821
822        
823        // Looking at the existance of TPC references
824
825        TParticle* ePosTPC;
826        mctruth->MCEvent()->GetParticleAndTR(posiGIndex,ePosTPC,clRefsP);
827
828        AliMCParticle *mcParticlePos = (AliMCParticle*) (mctruth->MCEvent()->GetTrack(posiGIndex));
829        if(!mcParticlePos) continue;
830
831        Int_t counter; 
832        Float_t tpcTrackLengthePos = mcParticlePos->GetTPCTrackLength(fESD->GetMagneticField(),0.05,counter,3.0); 
833
834  
835
836        int nPointsP =  clRefsP->GetEntries();
837
838        if (fRefTPC) delete fRefTPC;
839        fRefTPC = new TObjArray();
840
841        for(int iPoint=0; iPoint<nPointsP; iPoint++) {
842          AliTrackReference *ref = (AliTrackReference*)clRefsP->At(iPoint);
843          if (ref->DetectorId() == AliTrackReference::kTPC) fRefTPC->Add(new AliTrackReference(*ref));    
844        }
845
846        fRefTPC->Sort();
847
848        for(int i=0; i<fRefTPC->GetEntries(); i++) {
849          AliTrackReference *ref = (AliTrackReference*)(*fRefTPC)[i]; 
850          fLabelsTPC[i] = ref->GetTrack();
851        }
852        int labelPosRefs=0;
853        for(int iPoint=GetTPCReference(iTracks);iPoint<fRefTPC->GetEntries();iPoint++){
854          AliTrackReference* aRef = (AliTrackReference*)(*fRefTPC)[iPoint];
855          if (aRef->GetTrack() != posiGIndex ) break;                    
856          ++labelPosRefs;
857        }
858      
859
860
861
862        TParticle* eNegTPC;
863        mctruth->MCEvent()->GetParticleAndTR(elecGIndex,eNegTPC,clRefsN);
864
865        AliMCParticle *mcParticleNeg = (AliMCParticle*) (mctruth->MCEvent()->GetTrack(elecGIndex));
866        if(!mcParticleNeg) continue;
867
868        Int_t counterN; 
869        Float_t tpcTrackLengtheNeg = mcParticleNeg->GetTPCTrackLength(fESD->GetMagneticField(),0.05,counterN,3.0); 
870        int nPointsN =  clRefsN->GetEntries();
871
872        if (fRefTPC) delete fRefTPC;
873        fRefTPC = new TObjArray();
874
875        for(int iPoint=0; iPoint<nPointsN; iPoint++) {
876          AliTrackReference *ref = (AliTrackReference*)clRefsN->At(iPoint);
877          if (ref->DetectorId() == AliTrackReference::kTPC) fRefTPC->Add(new AliTrackReference(*ref));    
878        }
879
880        fRefTPC->Sort();
881
882        for(int i=0; i<fRefTPC->GetEntries(); i++) {
883          AliTrackReference *ref = (AliTrackReference*)(*fRefTPC)[i]; 
884          fLabelsTPC[i] = ref->GetTrack();
885        }
886        int labelNegRefs=0;
887        for(int iPoint=GetTPCReference(iTracks);iPoint<fRefTPC->GetEntries();iPoint++){
888          AliTrackReference* aRef = (AliTrackReference*)(*fRefTPC)[iPoint];
889          if (aRef->GetTrack() != elecGIndex ) break;                    
890          ++labelNegRefs;
891        }
892
893
894        
895        if ( labelNegRefs==0 || labelPosRefs==0) continue; // if e+/e- do not have a TPC ref continue;
896        ////////////////////////////////////////////////////////////////////
897        
898        
899        nConvGamGeant++;
900        gConvGamGeantIndex[nConvGamGeant]=iTracks;
901        eNegConvGamGeantIndex[nConvGamGeant] = elecGIndex;
902        ePosConvGamGeantIndex[nConvGamGeant] = posiGIndex;
903        
904        eNegConvGamGeantLength[nConvGamGeant] = tpcTrackLengtheNeg;
905        ePosConvGamGeantLength[nConvGamGeant] = tpcTrackLengthePos;
906
907      }
908
909      
910      TParticle* piPos = NULL;
911      TParticle* piNeg = NULL;
912      piPosK0Index=-1;
913      piNegK0Index=-1;
914
915      if (particle->GetPdgCode()== 310){          // k0short
916        if(particle->GetNDaughters() == 2){
917          for(Int_t daughterIndex=particle->GetFirstDaughter();daughterIndex<=particle->GetLastDaughter();daughterIndex++){
918            TParticle *tmpDaughter = stack->Particle(daughterIndex);
919            if(tmpDaughter->GetPdgCode() == 211){
920              piPos= tmpDaughter;
921              piPosK0Index=daughterIndex;
922            }
923            else if(tmpDaughter->GetPdgCode() == -211){
924              piNeg = tmpDaughter;
925              piNegK0Index=daughterIndex;
926            }
927          }
928        }
929
930        if(piPos == NULL || piNeg == NULL){ // means we do not have two daughters from K0short decay
931          continue;
932        }
933        
934        if(TMath::Abs(piPos->Eta())> MaxEta || TMath::Abs(piNeg->Eta())> MaxEta){
935          continue;
936        }        
937        
938        if(piPos->R()> MaxVertex ){
939          continue; // cuts on distance from collision point
940        }
941       
942       
943        if( (TMath::Abs(piPos->Vz()) * LineCutZRSlope - LineCutZValue)  > piPos->R() ){
944          continue;               // line cut to exclude regions where we do not reconstruct
945        }                
946
947       // Looking at the existance of TPC references
948
949        TParticle* ePosTPC;
950        mctruth->MCEvent()->GetParticleAndTR(piPosK0Index,ePosTPC,clRefsP);
951
952        AliMCParticle *mcParticlePos = (AliMCParticle*) (mctruth->MCEvent()->GetTrack(piPosK0Index));
953        if(!mcParticlePos) continue;
954
955        Int_t counter; 
956        Float_t tpcTrackLengthePos = mcParticlePos->GetTPCTrackLength(fESD->GetMagneticField(),0.05,counter,3.0); 
957  
958
959        int nPointsP =  clRefsP->GetEntries();
960        if (fRefTPC) delete fRefTPC;
961        fRefTPC = new TObjArray();
962
963        for(int iPoint=0; iPoint<nPointsP; iPoint++) {
964          AliTrackReference *ref = (AliTrackReference*)clRefsP->At(iPoint);
965          if (ref->DetectorId() == AliTrackReference::kTPC) fRefTPC->Add(new AliTrackReference(*ref));    
966        }
967
968        fRefTPC->Sort();
969
970        for(int i=0; i<fRefTPC->GetEntries(); i++) {
971          AliTrackReference *ref = (AliTrackReference*)(*fRefTPC)[i]; 
972          fLabelsTPC[i] = ref->GetTrack();
973        }
974        int labelPosRefs=0;
975        for(int iPoint=GetTPCReference(iTracks);iPoint<fRefTPC->GetEntries();iPoint++){
976          AliTrackReference* aRef = (AliTrackReference*)(*fRefTPC)[iPoint];
977          if (aRef->GetTrack() != piPosK0Index ) break;                  
978          ++labelPosRefs;
979        }
980      
981
982
983
984        TParticle* eNegTPC;
985        mctruth->MCEvent()->GetParticleAndTR(piNegK0Index,eNegTPC,clRefsN);
986
987        AliMCParticle *mcParticleNeg = (AliMCParticle*) (mctruth->MCEvent()->GetTrack(piNegK0Index));
988        if(!mcParticleNeg) continue;
989
990        Int_t counterN; 
991        Float_t tpcTrackLengtheNeg = mcParticleNeg->GetTPCTrackLength(fESD->GetMagneticField(),0.05,counterN,3.0); 
992
993        int nPointsN =  clRefsN->GetEntries();
994        if (fRefTPC) delete fRefTPC;
995        fRefTPC = new TObjArray();
996
997        for(int iPoint=0; iPoint<nPointsN; iPoint++) {
998          AliTrackReference *ref = (AliTrackReference*)clRefsN->At(iPoint);
999          if (ref->DetectorId() == AliTrackReference::kTPC) fRefTPC->Add(new AliTrackReference(*ref));    
1000        }
1001
1002        fRefTPC->Sort();
1003
1004        for(int i=0; i<fRefTPC->GetEntries(); i++) {
1005          AliTrackReference *ref = (AliTrackReference*)(*fRefTPC)[i]; 
1006          fLabelsTPC[i] = ref->GetTrack();
1007        }
1008        int labelNegRefs=0;
1009        for(int iPoint=GetTPCReference(iTracks);iPoint<fRefTPC->GetEntries();iPoint++){
1010          AliTrackReference* aRef = (AliTrackReference*)(*fRefTPC)[iPoint];
1011          if (aRef->GetTrack() != piNegK0Index ) break;                  
1012          ++labelNegRefs;
1013        }
1014        
1015        if ( labelNegRefs==0 || labelPosRefs==0) continue; // if pi+/pi- do not have a TPC ref continue;
1016        ////////////////////////////////////////////////////////////////////
1017        
1018        nDecayK0Geant++;
1019
1020        K0DecayK0GeantIndex[nDecayK0Geant]=iTracks;
1021        piNegDecayK0GeantIndex[nDecayK0Geant]=piNegK0Index;
1022        piPosDecayK0GeantIndex[nDecayK0Geant]=piPosK0Index;
1023        piNegDecayK0GeantLength[nDecayK0Geant]=tpcTrackLengtheNeg;
1024        piPosDecayK0GeantLength[nDecayK0Geant]=tpcTrackLengthePos;
1025       
1026      }    
1027
1028
1029      TParticle* pPos = NULL;
1030      TParticle* piNegL = NULL;
1031      pPosLIndex=-1;
1032      piNegLIndex=-1;
1033
1034  
1035      if (particle->GetPdgCode()== 3122){        //lambda
1036       
1037        if(particle->GetNDaughters() == 2){
1038          for(Int_t daughterIndex=particle->GetFirstDaughter();daughterIndex<=particle->GetLastDaughter();daughterIndex++){
1039            TParticle *tmpDaughter = stack->Particle(daughterIndex);
1040            if(tmpDaughter->GetPdgCode() == 2212){
1041              pPos= tmpDaughter;
1042              pPosLIndex=daughterIndex;
1043            }
1044            else if(tmpDaughter->GetPdgCode() == -211){
1045              piNegL = tmpDaughter;
1046              piNegLIndex=daughterIndex;
1047            }
1048          }
1049        }
1050
1051        if(pPos == NULL || piNegL == NULL){ // means we do not have two daughters from lambda decay
1052          continue;
1053        }
1054
1055        if(TMath::Abs(pPos->Eta())> MaxEta || TMath::Abs(piNegL->Eta())> MaxEta){
1056          continue;
1057        }        
1058        
1059        if(pPos->R()> MaxVertex ){
1060          continue; // cuts on distance from collision point
1061        }
1062       
1063       
1064        if( (TMath::Abs(pPos->Vz()) * LineCutZRSlope - LineCutZValue)  > pPos->R() ){
1065          continue;               // line cut to exclude regions where we do not reconstruct
1066        }
1067
1068
1069      // Looking at the existance of TPC references
1070
1071        TParticle* ePosTPC;
1072        mctruth->MCEvent()->GetParticleAndTR(pPosLIndex,ePosTPC,clRefsP);
1073
1074        AliMCParticle *mcParticlePos = (AliMCParticle*) (mctruth->MCEvent()->GetTrack(pPosLIndex));
1075        if(!mcParticlePos) continue;
1076
1077        Int_t counter; 
1078        Float_t tpcTrackLengthePos = mcParticlePos->GetTPCTrackLength(fESD->GetMagneticField(),0.05,counter,3.0); 
1079  
1080
1081        int nPointsP =  clRefsP->GetEntries();
1082        if (fRefTPC) delete fRefTPC;
1083        fRefTPC = new TObjArray();
1084
1085        for(int iPoint=0; iPoint<nPointsP; iPoint++) {
1086          AliTrackReference *ref = (AliTrackReference*)clRefsP->At(iPoint);
1087          if (ref->DetectorId() == AliTrackReference::kTPC) fRefTPC->Add(new AliTrackReference(*ref));    
1088        }
1089
1090        fRefTPC->Sort();
1091
1092        for(int i=0; i<fRefTPC->GetEntries(); i++) {
1093          AliTrackReference *ref = (AliTrackReference*)(*fRefTPC)[i]; 
1094          fLabelsTPC[i] = ref->GetTrack();
1095        }
1096        int labelPosRefs=0;
1097        for(int iPoint=GetTPCReference(iTracks);iPoint<fRefTPC->GetEntries();iPoint++){
1098          AliTrackReference* aRef = (AliTrackReference*)(*fRefTPC)[iPoint];
1099          if (aRef->GetTrack() != pPosLIndex ) break;                    
1100          ++labelPosRefs;
1101        }
1102      
1103
1104
1105
1106        TParticle* eNegTPC;
1107        mctruth->MCEvent()->GetParticleAndTR(piNegLIndex,eNegTPC,clRefsN);
1108
1109        AliMCParticle *mcParticleNeg = (AliMCParticle*) (mctruth->MCEvent()->GetTrack(piNegLIndex));
1110        if(!mcParticleNeg) continue;
1111
1112        Int_t counterN; 
1113        Float_t tpcTrackLengtheNeg = mcParticleNeg->GetTPCTrackLength(fESD->GetMagneticField(),0.05,counterN,3.0); 
1114
1115        int nPointsN =  clRefsN->GetEntries();
1116        if (fRefTPC) delete fRefTPC;
1117        fRefTPC = new TObjArray();
1118
1119        for(int iPoint=0; iPoint<nPointsN; iPoint++) {
1120          AliTrackReference *ref = (AliTrackReference*)clRefsN->At(iPoint);
1121          if (ref->DetectorId() == AliTrackReference::kTPC) fRefTPC->Add(new AliTrackReference(*ref));    
1122        }
1123
1124        fRefTPC->Sort();
1125
1126        for(int i=0; i<fRefTPC->GetEntries(); i++) {
1127          AliTrackReference *ref = (AliTrackReference*)(*fRefTPC)[i]; 
1128          fLabelsTPC[i] = ref->GetTrack();
1129        }
1130        int labelNegRefs=0;
1131        for(int iPoint=GetTPCReference(iTracks);iPoint<fRefTPC->GetEntries();iPoint++){
1132          AliTrackReference* aRef = (AliTrackReference*)(*fRefTPC)[iPoint];
1133          if (aRef->GetTrack() != piNegLIndex ) break;                   
1134          ++labelNegRefs;
1135        }
1136        
1137        if ( labelNegRefs==0 || labelPosRefs==0) continue; // if proton/pi- do not have a TPC ref continue;
1138        ////////////////////////////////////////////////////////////////////
1139      
1140        nDecayLGeant++;
1141
1142        lDecayLGeantIndex[nDecayLGeant]=iTracks;
1143      
1144        piNegDecayLGeantIndex[nDecayLGeant]=piNegLIndex;
1145        pPosDecayLGeantIndex[nDecayLGeant]=pPosLIndex;
1146
1147        piNegDecayLGeantLength[nDecayLGeant]=tpcTrackLengtheNeg;
1148        pPosDecayLGeantLength[nDecayLGeant]=tpcTrackLengthePos;
1149
1150      
1151      }
1152
1153      /////////////////////////////////////////////////////
1154
1155      // AntiLambda    
1156      TParticle* apNeg = NULL;
1157      TParticle* piPosAL = NULL;
1158      apNegALIndex=-1;
1159      piPosALIndex=-1;
1160
1161  
1162      if (particle->GetPdgCode()== -3122){        //antilambda
1163
1164        if(particle->GetNDaughters() == 2){
1165          for(Int_t daughterIndex=particle->GetFirstDaughter();daughterIndex<=particle->GetLastDaughter();daughterIndex++){
1166            TParticle *tmpDaughter = stack->Particle(daughterIndex);
1167            if(tmpDaughter->GetPdgCode() == -2212){
1168              apNeg= tmpDaughter;
1169              apNegALIndex=daughterIndex;
1170            }
1171            else if(tmpDaughter->GetPdgCode() == 211){
1172              piPosAL = tmpDaughter;
1173              piPosALIndex=daughterIndex;
1174            }
1175          }
1176        }
1177
1178        if(apNeg == NULL || piPosAL == NULL){ // means we do not have two daughters from antilambda decay
1179          continue;
1180        }
1181
1182        if(TMath::Abs(apNeg->Eta())> MaxEta || TMath::Abs(piPosAL->Eta())> MaxEta){
1183          continue;
1184        }        
1185        
1186        if(apNeg->R()> MaxVertex ){
1187          continue; // cuts on distance from collision point
1188        }
1189       
1190       
1191        if( (TMath::Abs(apNeg->Vz()) * LineCutZRSlope - LineCutZValue)  > apNeg->R() ){
1192          continue;               // line cut to exclude regions where we do not reconstruct
1193        }
1194
1195
1196      // Looking at the existance of TPC references
1197
1198        TParticle* ePosTPC;
1199        mctruth->MCEvent()->GetParticleAndTR(piPosALIndex,ePosTPC,clRefsP);
1200
1201        AliMCParticle *mcParticlePos = (AliMCParticle*) (mctruth->MCEvent()->GetTrack(piPosALIndex));
1202        if(!mcParticlePos) continue;
1203
1204        Int_t counter; 
1205        Float_t tpcTrackLengthePos = mcParticlePos->GetTPCTrackLength(fESD->GetMagneticField(),0.05,counter,3.0); 
1206
1207        int nPointsP =  clRefsP->GetEntries();
1208        if (fRefTPC) delete fRefTPC;
1209        fRefTPC = new TObjArray();
1210
1211        for(int iPoint=0; iPoint<nPointsP; iPoint++) {
1212          AliTrackReference *ref = (AliTrackReference*)clRefsP->At(iPoint);
1213          if (ref->DetectorId() == AliTrackReference::kTPC) fRefTPC->Add(new AliTrackReference(*ref));    
1214        }
1215
1216        fRefTPC->Sort();
1217
1218        for(int i=0; i<fRefTPC->GetEntries(); i++) {
1219          AliTrackReference *ref = (AliTrackReference*)(*fRefTPC)[i]; 
1220          fLabelsTPC[i] = ref->GetTrack();
1221        }
1222        int labelPosRefs=0;
1223        for(int iPoint=GetTPCReference(iTracks);iPoint<fRefTPC->GetEntries();iPoint++){
1224          AliTrackReference* aRef = (AliTrackReference*)(*fRefTPC)[iPoint];
1225          if (aRef->GetTrack() != piPosALIndex ) break;                  
1226          ++labelPosRefs;
1227        }
1228      
1229
1230        TParticle* eNegTPC;
1231        mctruth->MCEvent()->GetParticleAndTR(apNegALIndex,eNegTPC,clRefsN);
1232
1233        AliMCParticle *mcParticleNeg = (AliMCParticle*) (mctruth->MCEvent()->GetTrack(apNegALIndex));
1234        if(!mcParticleNeg) continue;
1235
1236        Int_t counterN; 
1237        Float_t tpcTrackLengtheNeg = mcParticleNeg->GetTPCTrackLength(fESD->GetMagneticField(),0.05,counterN,3.0); 
1238
1239        int nPointsN =  clRefsN->GetEntries();
1240        if (fRefTPC) delete fRefTPC;
1241        fRefTPC = new TObjArray();
1242
1243        for(int iPoint=0; iPoint<nPointsN; iPoint++) {
1244          AliTrackReference *ref = (AliTrackReference*)clRefsN->At(iPoint);
1245          if (ref->DetectorId() == AliTrackReference::kTPC) fRefTPC->Add(new AliTrackReference(*ref));    
1246        }
1247
1248        fRefTPC->Sort();
1249
1250        for(int i=0; i<fRefTPC->GetEntries(); i++) {
1251          AliTrackReference *ref = (AliTrackReference*)(*fRefTPC)[i]; 
1252          fLabelsTPC[i] = ref->GetTrack();
1253        }
1254        int labelNegRefs=0;
1255        for(int iPoint=GetTPCReference(iTracks);iPoint<fRefTPC->GetEntries();iPoint++){
1256          AliTrackReference* aRef = (AliTrackReference*)(*fRefTPC)[iPoint];
1257          if (aRef->GetTrack() != apNegALIndex ) break;                  
1258          ++labelNegRefs;
1259        }
1260
1261
1262        
1263        if ( labelNegRefs==0 || labelPosRefs==0) continue; // if proton/pi- do not have a TPC ref continue;
1264        ////////////////////////////////////////////////////////////////////
1265      
1266        nDecayALGeant++;
1267        alDecayALGeantIndex[nDecayALGeant]=iTracks;
1268      
1269        piPosDecayALGeantIndex[nDecayALGeant]=piPosALIndex;
1270        apNegDecayALGeantIndex[nDecayALGeant]=apNegALIndex;
1271
1272        piPosDecayALGeantLength[nDecayALGeant]=tpcTrackLengthePos;
1273        apNegDecayALGeantLength[nDecayALGeant]=tpcTrackLengtheNeg;
1274
1275      
1276      }  // AntiLambda    
1277
1278     } //track loop 
1279    }    
1280
1281   }
1282
1283
1284   AliKFParticle::SetField(fESD->GetMagneticField());
1285
1286   const AliESDVertex *pvertex = fESD->GetPrimaryVertex();
1287   Double_t xyzVtx[3];
1288   pvertex->GetXYZ(xyzVtx);
1289
1290   if(nTracksPrim>0) {
1291
1292     InspectListOfChargedParticles();
1293     InspectListOfV0s();
1294
1295
1296     if(nConvGamGeant>-1){
1297       FillHnSparseGamma();
1298     }
1299
1300     if(nDecayK0Geant>-1){
1301       FillHnSparseK0();
1302     }
1303
1304     if(nDecayLGeant>-1){
1305       FillHnSparseL();
1306     }
1307
1308     if(nDecayALGeant>-1){
1309       FillHnSparseAL();
1310     }
1311
1312   }
1313
1314  
1315   PostData(0,fOutputContainer );
1316   
1317
1318 }
1319
1320 void AliAnalysisTaskV0QA::Terminate(Option_t *) {
1321   // Draw some histogram at the end.
1322
1323 }
1324
1325
1326 Int_t AliAnalysisTaskV0QA::GetTPCReference(Int_t label) {
1327
1328   int start = TMath::BinarySearch(fRefTPC->GetEntries(), fLabelsTPC, label);
1329
1330   while (start >= 0) {
1331     AliTrackReference *ref = (AliTrackReference*)(*fRefTPC)[start];
1332     if (ref->GetTrack() != label) return start+1;
1333     start--;
1334   }
1335
1336   return 0;
1337 }
1338
1339
1340
1341
1342
1343 void AliAnalysisTaskV0QA::InspectListOfChargedParticles(){
1344
1345
1346   for(Int_t iTracks = 0; iTracks < fESD->GetNumberOfTracks(); iTracks++){
1347
1348     AliESDtrack* curTrack = fESD->GetTrack(iTracks);
1349
1350     if(!curTrack){
1351       continue;
1352     }
1353
1354
1355 //     if( !(curTrack->GetStatus() & AliESDtrack::kTPCrefit)){
1356 //       continue;
1357 //     }
1358
1359
1360     Int_t labelMC = TMath::Abs(curTrack->GetLabel());
1361
1362     if ( labelMC > stack->GetNtrack() ) continue;
1363
1364
1365     TParticle* curParticle = stack->Particle(labelMC);
1366     if(curParticle->GetMother(0)==-1){
1367       continue;
1368     }
1369
1370      
1371     if(TMath::Abs(curParticle->GetPdgCode()) == 11){ // e+/e-
1372       
1373       if( stack->Particle(curParticle->GetMother(0))->GetPdgCode()==22 ){  // e+/e- from gamma
1374         if( curParticle->GetUniqueID()!=5 ){ // e+/e- from gamma conversion
1375           continue;
1376         }
1377
1378         for(Int_t iGamConv=0;iGamConv<nConvGamGeant+1;iGamConv++ ){
1379           if(curTrack->GetSign()>0){
1380             if (labelMC== ePosConvGamGeantIndex[iGamConv]){
1381               ePosConvGamSingleRecIndex[iGamConv]=iTracks;
1382             }
1383           }else{
1384             if (labelMC== eNegConvGamGeantIndex[iGamConv]){
1385               eNegConvGamSingleRecIndex[iGamConv]=iTracks;
1386             }
1387           }
1388         } // loop over geant converted gammas
1389
1390       }
1391     } // condition to select reconstructed electrons
1392   
1393
1394
1395     if(TMath::Abs(curParticle->GetPdgCode()) == 211 || TMath::Abs(curParticle->GetPdgCode())==2212 ){ // pi+/pi-
1396       
1397       if( stack->Particle(curParticle->GetMother(0))->GetPdgCode()==310 || 
1398           stack->Particle(curParticle->GetMother(0))->GetPdgCode()==3122 || 
1399           stack->Particle(curParticle->GetMother(0))->GetPdgCode()==-3122 ){  // pi+/proton/pi- from K0/Lambda
1400
1401         for(Int_t iK0Dec=0;iK0Dec<nDecayK0Geant+1;iK0Dec++ ){
1402           if(curTrack->GetSign()>0){
1403             if (labelMC== piPosDecayK0GeantIndex[iK0Dec]){
1404               piPosDecayK0SingleRecIndex[iK0Dec]=iTracks;
1405             }
1406           }else{
1407             if (labelMC== piNegDecayK0GeantIndex[iK0Dec]){
1408               piNegDecayK0SingleRecIndex[iK0Dec]=iTracks;
1409             }
1410           }
1411         } // loop over geant decay K0
1412
1413         for(Int_t iLDec=0;iLDec<nDecayLGeant+1;iLDec++ ){
1414           if(curTrack->GetSign()>0){
1415             if (labelMC== pPosDecayLGeantIndex[iLDec]){
1416               pPosDecayLSingleRecIndex[iLDec]=iTracks;
1417             }
1418           }else{
1419             if (labelMC== piNegDecayLGeantIndex[iLDec]){
1420               piNegDecayLSingleRecIndex[iLDec]=iTracks;
1421             }
1422           }
1423         } // loop over geant decay Lambda
1424
1425         for(Int_t iALDec=0;iALDec<nDecayALGeant+1;iALDec++ ){
1426           if(curTrack->GetSign()<0){
1427             if (labelMC== apNegDecayALGeantIndex[iALDec]){
1428               apNegDecayALSingleRecIndex[iALDec]=iTracks;
1429             }
1430           }else{
1431             if (labelMC== piPosDecayALGeantIndex[iALDec]){
1432               piPosDecayALSingleRecIndex[iALDec]=iTracks;
1433             }
1434           }
1435         } // loop over geant decay antiLambda
1436       }
1437     } // condition to select reconstructed electrons
1438    } // all reconstructed track
1439
1440
1441 }
1442
1443 void AliAnalysisTaskV0QA::InspectListOfV0s(){
1444
1445   AliESDtrack* trackPos= NULL;
1446   AliESDtrack* trackNeg= NULL;
1447   Int_t grandMotherPos=-1;
1448   Int_t grandMotherNeg=-1;
1449   Int_t motherPos=-1;
1450   Int_t motherNeg=-1;
1451   Int_t pIndex=-1;
1452   Int_t nIndex=-1;
1453
1454   for(Int_t iV0MI = 0; iV0MI < fESD->GetNumberOfV0s(); iV0MI++) {
1455
1456     AliESDv0 * fV0MIs = fESD->GetV0(iV0MI);
1457     
1458
1459     if ( !fV0MIs->GetOnFlyStatus() ){
1460         continue;
1461     }
1462     if(nTracksPrim<=0) {
1463       continue;
1464     }
1465
1466
1467     AliESDtrack* trackPosTest = fESD->GetTrack(fV0MIs->GetPindex());
1468     AliESDtrack* trackNegTest = fESD->GetTrack(fV0MIs->GetNindex());
1469
1470
1471     if ( trackPosTest->GetSign() == trackNegTest->GetSign()){
1472      continue;
1473     }
1474
1475     //To avoid ghosts
1476
1477 //     if( !(trackPosTest->GetStatus() & AliESDtrack::kTPCrefit)){
1478 //       continue;
1479 //     }
1480
1481 //     if( !(trackNegTest->GetStatus() & AliESDtrack::kTPCrefit)){
1482 //       continue;
1483 //     }
1484
1485     if( trackPosTest->GetSign() ==1){
1486       trackPos =fESD->GetTrack(fV0MIs->GetPindex());
1487       trackNeg =fESD->GetTrack(fV0MIs->GetNindex());
1488       pIndex=fV0MIs->GetPindex();
1489       nIndex=fV0MIs->GetNindex();
1490     }
1491
1492     if( trackPosTest->GetSign() ==-1){
1493       trackPos =fESD->GetTrack(fV0MIs->GetNindex());
1494       trackNeg =fESD->GetTrack(fV0MIs->GetPindex());
1495       pIndex=fV0MIs->GetNindex();
1496       nIndex=fV0MIs->GetPindex();
1497
1498     }
1499
1500     Int_t labelNeg=TMath::Abs(trackNeg->GetLabel());
1501     if(labelNeg > stack->GetNtrack() ) continue;
1502     TParticle * particleNeg= stack->Particle(labelNeg);
1503
1504     Int_t labelPos=TMath::Abs(trackPos->GetLabel());
1505     if(labelPos > stack->GetNtrack() ) continue;
1506     TParticle * particlePos= stack->Particle(labelPos);
1507
1508
1509     if(particlePos->GetMother(0)>-1){
1510       grandMotherPos=stack->Particle(particlePos->GetMother(0))->GetMother(0);
1511       motherPos=particlePos->GetMother(0);
1512     }
1513     
1514     if(particleNeg->GetMother(0)>-1){
1515       grandMotherNeg=stack->Particle(particleNeg->GetMother(0))->GetMother(0);
1516       motherNeg=particleNeg->GetMother(0);
1517     }
1518
1519     if(motherPos == motherNeg &&  motherPos!=-1 ){
1520       if( particlePos->GetPdgCode() ==-11  &&   particleNeg->GetPdgCode()==11 ){
1521         for(Int_t iGamConv=0;iGamConv<nConvGamGeant+1;iGamConv++ ){
1522           if (labelPos== ePosConvGamGeantIndex[iGamConv]){
1523             ePosConvGamV0RecIndex[iGamConv]=pIndex;
1524             ConvGamV0RecIndexPos[iGamConv]=iV0MI;
1525           }
1526           if (labelNeg== eNegConvGamGeantIndex[iGamConv]){
1527             eNegConvGamV0RecIndex[iGamConv]=nIndex;
1528             ConvGamV0RecIndexNeg[iGamConv]=iV0MI;
1529           }
1530
1531         } // loop over geant converted gammas
1532       }
1533
1534       if( particlePos->GetPdgCode()==211  &&   particleNeg->GetPdgCode()==-211 ){
1535         for(Int_t iK0Dec=0;iK0Dec<nDecayK0Geant+1;iK0Dec++ ){
1536           if (labelPos== piPosDecayK0GeantIndex[iK0Dec]){
1537             piPosDecayK0V0RecIndex[iK0Dec]=pIndex;
1538             DecayK0V0RecIndexPos[iK0Dec]=iV0MI;
1539           }
1540           if (labelNeg== piNegDecayK0GeantIndex[iK0Dec]){
1541             piNegDecayK0V0RecIndex[iK0Dec]=nIndex;
1542             DecayK0V0RecIndexNeg[iK0Dec]=iV0MI;
1543           }
1544
1545         } // loop over geant K0
1546       }
1547
1548       if( particlePos->GetPdgCode()==2212  && particleNeg->GetPdgCode()==-211 ){
1549         for(Int_t iLDec=0;iLDec<nDecayLGeant+1;iLDec++ ){
1550           if (labelPos== pPosDecayLGeantIndex[iLDec]){
1551             pPosDecayLV0RecIndex[iLDec]=pIndex;
1552             DecayLV0RecIndexPos[iLDec]=iV0MI;
1553           }
1554           if (labelNeg== piNegDecayLGeantIndex[iLDec]){
1555             piNegDecayLV0RecIndex[iLDec]=nIndex;
1556             DecayLV0RecIndexNeg[iLDec]=iV0MI;
1557           }
1558
1559         } // loop over geant Lambda
1560       }
1561
1562       if( particleNeg->GetPdgCode()==-2212  && particlePos->GetPdgCode()==211 ){
1563         for(Int_t iALDec=0;iALDec<nDecayALGeant+1;iALDec++ ){
1564           if (labelNeg== apNegDecayALGeantIndex[iALDec]){
1565             apNegDecayALV0RecIndex[iALDec]=nIndex;
1566             DecayALV0RecIndexNeg[iALDec]=iV0MI;
1567           }
1568           if (labelPos== piPosDecayALGeantIndex[iALDec]){
1569             piPosDecayALV0RecIndex[iALDec]=pIndex;
1570             DecayALV0RecIndexPos[iALDec]=iV0MI;
1571           }
1572
1573         } // loop over geant antiLambda
1574       }
1575
1576
1577     }
1578     
1579   }
1580   for(Int_t iGamConv=0;iGamConv<nConvGamGeant+1;iGamConv++ ){
1581     if ( ConvGamV0RecIndexNeg[iGamConv]!=  ConvGamV0RecIndexPos[iGamConv]){
1582       ePosConvGamV0RecIndex[iGamConv]=-1;
1583       eNegConvGamV0RecIndex[iGamConv]=-1; 
1584       ConvGamV0RecIndexNeg[iGamConv]=-1;
1585       ConvGamV0RecIndexPos[iGamConv]=-1;
1586
1587     }
1588   }
1589
1590   for(Int_t iLDec=0;iLDec<nDecayLGeant+1;iLDec++ ){
1591     if(DecayLV0RecIndexPos[iLDec] !=  DecayLV0RecIndexNeg[iLDec]){
1592       piNegDecayLV0RecIndex[iLDec]=-1;
1593       pPosDecayLV0RecIndex[iLDec]=-1;
1594       DecayLV0RecIndexNeg[iLDec]=-1;
1595       DecayLV0RecIndexPos[iLDec]=-1;
1596     }
1597   }
1598
1599   for(Int_t iALDec=0;iALDec<nDecayALGeant+1;iALDec++ ){
1600     if(DecayALV0RecIndexPos[iALDec] !=  DecayALV0RecIndexNeg[iALDec]){
1601       piPosDecayALV0RecIndex[iALDec]=-1;
1602       apNegDecayALV0RecIndex[iALDec]=-1;
1603       DecayALV0RecIndexNeg[iALDec]=-1;
1604       DecayALV0RecIndexPos[iALDec]=-1;
1605     }
1606   }
1607
1608   for(Int_t iK0Dec=0;iK0Dec<nDecayK0Geant+1;iK0Dec++ ){
1609     if(DecayK0V0RecIndexPos[iK0Dec] !=  DecayK0V0RecIndexNeg[iK0Dec]){
1610       piNegDecayK0V0RecIndex[iK0Dec]=-1;
1611       piPosDecayK0V0RecIndex[iK0Dec]=-1;
1612       DecayK0V0RecIndexNeg[iK0Dec]=-1;
1613       DecayK0V0RecIndexPos[iK0Dec]=-1;
1614     }
1615   }
1616   
1617
1618 }
1619 void AliAnalysisTaskV0QA::FillHnSparseGamma()
1620 {
1621   Double_t massE=0.00051099892;
1622   Double_t ppSgl[3];
1623   Double_t pmSgl[3];
1624   Float_t bPosSgl[2];
1625   Float_t bNegSgl[2];
1626   Float_t bPosCov[3];
1627   Float_t bNegCov[3];
1628   
1629   Double_t ppV0[3];
1630   Double_t pmV0[3];
1631   Double_t xrG[3];
1632   
1633   TLorentzVector posSglTrack;
1634   TLorentzVector negSglTrack;
1635   Double_t posPt,posEta,posPhi;
1636   Double_t negPt,negEta,negPhi;
1637
1638   TLorentzVector posV0Track;
1639   TLorentzVector negV0Track;
1640   Double_t posV0Pt,posV0Eta,posV0Phi;
1641   Double_t negV0Pt,negV0Eta,negV0Phi;
1642   
1643   Float_t nClsITSPos=-1;
1644   Float_t nClsITSNeg=-1;
1645
1646   Float_t nClsTPCPos=-1;
1647   Float_t nClsTPCNeg=-1;
1648
1649   Int_t statusSingPos=-1;
1650   Int_t statusSingNeg=-1;
1651
1652   Int_t statusV0Pos=-1;
1653   Int_t statusV0Neg=-1;
1654
1655
1656   for(Int_t i=0;i<nConvGamGeant+1;i++){
1657     TParticle* gamPart = stack->Particle(gConvGamGeantIndex[i]);
1658     TParticle* ePosPart = stack->Particle(ePosConvGamGeantIndex[i]);
1659     TParticle* eNegPart = stack->Particle(eNegConvGamGeantIndex[i]);
1660     if (ePosConvGamSingleRecIndex[i]!=-1){
1661       AliESDtrack * ePosSglTrack = fESD->GetTrack(ePosConvGamSingleRecIndex[i]);
1662       ePosSglTrack->GetPxPyPz(ppSgl); 
1663       posSglTrack.SetXYZM(ppSgl[0],ppSgl[1],ppSgl[2],massE);
1664       posPt  = posSglTrack.Pt();
1665       posEta = posSglTrack.Eta();
1666       posPhi = posSglTrack.Phi();
1667       ePosSglTrack->GetImpactParameters(bPosSgl,bPosCov);
1668       nClsITSPos=ePosSglTrack->GetNcls(0);
1669       nClsTPCPos=ePosSglTrack->GetNcls(1);
1670       statusSingPos=1;
1671     }else{
1672       posPt  = 1000000;
1673       posEta = -2.;
1674       posPhi = -2*TMath::Pi();
1675       bPosSgl[0]=-100.;
1676       bPosSgl[1]=-100.;
1677       bPosCov[0]=-100;
1678       bPosCov[2]=-100;
1679       nClsITSPos=-1;
1680       nClsTPCPos=-1;
1681       statusSingPos=-1;
1682     }
1683     
1684     if (eNegConvGamSingleRecIndex[i]!=-1){
1685       AliESDtrack * eNegSglTrack = fESD->GetTrack(eNegConvGamSingleRecIndex[i]);
1686       eNegSglTrack->GetPxPyPz(pmSgl); 
1687       negSglTrack.SetXYZM(pmSgl[0],pmSgl[1],pmSgl[2],massE);
1688       negPt  = negSglTrack.Pt();
1689       negEta = negSglTrack.Eta();
1690       negPhi = negSglTrack.Phi();
1691       eNegSglTrack->GetImpactParameters(bNegSgl,bNegCov);
1692       nClsITSNeg=eNegSglTrack->GetNcls(0);
1693       nClsTPCNeg=eNegSglTrack->GetNcls(1);
1694       statusSingNeg=1;
1695     }else{
1696       negPt  = 1000000;
1697       negEta = -2.;
1698       negPhi = -2*TMath::Pi();
1699       bNegSgl[0]=-100.;
1700       bNegSgl[1]=-100.;
1701       bNegCov[0]=-100;
1702       bNegCov[2]=-100;
1703       nClsITSNeg=-1;
1704       nClsTPCNeg=-1;
1705       statusSingNeg=-1;
1706     }
1707
1708     posV0Pt  = 1000000;
1709     posV0Eta = -2.;
1710     posV0Phi = -2*TMath::Pi();
1711     negV0Pt  = 1000000;
1712     negV0Eta = -2.;
1713     negV0Phi = -2*TMath::Pi();
1714     
1715     if(ConvGamV0RecIndexPos[i]!=-1){
1716       AliESDv0 * fV0MIs = fESD->GetV0(ConvGamV0RecIndexPos[i]);
1717       AliESDtrack* trackPosTest = fESD->GetTrack(fV0MIs->GetPindex());
1718       AliESDtrack* trackNegTest = fESD->GetTrack(fV0MIs->GetNindex());
1719
1720       if (ePosConvGamV0RecIndex[i]!=-1 ){
1721         //      AliESDtrack * ePosV0Track = fESD->GetTrack(ePosConvGamV0RecIndex[i]);
1722         if ( trackPosTest->GetSign()==1 ) {
1723           fV0MIs->GetPPxPyPz(ppV0[0],ppV0[1],ppV0[2]);
1724         }else{
1725           fV0MIs->GetNPxPyPz(ppV0[0],ppV0[1],ppV0[2]);
1726         }
1727         posV0Track.SetXYZM(ppV0[0],ppV0[1],ppV0[2],massE);
1728
1729         posV0Pt  = posV0Track.Pt();
1730         posV0Eta = posV0Track.Eta();
1731         posV0Phi = posV0Track.Phi();
1732         statusV0Pos=1;
1733       }else{
1734         posV0Pt  = 1000000;
1735         posV0Eta = -2.;
1736         posV0Phi = -2*TMath::Pi();
1737         statusV0Pos=-1;
1738       }
1739       
1740       if (eNegConvGamV0RecIndex[i]!=-1 ){
1741         //      AliESDtrack * eNegV0Track = fESD->GetTrack(eNegConvGamV0RecIndex[i]);
1742         if ( trackNegTest->GetSign()==-1 ) {
1743           fV0MIs->GetNPxPyPz(pmV0[0],pmV0[1],pmV0[2]);
1744         }else{
1745           fV0MIs->GetPPxPyPz(pmV0[0],pmV0[1],pmV0[2]);
1746         }
1747         negV0Track.SetXYZM(pmV0[0],pmV0[1],pmV0[2],massE);
1748         
1749         negV0Pt  = negV0Track.Pt();
1750         negV0Eta = negV0Track.Eta();
1751         negV0Phi = negV0Track.Phi();
1752         statusV0Neg=1;
1753       }else{
1754         negV0Pt  = 1000000;
1755         negV0Eta = -2.;
1756         negV0Phi = -2*TMath::Pi();
1757         statusV0Neg=-1;
1758       }
1759     }
1760     
1761     xrG[0] = ePosPart->Vx();
1762     xrG[1] = ePosPart->Vy();
1763     xrG[2] = ePosPart->Vz();
1764     
1765     //--------- Geant variables ----------------------
1766     fValueV0[0] = 1./TMath::Sqrt(gamPart->Pt());
1767     fValueV0[1] = gamPart->Eta();
1768
1769     Double_t tmpGPhi=gamPart->Phi();
1770     if( gamPart->Phi()>TMath::Pi()){
1771       tmpGPhi=gamPart->Phi()-2*TMath::Pi();
1772     }
1773     fValueV0[2] = tmpGPhi;
1774
1775     fValueV0[3] = TMath::Sqrt(xrG[0]*xrG[0]+xrG[1]*xrG[1]);
1776     fValueV0[4] = xrG[2];
1777  
1778
1779     fValueV0[5] = 1./TMath::Sqrt(ePosPart->Pt());
1780     fValueV0[6] = ePosPart->Eta();
1781
1782     Double_t tmpPPhi=ePosPart->Phi();
1783     if( ePosPart->Phi()>TMath::Pi()){
1784       tmpPPhi = ePosPart->Phi()-2*TMath::Pi();
1785     }
1786      fValueV0[7] = tmpPPhi;
1787      fValueV0[8] = ePosConvGamGeantLength[i];
1788
1789     fValueV0[9] = 1./TMath::Sqrt(eNegPart->Pt());
1790     fValueV0[10] = eNegPart->Eta();
1791
1792     Double_t tmpNPhi=eNegPart->Phi();
1793     if( eNegPart->Phi()>TMath::Pi()){
1794       tmpNPhi = eNegPart->Phi()-2*TMath::Pi();
1795     }
1796     fValueV0[11] = tmpNPhi;
1797     fValueV0[12] = eNegConvGamGeantLength[i];    
1798
1799     //---- Single track variables----------------------
1800
1801     fValueV0[13] = (posPt-ePosPart->Pt())/ePosPart->Pt();
1802     fValueV0[14] = posEta;
1803     fValueV0[15] = posPhi;
1804     fValueV0[16] = TMath::Sqrt( bPosSgl[0]* bPosSgl[0] +  bPosSgl[1]* bPosSgl[1] );
1805     fValueV0[17] = TMath::Sqrt( bPosSgl[0]* bPosSgl[0] +  bPosSgl[1]* bPosSgl[1] )/TMath::Sqrt(bPosCov[0]*bPosCov[0]+bPosCov[2]*bPosCov[2]);
1806     fValueV0[18] = nClsITSPos;
1807     fValueV0[19] = nClsTPCPos;
1808     fValueV0[20] = statusSingPos;
1809
1810
1811     fValueV0[21] = (negPt-eNegPart->Pt())/eNegPart->Pt();
1812     fValueV0[22] = negEta;
1813     fValueV0[23] = negPhi;
1814     fValueV0[24] = TMath::Sqrt( bNegSgl[0]* bNegSgl[0] +  bNegSgl[1]* bNegSgl[1] );
1815     fValueV0[25] = TMath::Sqrt( bNegSgl[0]* bNegSgl[0] +  bNegSgl[1]* bNegSgl[1] )/TMath::Sqrt(bNegCov[0]*bNegCov[0]+bNegCov[2]*bNegCov[2]);
1816
1817     fValueV0[26] = nClsITSNeg;
1818     fValueV0[27] = nClsTPCNeg;
1819     fValueV0[28] = statusSingNeg;
1820
1821     
1822     //---- V0 track variables----------------------
1823
1824     fValueV0[29] = (posV0Pt-ePosPart->Pt())/ePosPart->Pt();
1825     fValueV0[30] = posV0Eta;
1826     fValueV0[31] = posV0Phi;
1827     fValueV0[32] = statusV0Pos;
1828
1829     fValueV0[33] = (negV0Pt-eNegPart->Pt())/eNegPart->Pt();
1830     fValueV0[34] = negV0Eta;
1831     fValueV0[35] = negV0Phi;
1832     fValueV0[36] = statusV0Neg;
1833
1834     fSparseV0->Fill(fValueV0);
1835   }
1836
1837
1838 }
1839
1840 void AliAnalysisTaskV0QA::FillHnSparseK0()
1841 {
1842
1843   Double_t massPi=0.13957018;
1844   Double_t ppSgl[3];
1845   Double_t pmSgl[3];
1846   Float_t bPosSgl[2];
1847   Float_t bNegSgl[2];
1848   Float_t bPosCov[3];
1849   Float_t bNegCov[3];
1850   
1851   Double_t ppV0[3];
1852   Double_t pmV0[3];
1853   Double_t xrG[3];
1854   
1855   TLorentzVector posSglTrack;
1856   TLorentzVector negSglTrack;
1857   Double_t posPt,posEta,posPhi;
1858   Double_t negPt,negEta,negPhi;
1859
1860   TLorentzVector posV0Track;
1861   TLorentzVector negV0Track;
1862   Double_t posV0Pt,posV0Eta,posV0Phi;
1863   Double_t negV0Pt,negV0Eta,negV0Phi;
1864   
1865   Float_t nClsITSPos=-1;
1866   Float_t nClsITSNeg=-1;
1867
1868   Float_t nClsTPCPos=-1;
1869   Float_t nClsTPCNeg=-1;
1870
1871   Int_t statusSingPos=-1;
1872   Int_t statusSingNeg=-1;
1873
1874   Int_t statusV0Pos=-1;
1875   Int_t statusV0Neg=-1;
1876   
1877   for(Int_t i=0;i<nDecayK0Geant+1;i++){
1878     TParticle* K0Part = stack->Particle(K0DecayK0GeantIndex[i]);
1879     TParticle* ePosPart = stack->Particle(piPosDecayK0GeantIndex[i]);
1880     TParticle* eNegPart = stack->Particle(piNegDecayK0GeantIndex[i]);
1881     if (piPosDecayK0SingleRecIndex[i]!=-1){
1882       AliESDtrack * ePosSglTrack = fESD->GetTrack(piPosDecayK0SingleRecIndex[i]);
1883       ePosSglTrack->GetPxPyPz(ppSgl); 
1884       posSglTrack.SetXYZM(ppSgl[0],ppSgl[1],ppSgl[2],massPi);
1885       posPt  = posSglTrack.Pt();
1886       posEta = posSglTrack.Eta();
1887       posPhi = posSglTrack.Phi();
1888       ePosSglTrack->GetImpactParameters(bPosSgl,bPosCov);
1889       nClsITSPos=ePosSglTrack->GetNcls(0);
1890       nClsTPCPos=ePosSglTrack->GetNcls(1);
1891       statusSingPos=1;
1892     }else{
1893       posPt  = 1000000;
1894       posEta = -2.;
1895       posPhi = -2*TMath::Pi();
1896       bPosSgl[0]=-100.;
1897       bPosSgl[1]=-100.;
1898       bPosCov[0]=-100;
1899       bPosCov[2]=-100;
1900
1901       nClsITSPos=-1;
1902       nClsTPCPos=-1;
1903       statusSingPos=-1;
1904     }
1905
1906     if (piNegDecayK0SingleRecIndex[i]!=-1){
1907       AliESDtrack * eNegSglTrack = fESD->GetTrack(piNegDecayK0SingleRecIndex[i]);
1908       eNegSglTrack->GetPxPyPz(pmSgl); 
1909       negSglTrack.SetXYZM(pmSgl[0],pmSgl[1],pmSgl[2],massPi);
1910       negPt  = negSglTrack.Pt();
1911       negEta = negSglTrack.Eta();
1912       negPhi = negSglTrack.Phi();
1913       eNegSglTrack->GetImpactParameters(bNegSgl,bNegCov);
1914       nClsITSNeg=eNegSglTrack->GetNcls(0);
1915       nClsTPCNeg=eNegSglTrack->GetNcls(1);
1916       statusSingNeg=1;
1917     }else{
1918       negPt  = 1000000;
1919       negEta = -2.;
1920       negPhi = -2*TMath::Pi();
1921       bNegSgl[0]=-100.;
1922       bNegSgl[1]=-100.;
1923       bNegCov[0]=-100;
1924       bNegCov[2]=-100;
1925       nClsITSNeg=-1;
1926       nClsTPCNeg=-1;
1927       statusSingNeg=-1;
1928     }
1929
1930     posV0Pt  = 1000000;
1931     posV0Eta = -2.;
1932     posV0Phi = -2*TMath::Pi();
1933     negV0Pt  = 1000000;
1934     negV0Eta = -2.;
1935     negV0Phi = -2*TMath::Pi();
1936
1937     if(DecayK0V0RecIndexPos[i]!=-1){
1938       AliESDv0 * fV0MIs = fESD->GetV0(DecayK0V0RecIndexPos[i]);
1939       AliESDtrack* trackPosTest = fESD->GetTrack(fV0MIs->GetPindex());
1940       AliESDtrack* trackNegTest = fESD->GetTrack(fV0MIs->GetNindex());
1941
1942       if (piPosDecayK0V0RecIndex[i]!=-1 ){
1943         //      AliESDtrack * ePosV0Track = fESD->GetTrack(piPosDecayK0V0RecIndex[i]);
1944         if ( trackPosTest->GetSign()==1 ) {
1945           fV0MIs->GetPPxPyPz(ppV0[0],ppV0[1],ppV0[2]);
1946         }else{
1947           fV0MIs->GetNPxPyPz(ppV0[0],ppV0[1],ppV0[2]);
1948         }
1949         posV0Track.SetXYZM(ppV0[0],ppV0[1],ppV0[2],massPi);
1950
1951         posV0Pt  = posV0Track.Pt();
1952         posV0Eta = posV0Track.Eta();
1953         posV0Phi = posV0Track.Phi();
1954         statusV0Pos=1;
1955       }else{
1956         posV0Pt  = 1000000;
1957         posV0Eta = -2.;
1958         posV0Phi = -2*TMath::Pi();
1959         statusV0Pos=-1;
1960       }
1961       
1962       if (piNegDecayK0V0RecIndex[i]!=-1 ){
1963         //      AliESDtrack * eNegV0Track = fESD->GetTrack(piNegDecayK0V0RecIndex[i]);
1964         if ( trackNegTest->GetSign()==-1 ) {
1965           fV0MIs->GetNPxPyPz(pmV0[0],pmV0[1],pmV0[2]);
1966         }else{
1967           fV0MIs->GetPPxPyPz(pmV0[0],pmV0[1],pmV0[2]);
1968         }
1969         negV0Track.SetXYZM(pmV0[0],pmV0[1],pmV0[2],massPi);
1970         
1971         negV0Pt  = negV0Track.Pt();
1972         negV0Eta = negV0Track.Eta();
1973         negV0Phi = negV0Track.Phi();
1974         statusV0Neg=1;
1975       }else{
1976         negV0Pt  = 1000000;
1977         negV0Eta = -2.;
1978         negV0Phi = -2*TMath::Pi();
1979         statusV0Neg=-1;
1980       }
1981     }
1982
1983     xrG[0] = ePosPart->Vx();
1984     xrG[1] = ePosPart->Vy();
1985     xrG[2] = ePosPart->Vz();
1986     
1987     
1988     //--------- Geant variables ----------------------
1989     fValueK0[0] = 1./TMath::Sqrt(K0Part->Pt());
1990     fValueK0[1] = K0Part->Eta();
1991
1992     Double_t tmpGPhi=K0Part->Phi();
1993     if( K0Part->Phi()>TMath::Pi()){
1994       tmpGPhi=K0Part->Phi()-2*TMath::Pi();
1995     }
1996     fValueK0[2] = tmpGPhi;
1997
1998     fValueK0[3] = TMath::Sqrt(xrG[0]*xrG[0]+xrG[1]*xrG[1]);
1999     fValueK0[4] = xrG[2];
2000  
2001
2002     fValueK0[5] = 1./TMath::Sqrt(ePosPart->Pt());
2003     fValueK0[6] = ePosPart->Eta();
2004
2005     Double_t tmpPPhi=ePosPart->Phi();
2006     if( ePosPart->Phi()>TMath::Pi()){
2007       tmpPPhi = ePosPart->Phi()-2*TMath::Pi();
2008     }
2009      fValueK0[7] = tmpPPhi;
2010      fValueK0[8] = piPosDecayK0GeantLength[i];
2011
2012     fValueK0[9] = 1./TMath::Sqrt(eNegPart->Pt());
2013     fValueK0[10] = eNegPart->Eta();
2014
2015     Double_t tmpNPhi=eNegPart->Phi();
2016     if( eNegPart->Phi()>TMath::Pi()){
2017       tmpNPhi = eNegPart->Phi()-2*TMath::Pi();
2018     }
2019     fValueK0[11] = tmpNPhi;
2020     fValueK0[12] = piNegDecayK0GeantLength[i];    
2021     //---- Single track variables----------------------
2022
2023     fValueK0[13] = (posPt-ePosPart->Pt())/ePosPart->Pt() ;
2024     fValueK0[14] = posEta;
2025     fValueK0[15] = posPhi;
2026     fValueK0[16] = TMath::Sqrt( bPosSgl[0]* bPosSgl[0] +  bPosSgl[1]* bPosSgl[1] );
2027     fValueK0[17] = TMath::Sqrt( bPosSgl[0]* bPosSgl[0] +  bPosSgl[1]* bPosSgl[1] )/TMath::Sqrt(bPosCov[0]*bPosCov[0]+bPosCov[2]*bPosCov[2]);
2028    
2029     fValueK0[18] = nClsITSPos;
2030     fValueK0[19] = nClsTPCPos;
2031     fValueK0[20] = statusSingPos;
2032
2033     fValueK0[21] = (negPt-eNegPart->Pt())/eNegPart->Pt();
2034     fValueK0[22] = negEta;
2035     fValueK0[23] = negPhi;
2036     fValueK0[24] = TMath::Sqrt( bNegSgl[0]* bNegSgl[0]+  bNegSgl[1]* bNegSgl[1] );
2037     fValueK0[25] = TMath::Sqrt( bNegSgl[0]* bNegSgl[0]+  bNegSgl[1]* bNegSgl[1] )/TMath::Sqrt(bNegCov[0]*bNegCov[0]+bNegCov[2]*bNegCov[2]);
2038     fValueK0[26] = nClsITSNeg;
2039     fValueK0[27] = nClsTPCNeg;
2040     fValueK0[28] = statusSingNeg;
2041
2042     
2043     //---- V0 track variables----------------------
2044
2045     fValueK0[29] = (posV0Pt-ePosPart->Pt())/ePosPart->Pt();
2046     fValueK0[30] = posV0Eta;
2047     fValueK0[31] = posV0Phi;
2048     fValueK0[32] = statusV0Pos;
2049
2050     fValueK0[33] = (negV0Pt-eNegPart->Pt())/eNegPart->Pt();
2051     fValueK0[34] = negV0Eta;
2052     fValueK0[35] = negV0Phi;
2053     fValueK0[36] = statusV0Neg;
2054
2055     fSparseK0->Fill(fValueK0);
2056   }
2057   
2058
2059 }
2060 void AliAnalysisTaskV0QA::FillHnSparseL()
2061 {
2062
2063   Double_t massPi=0.13957018;
2064   Double_t massP=0.93827203;
2065
2066
2067   Double_t ppSgl[3];
2068   Double_t pmSgl[3];
2069   Float_t bPosSgl[2];
2070   Float_t bNegSgl[2];
2071   Float_t bPosCov[3];
2072   Float_t bNegCov[3];
2073   
2074   Double_t ppV0[3];
2075   Double_t pmV0[3];
2076   Double_t xrG[3];
2077   
2078   TLorentzVector posSglTrack;
2079   TLorentzVector negSglTrack;
2080   Double_t posPt,posEta,posPhi;
2081   Double_t negPt,negEta,negPhi;
2082
2083   TLorentzVector posV0Track;
2084   TLorentzVector negV0Track;
2085   Double_t posV0Pt,posV0Eta,posV0Phi;
2086   Double_t negV0Pt,negV0Eta,negV0Phi;
2087   
2088   Float_t nClsITSPos=-1;
2089   Float_t nClsITSNeg=-1;
2090
2091   Float_t nClsTPCPos=-1;
2092   Float_t nClsTPCNeg=-1;
2093
2094   Int_t statusSingPos=-1;
2095   Int_t statusSingNeg=-1;
2096
2097   Int_t statusV0Pos=-1;
2098   Int_t statusV0Neg=-1;
2099
2100   for(Int_t i=0;i<nDecayLGeant+1;i++){
2101     TParticle* lPart = stack->Particle(lDecayLGeantIndex[i]);
2102     TParticle* ePosPart = stack->Particle(pPosDecayLGeantIndex[i]);
2103     TParticle* eNegPart = stack->Particle(piNegDecayLGeantIndex[i]);
2104     if (pPosDecayLSingleRecIndex[i]!=-1){
2105       AliESDtrack * ePosSglTrack = fESD->GetTrack(pPosDecayLSingleRecIndex[i]);
2106       ePosSglTrack->GetPxPyPz(ppSgl); 
2107       posSglTrack.SetXYZM(ppSgl[0],ppSgl[1],ppSgl[2],massP);
2108       posPt  = posSglTrack.Pt();
2109       posEta = posSglTrack.Eta();
2110       posPhi = posSglTrack.Phi();
2111       ePosSglTrack->GetImpactParameters(bPosSgl,bPosCov);
2112       nClsITSPos=ePosSglTrack->GetNcls(0);
2113       nClsTPCPos=ePosSglTrack->GetNcls(1);
2114       statusSingPos=1;
2115     }else{
2116       posPt  = 1000000;
2117       posEta = -2.;
2118       posPhi = -2*TMath::Pi();
2119       bPosSgl[0]=-100.;
2120       bPosSgl[1]=-100.;
2121       bPosCov[0]=-100;
2122       bPosCov[2]=-100;
2123       nClsITSPos=-1;
2124       nClsTPCPos=-1;
2125       statusSingPos=-1;
2126     }
2127     
2128     if (piNegDecayLSingleRecIndex[i]!=-1){
2129       AliESDtrack * eNegSglTrack = fESD->GetTrack(piNegDecayLSingleRecIndex[i]);
2130       eNegSglTrack->GetPxPyPz(pmSgl); 
2131       negSglTrack.SetXYZM(pmSgl[0],pmSgl[1],pmSgl[2],massPi);
2132       negPt  = negSglTrack.Pt();
2133       negEta = negSglTrack.Eta();
2134       negPhi = negSglTrack.Phi();
2135       eNegSglTrack->GetImpactParameters(bNegSgl,bNegCov);
2136       nClsITSNeg=eNegSglTrack->GetNcls(0);
2137       nClsTPCNeg=eNegSglTrack->GetNcls(1);
2138       statusSingNeg=1;
2139     }else{
2140       negPt  = 1000000;
2141       negEta = -2.;
2142       negPhi = -2*TMath::Pi();
2143       bNegSgl[0]=-100.;
2144       bNegSgl[1]=-100.;
2145       bNegCov[0]=-100;
2146       bNegCov[2]=-100;
2147       nClsITSNeg=-1;
2148       nClsTPCNeg=-1;
2149       statusSingNeg=-1;
2150     }
2151
2152     posV0Pt  = 1000000;
2153     posV0Eta = -2.;
2154     posV0Phi = -2*TMath::Pi();
2155     negV0Pt  = 1000000;
2156     negV0Eta = -2.;
2157     negV0Phi = -2*TMath::Pi();
2158     
2159     if(DecayLV0RecIndexPos[i]!=-1){
2160       AliESDv0 * fV0MIs = fESD->GetV0(DecayLV0RecIndexPos[i]);
2161       AliESDtrack* trackPosTest = fESD->GetTrack(fV0MIs->GetPindex());
2162       AliESDtrack* trackNegTest = fESD->GetTrack(fV0MIs->GetNindex());
2163
2164       if (pPosDecayLV0RecIndex[i]!=-1 ){
2165         //      AliESDtrack * ePosV0Track = fESD->GetTrack(pPosDecayLV0RecIndex[i]);
2166         if ( trackPosTest->GetSign()==1 ) {
2167           fV0MIs->GetPPxPyPz(ppV0[0],ppV0[1],ppV0[2]);
2168         }else{
2169           fV0MIs->GetNPxPyPz(ppV0[0],ppV0[1],ppV0[2]);
2170         }
2171         posV0Track.SetXYZM(ppV0[0],ppV0[1],ppV0[2],massP);
2172
2173         posV0Pt  = posV0Track.Pt();
2174         posV0Eta = posV0Track.Eta();
2175         posV0Phi = posV0Track.Phi();
2176         statusV0Pos=1;
2177       }else{
2178         posV0Pt  = 1000000;
2179         posV0Eta = -2.;
2180         posV0Phi = -2*TMath::Pi();
2181         statusV0Pos=-1;
2182       }
2183       
2184       if (piNegDecayLV0RecIndex[i]!=-1 ){
2185         //      AliESDtrack * eNegV0Track = fESD->GetTrack(piNegDecayLV0RecIndex[i]);
2186         if ( trackNegTest->GetSign()==-1 ) {
2187           fV0MIs->GetNPxPyPz(pmV0[0],pmV0[1],pmV0[2]);
2188         }else{
2189           fV0MIs->GetPPxPyPz(pmV0[0],pmV0[1],pmV0[2]);
2190         }
2191         negV0Track.SetXYZM(pmV0[0],pmV0[1],pmV0[2],massPi);
2192         
2193         negV0Pt  = negV0Track.Pt();
2194         negV0Eta = negV0Track.Eta();
2195         negV0Phi = negV0Track.Phi();
2196         statusV0Neg=1;
2197       }else{
2198         negV0Pt  = 1000000;
2199         negV0Eta = -2.;
2200         negV0Phi = -2*TMath::Pi();
2201         statusV0Neg=-1;
2202       }
2203     }
2204     
2205     xrG[0] = ePosPart->Vx();
2206     xrG[1] = ePosPart->Vy();
2207     xrG[2] = ePosPart->Vz();
2208     
2209     //--------- Geant variables ----------------------
2210     fValueL[0] = 1./TMath::Sqrt(lPart->Pt());
2211     fValueL[1] = lPart->Eta();
2212
2213     Double_t tmpGPhi=lPart->Phi();
2214     if( lPart->Phi()>TMath::Pi()){
2215       tmpGPhi=lPart->Phi()-2*TMath::Pi();
2216     }
2217     fValueL[2] = tmpGPhi;
2218
2219     fValueL[3] = TMath::Sqrt(xrG[0]*xrG[0]+xrG[1]*xrG[1]);
2220     fValueL[4] = xrG[2];
2221  
2222
2223     fValueL[5] = 1./TMath::Sqrt(ePosPart->Pt());
2224     fValueL[6] = ePosPart->Eta();
2225
2226     Double_t tmpPPhi=ePosPart->Phi();
2227     if( ePosPart->Phi()>TMath::Pi()){
2228       tmpPPhi = ePosPart->Phi()-2*TMath::Pi();
2229     }
2230      fValueL[7] = tmpPPhi;
2231      fValueL[8] = pPosDecayLGeantLength[i];
2232
2233     fValueL[9] = 1./TMath::Sqrt(eNegPart->Pt());
2234     fValueL[10] = eNegPart->Eta();
2235
2236     Double_t tmpNPhi=eNegPart->Phi();
2237     if( eNegPart->Phi()>TMath::Pi()){
2238       tmpNPhi = eNegPart->Phi()-2*TMath::Pi();
2239     }
2240     fValueL[11] = tmpNPhi;
2241     fValueL[12] = piNegDecayLGeantLength[i];    
2242     //---- Single track variables----------------------
2243
2244     fValueL[13] = (posPt-ePosPart->Pt())/ePosPart->Pt();
2245     fValueL[14] = posEta;
2246     fValueL[15] = posPhi;
2247     fValueL[16] = TMath::Sqrt( bPosSgl[0]* bPosSgl[0] +  bPosSgl[1]* bPosSgl[1]);
2248     fValueL[17] = TMath::Sqrt( bPosSgl[0]* bPosSgl[0] +  bPosSgl[1]* bPosSgl[1] )/TMath::Sqrt(bPosCov[0]*bPosCov[0]+bPosCov[2]*bPosCov[2]);   
2249     fValueL[18] = nClsITSPos;
2250     fValueL[19] = nClsTPCPos;
2251     fValueL[20] = statusSingPos;
2252
2253     fValueL[21] = (negPt-eNegPart->Pt())/eNegPart->Pt() ;
2254     fValueL[22] = negEta;
2255     fValueL[23] = negPhi;
2256     fValueL[24] = TMath::Sqrt( bNegSgl[0]* bNegSgl[0] +  bNegSgl[1]* bNegSgl[1] );
2257     fValueL[25] = TMath::Sqrt( bNegSgl[0]* bNegSgl[0] +  bNegSgl[1]* bNegSgl[1] )/TMath::Sqrt(bNegCov[0]*bNegCov[0]+bNegCov[2]*bNegCov[2]);
2258     fValueL[26] = nClsITSNeg;
2259     fValueL[27] = nClsTPCNeg;
2260     fValueL[28] = statusSingNeg;
2261
2262
2263     
2264     //---- V0 track variables----------------------
2265
2266     fValueL[29] = (posV0Pt-ePosPart->Pt())/ePosPart->Pt();
2267     fValueL[30] = posV0Eta;
2268     fValueL[31] = posV0Phi;
2269     fValueL[32] = statusV0Pos;
2270
2271
2272     fValueL[33] = (negV0Pt-eNegPart->Pt())/eNegPart->Pt();
2273     fValueL[34] = negV0Eta;
2274     fValueL[35] = negV0Phi;
2275     fValueL[36] = statusV0Neg;
2276
2277     fSparseL->Fill(fValueL);
2278   }
2279
2280
2281 }
2282
2283 void AliAnalysisTaskV0QA::FillHnSparseAL()
2284 {
2285
2286   Double_t massPi=0.13957018;
2287   Double_t massP=0.93827203;
2288
2289
2290   Double_t ppSgl[3];
2291   Double_t pmSgl[3];
2292   Float_t bPosSgl[2];
2293   Float_t bNegSgl[2];
2294   Float_t bPosCov[3];
2295   Float_t bNegCov[3];
2296   
2297   Double_t ppV0[3];
2298   Double_t pmV0[3];
2299   Double_t xrG[3];
2300   
2301   TLorentzVector posSglTrack;
2302   TLorentzVector negSglTrack;
2303   Double_t posPt,posEta,posPhi;
2304   Double_t negPt,negEta,negPhi;
2305
2306   TLorentzVector posV0Track;
2307   TLorentzVector negV0Track;
2308   Double_t posV0Pt,posV0Eta,posV0Phi;
2309   Double_t negV0Pt,negV0Eta,negV0Phi;
2310   
2311   Float_t nClsITSPos=-1;
2312   Float_t nClsITSNeg=-1;
2313
2314   Float_t nClsTPCPos=-1;
2315   Float_t nClsTPCNeg=-1;
2316
2317   Int_t statusSingPos=-1;
2318   Int_t statusSingNeg=-1;
2319
2320   Int_t statusV0Pos=-1;
2321   Int_t statusV0Neg=-1;
2322
2323
2324   for(Int_t i=0;i<nDecayALGeant+1;i++){
2325     TParticle* alPart = stack->Particle(alDecayALGeantIndex[i]);
2326     TParticle* eNegPart = stack->Particle(apNegDecayALGeantIndex[i]);
2327     TParticle* ePosPart = stack->Particle(piPosDecayALGeantIndex[i]);
2328     if (piPosDecayALSingleRecIndex[i]!=-1){
2329       AliESDtrack * ePosSglTrack = fESD->GetTrack(piPosDecayALSingleRecIndex[i]);
2330       ePosSglTrack->GetPxPyPz(ppSgl); 
2331       posSglTrack.SetXYZM(ppSgl[0],ppSgl[1],ppSgl[2],massPi);
2332       posPt  = posSglTrack.Pt();
2333       posEta = posSglTrack.Eta();
2334       posPhi = posSglTrack.Phi();
2335       ePosSglTrack->GetImpactParameters(bPosSgl,bPosCov);
2336       nClsITSPos=ePosSglTrack->GetNcls(0);
2337       nClsTPCPos=ePosSglTrack->GetNcls(1);
2338       statusSingPos=1;
2339     }else{
2340       posPt  = 1000000;
2341       posEta = -2.;
2342       posPhi = -2*TMath::Pi();
2343       bPosSgl[0]=-100.;
2344       bPosSgl[1]=-100.;
2345       bPosCov[0]=-100;
2346       bPosCov[2]=-100;
2347       nClsITSPos=-1;
2348       nClsTPCPos=-1;
2349       statusSingPos=-1;
2350     }
2351     
2352     if (apNegDecayALSingleRecIndex[i]!=-1){
2353       AliESDtrack * eNegSglTrack = fESD->GetTrack(apNegDecayALSingleRecIndex[i]);
2354       eNegSglTrack->GetPxPyPz(pmSgl); 
2355       negSglTrack.SetXYZM(pmSgl[0],pmSgl[1],pmSgl[2],massP);
2356       negPt  = negSglTrack.Pt();
2357       negEta = negSglTrack.Eta();
2358       negPhi = negSglTrack.Phi();
2359       eNegSglTrack->GetImpactParameters(bNegSgl,bNegCov);
2360       nClsITSNeg=eNegSglTrack->GetNcls(0);
2361       nClsTPCNeg=eNegSglTrack->GetNcls(1);
2362       statusSingNeg=1;
2363     }else{
2364       negPt  = 1000000;
2365       negEta = -2.;
2366       negPhi = -2*TMath::Pi();
2367       bNegSgl[0]=-100.;
2368       bNegSgl[1]=-100.;
2369       bNegCov[0]=-100;
2370       bNegCov[2]=-100;
2371       nClsITSNeg=-1;
2372       nClsTPCNeg=-1;
2373       statusSingNeg=-1;
2374     }
2375
2376     posV0Pt  = 1000000;
2377     posV0Eta = -2.;
2378     posV0Phi = -2*TMath::Pi();
2379     negV0Pt  = 1000000;
2380     negV0Eta = -2.;
2381     negV0Phi = -2*TMath::Pi();
2382     
2383     if(DecayALV0RecIndexPos[i]!=-1){
2384       AliESDv0 * fV0MIs = fESD->GetV0(DecayALV0RecIndexPos[i]);
2385       AliESDtrack* trackPosTest = fESD->GetTrack(fV0MIs->GetPindex());
2386       AliESDtrack* trackNegTest = fESD->GetTrack(fV0MIs->GetNindex());
2387
2388       if (piPosDecayALV0RecIndex[i]!=-1 ){
2389         //      AliESDtrack * ePosV0Track = fESD->GetTrack(piPosDecayALV0RecIndex[i]);
2390         if ( trackPosTest->GetSign()==1 ) {
2391           fV0MIs->GetPPxPyPz(ppV0[0],ppV0[1],ppV0[2]);
2392         }else{
2393           fV0MIs->GetNPxPyPz(ppV0[0],ppV0[1],ppV0[2]);
2394         }
2395         posV0Track.SetXYZM(ppV0[0],ppV0[1],ppV0[2],massPi);
2396
2397         posV0Pt  = posV0Track.Pt();
2398         posV0Eta = posV0Track.Eta();
2399         posV0Phi = posV0Track.Phi();
2400         statusV0Pos=1;
2401       }else{
2402         posV0Pt  = 1000000;
2403         posV0Eta = -2.;
2404         posV0Phi = -2*TMath::Pi();
2405         statusV0Pos=-1;
2406       }
2407       
2408       if (apNegDecayALV0RecIndex[i]!=-1 ){
2409         //      AliESDtrack * eNegV0Track = fESD->GetTrack(apNegDecayALV0RecIndex[i]);
2410         if ( trackNegTest->GetSign()==-1 ) {
2411           fV0MIs->GetNPxPyPz(pmV0[0],pmV0[1],pmV0[2]);
2412         }else{
2413           fV0MIs->GetPPxPyPz(pmV0[0],pmV0[1],pmV0[2]);
2414         }
2415         negV0Track.SetXYZM(pmV0[0],pmV0[1],pmV0[2],massP);
2416         
2417         negV0Pt  = negV0Track.Pt();
2418         negV0Eta = negV0Track.Eta();
2419         negV0Phi = negV0Track.Phi();
2420         statusV0Neg=1;
2421       }else{
2422         negV0Pt  = 1000000;
2423         negV0Eta = -2.;
2424         negV0Phi = -2*TMath::Pi();
2425         statusV0Neg=-1;
2426       }
2427     }
2428     
2429     xrG[0] = ePosPart->Vx();
2430     xrG[1] = ePosPart->Vy();
2431     xrG[2] = ePosPart->Vz();
2432     
2433     //--------- Geant variables ----------------------
2434     fValueAL[0] = 1./TMath::Sqrt(alPart->Pt());
2435     fValueAL[1] = alPart->Eta();
2436
2437     Double_t tmpGPhi=alPart->Phi();
2438     if( alPart->Phi()>TMath::Pi()){
2439       tmpGPhi=alPart->Phi()-2*TMath::Pi();
2440     }
2441     fValueAL[2] = tmpGPhi;
2442
2443     fValueAL[3] = TMath::Sqrt(xrG[0]*xrG[0]+xrG[1]*xrG[1]);
2444     fValueAL[4] = xrG[2];
2445  
2446
2447     fValueAL[5] = 1./TMath::Sqrt(ePosPart->Pt());
2448     fValueAL[6] = ePosPart->Eta();
2449
2450     Double_t tmpPPhi=ePosPart->Phi();
2451     if( ePosPart->Phi()>TMath::Pi()){
2452       tmpPPhi = ePosPart->Phi()-2*TMath::Pi();
2453     }
2454      fValueAL[7] = tmpPPhi;
2455      fValueAL[8] = piPosDecayALGeantLength[i];
2456
2457     fValueAL[9] = 1./TMath::Sqrt(eNegPart->Pt());
2458     fValueAL[10] = eNegPart->Eta();
2459
2460     Double_t tmpNPhi=eNegPart->Phi();
2461     if( eNegPart->Phi()>TMath::Pi()){
2462       tmpNPhi = eNegPart->Phi()-2*TMath::Pi();
2463     }
2464     fValueAL[11] = tmpNPhi;
2465     fValueAL[12] = apNegDecayALGeantLength[i];    
2466     //---- Single track variables----------------------
2467
2468     fValueAL[13] = (posPt-ePosPart->Pt())/ePosPart->Pt();
2469     fValueAL[14] = posEta;
2470     fValueAL[15] = posPhi;
2471     fValueAL[16] = TMath::Sqrt( bPosSgl[0]* bPosSgl[0] +  bPosSgl[1]* bPosSgl[1]);
2472     fValueAL[17] = TMath::Sqrt( bPosSgl[0]* bPosSgl[0] +  bPosSgl[1]* bPosSgl[1] )/TMath::Sqrt(bPosCov[0]*bPosCov[0]+bPosCov[2]*bPosCov[2]);   
2473     fValueAL[18] = nClsITSPos;
2474     fValueAL[19] = nClsTPCPos;
2475     fValueAL[20] = statusSingPos;
2476
2477     fValueAL[21] = (negPt-eNegPart->Pt())/eNegPart->Pt() ;
2478     fValueAL[22] = negEta;
2479     fValueAL[23] = negPhi;
2480     fValueAL[24] = TMath::Sqrt( bNegSgl[0]* bNegSgl[0] +  bNegSgl[1]* bNegSgl[1] );
2481     fValueAL[25] = TMath::Sqrt( bNegSgl[0]* bNegSgl[0] +  bNegSgl[1]* bNegSgl[1] )/TMath::Sqrt(bNegCov[0]*bNegCov[0]+bNegCov[2]*bNegCov[2]);
2482     fValueAL[26] = nClsITSNeg;
2483     fValueAL[27] = nClsTPCNeg;
2484     fValueAL[28] = statusSingNeg;
2485
2486
2487     
2488     //---- V0 track variables----------------------
2489
2490     fValueAL[29] = (posV0Pt-ePosPart->Pt())/ePosPart->Pt();
2491     fValueAL[30] = posV0Eta;
2492     fValueAL[31] = posV0Phi;
2493     fValueAL[32] = statusV0Pos;
2494
2495
2496     fValueAL[33] = (negV0Pt-eNegPart->Pt())/eNegPart->Pt();
2497     fValueAL[34] = negV0Eta;
2498     fValueAL[35] = negV0Phi;
2499     fValueAL[36] = statusV0Neg;
2500
2501     fSparseAL->Fill(fValueAL);
2502   }
2503 }
2504
2505
2506 // void AliAnalysisTaskV0QA::SetESDtrackCuts()
2507 // {
2508
2509 //   fEsdTrackCuts = new AliESDtrackCuts("AliESDtrackCuts");
2510
2511 //   fEsdTrackCuts->SetRequireTPCRefit(kTRUE);
2512 //   fEsdTrackCuts->SetAcceptKinkDaughters(kFALSE);
2513
2514
2515
2516 // }