]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/hfe/AliHFEpidQA.cxx
Various updates, including corrections for code rule violations
[u/mrichter/AliRoot.git] / PWG3 / hfe / AliHFEpidQA.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15 //
16 // Class for PID QA
17 // Several studies done on clean samples of electrons, pions and kaons
18 // coming from V0 PID
19 // Compatible with both ESDs and AODs
20 //
21 // Autors:
22 //    Matus Kalisky <matus.kalisky@cern.ch>
23 //    Markus Heide <mheide@uni-muenster.de>
24 //    Markus Fasel <M.Fasel@gsi.de>
25 //
26
27
28 #include <TMath.h>
29 #include <TObjArray.h>
30 #include <TPDGCode.h>
31 #include <TString.h>
32 #include <TMultiLayerPerceptron.h>
33 #include <TFile.h>
34
35 #include "AliAODMCParticle.h"
36 #include "AliAODEvent.h"
37 #include "AliAODTrack.h"
38 #include "AliESDtrack.h"
39 #include "AliESDEvent.h"
40 #include "AliMCEvent.h"
41 #include "AliMCParticle.h"
42 #include "AliPID.h"
43 #include "AliESDpid.h"
44 #include "AliVParticle.h"
45
46
47 #include "AliHFEcollection.h"
48 #include "AliHFEpidQA.h"
49 #include "AliHFEV0info.h"
50 #include "AliHFEV0pid.h"
51 #include "AliHFEV0pidMC.h"
52 #include "AliHFEV0cuts.h"
53 #include "AliHFEtrdPIDqa.h"
54
55
56 ClassImp(AliHFEpidQA)
57
58   //__________________________________________
59   AliHFEpidQA::AliHFEpidQA():
60     fEvent(NULL),
61     fMC(NULL), 
62     fV0pid(NULL), 
63     fV0pidMC(NULL), 
64     fTRDpidQA(NULL), 
65     fOutput(NULL), 
66     fESDpid(NULL),
67     fNNref(NULL)
68 {
69   //
70   // Default constructor
71   //
72   for(Int_t mom = 0; mom < 11; mom++){
73     fNet[mom] = NULL;
74   }
75 }
76
77 //__________________________________________
78 AliHFEpidQA::AliHFEpidQA(const AliHFEpidQA &ref):
79   TObject(ref),
80   fEvent(NULL), 
81   fMC(NULL),
82   fV0pid(NULL),
83   fV0pidMC(NULL), 
84   fTRDpidQA(NULL),
85   fOutput(NULL), 
86   fESDpid(NULL),
87   fNNref(NULL)
88 {
89   //
90   // Copy constructor
91   //
92   for(Int_t mom = 0; mom < 11; mom++){
93     fNet[mom] = NULL;
94   }
95   ref.Copy(*this);
96 }
97
98 //__________________________________________
99 AliHFEpidQA &AliHFEpidQA::operator=(const AliHFEpidQA &ref){
100   //
101   // Assignment operator
102   //
103   if(this != &ref)
104     ref.Copy(*this);
105   return *this;
106 }
107
108 //__________________________________________
109 AliHFEpidQA::~AliHFEpidQA(){
110   //
111   // Destructor
112   //
113   if(fV0pid) delete fV0pid;
114   if(fV0pidMC) delete fV0pidMC;
115   if(fTRDpidQA) delete fTRDpidQA;
116   if(fOutput) delete fOutput;
117   //  if(fTRDpidResponse) delete fTRDpidResponse; 
118 }
119
120 //__________________________________________
121 void AliHFEpidQA::Copy(TObject &o) const {
122   //
123   // Copy function
124   //
125   
126   TObject::Copy(o);
127
128   AliHFEpidQA &target = dynamic_cast<AliHFEpidQA &>(o);
129   target.fMC = fMC;
130
131   if(target.fESDpid) delete target.fESDpid;
132   target.fESDpid = new AliESDpid;
133   if(target.fV0pid) delete target.fV0pid;
134   if(fV0pid)
135     target.fV0pid = dynamic_cast<AliHFEV0pid *>(fV0pid->Clone());
136   else
137     target.fV0pid = NULL;
138   if(target.fV0pidMC) delete target.fV0pidMC;
139   if(fV0pidMC) 
140     target.fV0pidMC = dynamic_cast<AliHFEV0pidMC *>(fV0pidMC->Clone());
141   else
142     target.fV0pidMC = NULL;
143   if(target.fTRDpidQA) delete target.fTRDpidQA;
144   if(fTRDpidQA)
145     target.fTRDpidQA = dynamic_cast<AliHFEtrdPIDqa *>(fTRDpidQA->Clone());
146   else
147     target.fTRDpidQA = NULL;
148   if(target.fOutput) delete target.fOutput;
149   if(fOutput)
150     target.fOutput = dynamic_cast<AliHFEcollection *>(fOutput->Clone());
151 }
152
153 //__________________________________________
154 void AliHFEpidQA::Init(){
155   //
156   // Prepare task output
157   //
158
159   // load networks
160   if(fNNref){
161     for(Int_t mom = 0; mom < 11; mom++){                      
162       fNet[mom] = (TMultiLayerPerceptron*) fNNref->Get(Form("NN_Mom%d", mom));
163       if(!fNet[mom]){
164         AliError(Form("No reference network for momentum bin %d!", mom));
165       }
166     }
167   }
168   
169   fV0pid = new AliHFEV0pid();
170   if(HasV0pidQA()) fV0pid->InitQA();
171   fV0pidMC = new AliHFEV0pidMC();
172   fV0pidMC->Init();
173
174   fTRDpidQA = new AliHFEtrdPIDqa;
175   fTRDpidQA->Init();
176
177   fOutput = new AliHFEcollection("pidQA", "PID QA output");
178
179   const char *name[4] = {"Electron", "PionK0", "PionL", "Proton"};
180   const char *title[4] = {"Electron", "K0 Pion", "Lambda Pion", "Proton"};
181   const char *det[4] = {"ITS", "TPC", "TRD", "TOF"};
182   for(Int_t i = 0; i < 4; i++){
183     fOutput->CreateTH2F(Form("purity%s", name[i]), Form("%s Purity", title[i]), 2, -0.5, 1.5, 20, 0.1, 10, 1);
184
185     for(Int_t idet = 0; idet < 4; idet++){
186       // create all the histograms which all the detectors have in common
187       if(idet != 2){ // no nSigma histogram for TRD
188         fOutput->CreateTH2F(Form("h%s_nSigma_%s", det[idet], name[i]), Form("%s number of sigmas for %ss; p (GeV/c); number of sigmas", det[idet], title[i]), 20, 0.1, 10, 100, -7, 7, 0);
189       }
190       fOutput->CreateTH2F(Form("h%s_PID_p_%s", det[idet], name[i]), Form("%s PID for %ss; p (GeV/c); ITS PID", det[idet], title[i]), 100, 0.1, 10, 5, -0.5, 4.5, 0);
191       fOutput->CreateTH2F(Form("h%s_El_like_%s", det[idet], name[i]), Form("%s Electron Likelihoods for %ss; p (GeV/c); likelihood", det[idet], title[i]), 25, 0.1, 10, 1000, 0., 1., 0);
192       fOutput->CreateTH2F(Form("h%s_El_like_MC_%s", det[idet], name[i]), Form("%s Electron Likelihoods for MC %ss; p (GeV/c); likelihood", det[idet], title[i]), 25, 0.1, 10, 1000, 0., 1., 0);
193     }
194
195     //
196     // ITS pid response
197     //
198     fOutput->CreateTH2F(Form("hITS_Signal_%s", name[i]), Form("ITS Signal org. for %ss", title[i]), 40, 0.1, 10, 400, 0, 1000, 0);
199     fOutput->CreateTH2Fvector1(5, Form("hITS_dEdx_%s", name[i]), Form("ITS dEdx for %ss; p (GeV/c); dEdx (a.u.)", title[i]), 40, 0.1, 10, 400, 0, 1000, 0);
200     
201     //
202     // TPC pid response
203     //
204     fOutput->CreateTH2F(Form("hTPC_dEdx_%s", name[i]), Form("TPC dEdx for %ss; p (GeV/c); dEdx (a.u.)", title[i]), 20, 0.1, 10, 200, 0, 200, 0);
205
206     //
207     // TRD pid response 
208     //
209     fOutput->CreateTH2F(Form("hTRD_trk_%s", name[i]), Form("%s PID tracklets; p (GeV/c); N TRD tracklets", title[i]), 100, 0.1, 10, 7, -0.5, 6.5, 0);
210     // number of the non 0 slices per tracklet
211     fOutput->CreateTH2F(Form("hTRD_Nslc_%s", name[i]), Form("%s PID slices > 0; p (GeV/c); N slc > 0", title[i]), 100, 0.1, 10, 9, -0.5, 8.5, 0);
212     // location of the slices > 0 - where are the emtpy slices located ?
213     fOutput->CreateTH2F(Form("hTRD_slc_%s", name[i]), Form("%s PID slices > 0 position; p (GeV/c); slice", title[i]), 100, 0.1, 10, 9, -0.5, 8.5, 0);
214     fOutput->CreateTH2F(Form("hTRD_cls_%s", name[i]), Form("TRD clusters for %s candidates; p (GeV/c); N cls", title[i]), 25, 0.1, 10, 1000, 0, 1000);
215     fOutput->CreateTH2F(Form("hTRD_dEdx_%s", name[i]), Form("TRD dEdx (trk) for %s candidates; p (GeV/c); tracklet dEdx (a.u.)", title[i]), 25, 0.1, 10, 1000, 0, 100000, 0);
216  
217     //
218     // TOF pid response
219     //
220     fOutput->CreateTH2F(Form("hTOF_beta_%s", name[i]), Form("TOF beta for %s; p (GeV/c); #beta", title[i]),  50, 0.1, 5, 350, 0.4, 1.1, 0);
221   }//.. loop over identified particle species
222
223   // Global histograms
224   fOutput->CreateTH1F("hITS_dEdx_nSamples", "ITS - number of non 0 dEdx samples", 4, 0.5, 4.5);
225   fOutput->CreateTH2F("hTPC_PID", "TPC pid all tracks; tpc pid probability; species",100, 0, 1, 5, -0.5, 4.5 );
226   fOutput->CreateTH2F("hTOF_PID", "TOF pid all tracks; tof pid probability; species",100, 0, 1,5,  -0.5, 4.5 );
227   fOutput->CreateTH2F("hTOF_beta_all", "TOF beta for all nice single tracks; p (GeV/c); #beta", 100, 0.1, 10, 480, 0, 1.2, 0);
228   fOutput->CreateTH2F("hTOF_qa_TmT0mT", "TOF (t - t0 - t[pion]) qa verus run",  10000, 114000, 124000, 200, -200, 200);
229   fOutput->CreateTH2F("hTOF_qa_length", "TOF track length verus run",  10000, 114000, 112400, 200, 0, 1000);
230   
231   //
232   // debug histograms
233   //
234   fOutput->CreateTH2F("hITS_kFlags", "ITS flags; flag; V0 candidates", 5, 0.5, 5.5, 5, -0.5, 4.5);  
235   fOutput->CreateTH2F("hTPC_kFlags", "TPC flags; flag; V0 candidates", 5, 0.5, 5.5, 5, -0.5, 4.5);
236   fOutput->CreateTH2F("hTRD_kFlags", "TRD flags; flag; V0 candidates", 5, 0.5, 5.5, 5, -0.5, 4.5);
237   fOutput->CreateTH2F("hTOF_kFlags", "TOF flags; flag; V0 candidates", 5, 0.5, 5.5, 5, -0.5, 4.5);
238   
239
240   //
241   // event based histograms
242   //
243   Int_t nT0[2] = {10000, 100};
244   Double_t minT0[2] = {114500, -1000};
245   Double_t maxT0[2] = {124500, 1000};
246   fOutput->CreateTHnSparse("hEvent_T0", "T0 as a function of run number; run number; T0 (ns)", 2, nT0, minT0, maxT0);
247
248   //
249   // test the tender V0 supply
250   //
251   fOutput->CreateTH2F("h_tender_check_01", "tender -vs- HFEpidQA; HFEpidQA V0 candiadates; tender V0 candidates", 4, -1.5, 2.5, 4, -1.5, 2.5);
252   fOutput->CreateTH2Fvector1(3, "h_tender_check_02", "tender -vs- HFEpidQA per type; AliHFEpidQA V0 ; tender V0", 4, -1.5, 2.5, 4, -1.5, 2.5);
253
254
255   //
256   // THnSpasre objects
257   //
258
259   // Create Illumination Plot
260   // bins: species, pt, eta, phi, TPC status, TRD status
261   {
262     Int_t nbins[6] = {4, 40, 16, 72, 2, 2};
263     Double_t min[6] = { 0, 0.1, -0.8, 0., 0., 0.};
264     Double_t max[6] = { 4., 20., 0.8, 2*TMath::Pi(), 2., 2.};
265     fOutput->CreateTHnSparse("hIllumination", "Illumination", 6, nbins, min, max);
266     fOutput->BinLogAxis("hIllumination", 1);
267   }
268
269   // TPC THnSparse
270   // bins: species, pt, n TPC clusters., TPC electron likelihood, TPC n sigmas, TPC signal
271   {
272     Int_t nbins[6] = { 5, 40, 20, 100, 100, 200};
273     Double_t min[6] = { -0.5, 0.1, 0., 0., -5., 0.};
274     Double_t max[6] = { 4.5, 20., 200, 1., 5., 120.};
275     TString htitle = "TPC info sparse; VO identified species; p (GeV/c); n TPC clusters; TPC N sigma; TPC signal";
276     fOutput->CreateTHnSparse("hTPCclusters",htitle,6, nbins, min, max);
277     fOutput->BinLogAxis("hTPCclusters", 1);
278   }
279   // TRD THnSparse - entries per tracklet
280   // species, p, tracklet position, number of PID tracklets, number of slices (non 0), number of clusters, electron likelihood, 
281
282   {
283     Int_t nbins[7] = {5, 20, 6, 7, 9, 45, 100};
284     Double_t min[7] = {-0.5, 0.5, -0.5, -0.5, -0.5, -0.5, 0.};
285     Double_t max[7] = {4.5, 10, 5.5, 6.5, 8.5, 179.5, 1.};
286     TString htitle = "TRD tracklets; VO identified species; p (GeV/c); tracklet position; No. PID tacklets; No. slices; No. clusters; El. likelihood";
287     fOutput->CreateTHnSparse("hTRDtracklets",htitle,7, nbins, min, max);
288     fOutput->BinLogAxis("hTRDtracklets", 1);
289   }
290   
291
292 }
293 //__________________________________________
294 void AliHFEpidQA::Process(){
295   //
296   // Run PID QA
297   //
298
299   if(!fV0pid){
300     AliError("V0pid not available! Forgotten to initialize?");
301     return;
302   }
303   if(!fESDpid){
304     AliError("fESDpid not initialized, I am leaving this!");
305     return;
306   }
307
308   // to be udpated to AOD save mdoe
309   if(!fEvent){
310     AliError("AliVEvent not available, returning");
311   }
312
313
314   if(fMC) fV0pidMC->SetMCEvent(fMC);
315   if(fMC) fV0pid->SetMCEvent(fMC);
316
317   fV0pid->Process(fEvent);
318   TObjArray *hfeelectrons = fV0pid->GetListOfElectrons();
319   TObjArray *hfepionsK0 = fV0pid->GetListOfPionsK0();
320   TObjArray *hfepionsL = fV0pid->GetListOfPionsL();
321   TObjArray *hfeprotons = fV0pid->GetListOfProtons();
322
323   // Get Track list for normal purpose
324   TObjArray *electrons = MakeTrackList(hfeelectrons);
325   TObjArray *pionsK0 = MakeTrackList(hfepionsK0);
326   TObjArray *pionsL = MakeTrackList(hfepionsL);
327   TObjArray *protons = MakeTrackList(hfeprotons);
328   TObjArray *cleanElectrons = MakeCleanListElectrons(hfeelectrons);
329
330   if(fMC){
331     fV0pidMC->Process(electrons, AliHFEV0cuts::kRecoElectron);
332     fV0pidMC->Process(pionsK0, AliHFEV0cuts::kRecoPionK0);
333     fV0pidMC->Process(pionsL, AliHFEV0cuts::kRecoPionL);
334     fV0pidMC->Process(protons, AliHFEV0cuts::kRecoProton);
335   }
336
337   AliDebug(2, Form("Number of Electrons      : %d", electrons->GetEntries()));
338   AliDebug(2, Form("Number of K0 Pions       : %d", pionsK0->GetEntries()));
339   AliDebug(2, Form("Number of Lambda Pions   : %d", pionsL->GetEntries()));
340   AliDebug(2, Form("Number of Protons        : %d", protons->GetEntries()));
341   if(fMC){
342     AliDebug(2, "MC Information available. Doing Purity checks...");
343     // Calculate the purity of the clean samples using MC 
344     MakePurity(electrons, AliHFEV0cuts::kRecoElectron);
345     MakePurity(pionsK0,  AliHFEV0cuts::kRecoPionK0);
346     MakePurity(pionsL,  AliHFEV0cuts::kRecoPionL);
347     MakePurity(protons,  AliHFEV0cuts::kRecoProton);
348   }
349
350   // some event wise checks
351   CheckEvent();
352
353   // Make Illumination Plot
354   FillIllumination(electrons, AliHFEV0cuts::kRecoElectron);
355   FillIllumination(pionsK0, AliHFEV0cuts::kRecoPionK0);
356   FillIllumination(pionsL, AliHFEV0cuts::kRecoPionL);
357   FillIllumination(protons, AliHFEV0cuts::kRecoProton);
358
359   // Now we can do studies on the PID itself
360   // For TRD use the TRD PID QA object
361   fTRDpidQA->ProcessTracks(cleanElectrons, AliPID::kElectron);
362   fTRDpidQA->ProcessTracks(pionsK0, AliPID::kPion);
363   fTRDpidQA->ProcessTracks(pionsL, AliPID::kPion);
364   fTRDpidQA->ProcessTracks(protons, AliPID::kProton);
365
366   FillElectronLikelihoods(electrons,  AliHFEV0cuts::kRecoElectron); 
367   FillElectronLikelihoods(pionsK0,  AliHFEV0cuts::kRecoPionK0); 
368   FillElectronLikelihoods(pionsL,  AliHFEV0cuts::kRecoPionL); 
369   FillElectronLikelihoods(protons,  AliHFEV0cuts::kRecoProton); 
370   
371   FillPIDresponse(electrons, AliHFEV0cuts::kRecoElectron);
372   FillPIDresponse(pionsK0, AliHFEV0cuts::kRecoPionK0);
373   FillPIDresponse(pionsL, AliHFEV0cuts::kRecoPionL);
374   FillPIDresponse(protons, AliHFEV0cuts::kRecoProton);
375
376   // check the tender V0s
377   CheckTenderV0pid(electrons, AliHFEV0cuts::kRecoElectron);
378   CheckTenderV0pid(pionsK0, AliHFEV0cuts::kRecoPionK0);
379   CheckTenderV0pid(pionsL, AliHFEV0cuts::kRecoPionL);
380   CheckTenderV0pid(protons, AliHFEV0cuts::kRecoProton);
381
382   // Analysis done, flush the containers
383   fV0pid->Flush();
384
385   delete electrons;
386   delete pionsL;
387   delete pionsK0;
388   delete protons;
389   delete cleanElectrons;
390 }
391
392 //__________________________________________
393 void AliHFEpidQA::FillIllumination(const TObjArray * const tracks, Int_t species){
394   //
395   // Fill Illumination Plot
396   //
397   THnSparseF *hIllumination = dynamic_cast<THnSparseF *>(fOutput->Get("hIllumination"));
398   if(!hIllumination) return;
399
400   Double_t quantities[6]; memset(quantities, 0, sizeof(Double_t) *6);
401   TIter trackIter(tracks);
402
403   quantities[0] = species;
404   TObject *o = NULL; AliESDtrack *esdtrack = NULL;
405   while((o = trackIter())){
406     if(!TString(o->IsA()->GetName()).CompareTo("AliESDtrack")){
407       // work on local copy in order to not spoil others
408       esdtrack = new AliESDtrack(*(dynamic_cast<AliESDtrack *>(o)));
409     } else if(!TString(o->IsA()->GetName()).CompareTo("AliAODrack")){
410       // Bad hack: Fill ESD track with AOD information
411       esdtrack = new AliESDtrack(dynamic_cast<AliAODTrack *>(o));
412     } else {
413       // Non usable
414       continue;
415     }
416
417     // Propagate to the entrance of the TRD
418     esdtrack->PropagateTo(300, fEvent->GetMagneticField());
419     quantities[1] = esdtrack->Pt();
420     quantities[2] = esdtrack->Eta();
421     quantities[3] = esdtrack->Phi();
422     quantities[4] = (esdtrack->GetStatus() & AliESDtrack::kTPCrefit) ? 1 : 0;
423     quantities[5] = (esdtrack->GetStatus() & AliESDtrack::kTRDout) ? 1. : 0.;
424     hIllumination->Fill(quantities);
425
426     delete esdtrack;
427   }
428 }
429 //__________________________________________
430 void AliHFEpidQA::FillTPCinfo(AliESDtrack *const esdtrack, Int_t species){
431   //
432   // Fill TPC Cluster Plots
433   //
434   THnSparseF *hTPCclusters = dynamic_cast<THnSparseF *>(fOutput->Get("hTPCclusters"));
435   if(!hTPCclusters) return;
436
437   Double_t quantities[6]; memset(quantities, 0, sizeof(Double_t) *6);
438   
439   Double_t pidProbs[5];
440   const Int_t typePID[5] = {0, 2, 2, 3, 4};
441
442
443   quantities[0] = species;
444
445   
446   esdtrack->GetTPCpid(pidProbs);   
447     
448   quantities[1] = esdtrack->P();
449   quantities[2] = esdtrack->GetTPCNcls();
450   quantities[3] = pidProbs[0];
451   quantities[4] = fESDpid->NumberOfSigmasTPC(esdtrack,(AliPID::EParticleType)typePID[species]);
452   quantities[5] = esdtrack->GetTPCsignal();
453   hTPCclusters->Fill(quantities);
454
455 }
456
457 //__________________________________________
458 void AliHFEpidQA::MakePurity(const TObjArray *tracks, Int_t species){
459   //
460   // Fill the QA histos for a given species
461   //
462   if(!fMC) return;
463   AliDebug(3, Form("Doing Purity checks for species %d", species));
464   Int_t pdg = GetPDG(species);
465   Char_t hname[256];
466
467   switch(species){
468   case  AliHFEV0cuts::kRecoElectron:
469     sprintf(hname, "purityElectron");
470     break;
471   case  AliHFEV0cuts::kRecoPionK0:
472     sprintf(hname, "purityPionK0");
473     break;
474   case  AliHFEV0cuts::kRecoPionL:
475     sprintf(hname, "purityPionL");
476     break;
477   case  AliHFEV0cuts::kRecoProton:
478     sprintf(hname, "purityProton");
479     break;
480   default:  // non investigated species
481     AliDebug(3, "Species not investigated");
482     return;
483   }  
484
485   AliDebug(3, Form("Number of tracks: %d", tracks->GetEntries()));
486   TIter trackIter(tracks);
487   AliVParticle *recTrack = NULL, *mcTrack = NULL;
488   while((recTrack = dynamic_cast<AliVParticle *>(trackIter()))){
489     Int_t label = recTrack->GetLabel();
490     AliDebug(4, Form("MC Label %d", label));
491     mcTrack =fMC->GetTrack(TMath::Abs(label));
492     if(!mcTrack){
493       AliDebug(4, "MC track not available");
494       continue; // we don't know
495     }
496
497     // Get the pdg code
498     Int_t trackPdg = 0;
499     if(!TString(mcTrack->IsA()->GetName()).CompareTo("AliMCParticle")){
500       // case ESD
501       AliMCParticle *mcp = dynamic_cast<AliMCParticle *>(mcTrack);
502       trackPdg = TMath::Abs(mcp->Particle()->GetPdgCode());
503     } else {
504       // case AOD
505       AliAODMCParticle *aodmcp = dynamic_cast<AliAODMCParticle *>(mcTrack);
506       trackPdg = TMath::Abs(aodmcp->GetPdgCode());
507     }
508     if(trackPdg == pdg)    // Correct identification
509       {
510         fOutput->Fill(hname, 0., recTrack->Pt());
511       }
512     else  // Wrong identification
513       fOutput->Fill(hname, 1., recTrack->Pt());
514   }
515 }
516
517 //__________________________________________
518 void AliHFEpidQA::FillElectronLikelihoods(const TObjArray * const particles, Int_t species){
519   //
520   // Fill electron Likelihoods for the ITS, TPC and TOF
521   // Required for the calculation of the electron efficiency, 
522   // pion and proton efficiency and the thresholds
523   //
524   Long_t status = 0;
525   const Char_t *detname[4] = {"ITS", "TPC", "TRD", "TOF"};
526   Char_t specname[256];
527
528   switch(species){
529   case  AliHFEV0cuts::kRecoElectron:
530     sprintf(specname, "Electron");
531     break;
532   case  AliHFEV0cuts::kRecoPionK0:
533     sprintf(specname, "PionK0");
534     break;
535   case  AliHFEV0cuts::kRecoPionL:
536     sprintf(specname, "PionL");
537     break;
538   case  AliHFEV0cuts::kRecoProton:
539     sprintf(specname, "Proton");
540     break;
541   default:
542     AliDebug(2, Form("Species %d not investigated", species));
543     return;
544   };
545   AliVParticle *recTrack = NULL;
546   //  mcTrack =fMC->GetTrack(TMath::Abs(label));
547   //   if(!mcTrack){
548   //     AliDebug(4, "MC track not available");
549   //     continue; // we don't know
550   //   }
551
552   TIter trackIter(particles);
553
554   Double_t quantities[2];
555   Double_t pidProbs[5];
556
557   while((recTrack = dynamic_cast<AliVParticle *>(trackIter()))){
558     if(!TString(recTrack->IsA()->GetName()).CompareTo("AliESDtrack")){
559       // case ESD
560       AliESDtrack *esdTrack = dynamic_cast<AliESDtrack *>(recTrack);
561       status = esdTrack->GetStatus();
562
563       //TPC momentum and likelihoods
564       Double_t pTPC = 0.;
565       pTPC = esdTrack->GetInnerParam() ? esdTrack->GetInnerParam()->P() : esdTrack->P();
566       Bool_t mcFound = kFALSE;
567       if(fMC){
568         Int_t label = esdTrack->GetLabel();
569         AliMCParticle *mcTrack = dynamic_cast<AliMCParticle *>(fMC->GetTrack(label));
570         Int_t pdg = GetPDG(species);
571         Int_t trackPdg = 0; 
572         if(mcTrack){
573           trackPdg = TMath::Abs(mcTrack->Particle()->GetPdgCode());
574         }
575         if(pdg == trackPdg) mcFound = kTRUE;    
576       }
577       quantities[0] = pTPC;
578       Bool_t detFlagSet = kFALSE;
579       for(Int_t idet = 0; idet < 4; idet++){
580         Char_t histname[256], histnameMC[256];
581         sprintf(histname, "h%s_El_like_%s", detname[idet], specname);
582         sprintf(histnameMC, "h%s_El_like_MC_%s", detname[idet], specname);
583         switch(idet){
584           case kITS:  esdTrack->GetITSpid(pidProbs);
585                       detFlagSet = status & AliESDtrack::kITSpid;
586                       break;
587           case kTPC:  esdTrack->GetTPCpid(pidProbs);
588                       detFlagSet = status & AliESDtrack::kTPCpid;
589                       break;
590           case kTRD:  esdTrack->GetTRDpid(pidProbs);
591                       detFlagSet = status & AliESDtrack::kTRDpid;
592                       break;
593           case kTOF:  esdTrack->GetTOFpid(pidProbs);  
594                       detFlagSet = status & AliESDtrack::kTOFpid;
595                       break;
596         };
597         quantities[1] = pidProbs[AliPID::kElectron];
598         // in case of TRD require 6 PID tracklets
599         if(kTRD == idet && esdTrack->GetTRDntrackletsPID() != 6) continue;
600         if(detFlagSet){
601           fOutput->Fill(histname, quantities[0], quantities[1]);
602           if(mcFound)
603             fOutput->Fill(histnameMC, quantities[0], quantities[1]);
604         }
605       }
606     }//.. ESD
607     else{
608       //AOD
609     }//.. aod 
610   }//.. while tracks 
611 }
612 //__________________________________________
613 void AliHFEpidQA::FillPIDresponse(const TObjArray * const particles, Int_t species){
614   //
615   // Fill the PID response of different detectors to V0 daughter particles
616   //
617   Char_t hname[256] = "";
618   Char_t hname2[256] = "";
619   Char_t hname3[256] = "";
620    
621   const Char_t *typeName[5] = {"Electron", "PionK0", "PionL", "Kaon", "Proton"};
622   const Int_t typePID[5] = {0, 2, 2, 3, 4};
623
624   // PID THnSparse
625   // axes:
626   // 0) species,  1) momentum, 2) DCA xy, 3) DCA z 
627   // 4) ITS signal 
628   // 5) TPC Ncls 6) TPC signal 7) TPC nSigma, 
629   // 8) TRD Ntrk, 9) TRD Ncls, 10) TRD dEdx, 
630
631   Double_t data[12];
632   memset(data, -99, sizeof(Double_t) *12);
633
634   Int_t run = fEvent->GetRunNumber();
635     
636   AliVParticle *recTrack = NULL;
637   TIter trackIter(particles); 
638   while((recTrack = dynamic_cast<AliVParticle *>(trackIter()))){
639     memset(data, -99, sizeof(Double_t) *10);
640     // ESD
641     if(!TString(recTrack->IsA()->GetName()).CompareTo("AliESDtrack")){
642       // case ESD
643       AliESDtrack *esdTrack = dynamic_cast<AliESDtrack *>(recTrack);
644       if(!esdTrack) continue;
645       
646       // for the PID THnSparse
647       data[0] = species;
648       data[1] = esdTrack->P();
649       Float_t impactR = -1.;
650       Float_t impactZ = -1.;
651       esdTrack->GetImpactParameters(impactR, impactZ);
652       data[2] = impactR;
653       data[3] = impactZ;
654       data[11] = 0; // initialize the TOF pid cut on elecgrons to false
655       // use ONLY tracks with PID flag TRUE
656       ULong_t status = 0;
657       status = esdTrack->GetStatus();
658
659       //
660       // DEBUG 
661       //
662       
663       fOutput->Fill("hITS_kFlags", 5., species);
664       if(status & AliESDtrack::kITSin)    fOutput->Fill("hITS_kFlags", 1., species);
665       if(status & AliESDtrack::kITSout)   fOutput->Fill("hITS_kFlags", 2., species);
666       if(status & AliESDtrack::kITSrefit) fOutput->Fill("hITS_kFlags", 3., species);
667       if(status & AliESDtrack::kITSpid)   fOutput->Fill("hITS_kFlags", 4., species);
668       
669       fOutput->Fill("hTPC_kFlags", 5., species);
670       if(status & AliESDtrack::kTPCin)    fOutput->Fill("hTPC_kFlags", 1., species);
671       if(status & AliESDtrack::kTPCout)   fOutput->Fill("hTPC_kFlags", 2., species);
672       if(status & AliESDtrack::kTPCrefit) fOutput->Fill("hTPC_kFlags", 3., species);
673       if(status & AliESDtrack::kTPCpid)   fOutput->Fill("hTPC_kFlags", 4., species);
674
675       fOutput->Fill("hTRD_kFlags", 5., species);
676       if(status & AliESDtrack::kTRDin)    fOutput->Fill("hTRD_kFlags", 1., species);
677       if(status & AliESDtrack::kTRDout)   fOutput->Fill("hTRD_kFlags", 2., species);
678       if(status & AliESDtrack::kTRDrefit) fOutput->Fill("hTRD_kFlags", 3., species);
679       if(status & AliESDtrack::kTRDpid)   fOutput->Fill("hTRD_kFlags", 4., species);
680
681       fOutput->Fill("hTOF_kFlags", 5., species);
682       if(status & AliESDtrack::kTOFin)    fOutput->Fill("hTOF_kFlags", 1., species);
683       if(status & AliESDtrack::kTOFout)   fOutput->Fill("hTOF_kFlags", 2., species);
684       if(status & AliESDtrack::kTOFrefit) fOutput->Fill("hTOF_kFlags", 3., species);
685       if(status & AliESDtrack::kTOFpid)   fOutput->Fill("hTOF_kFlags", 4., species);
686
687       
688       //
689       // ITS - 
690       //
691       if(status & AliESDtrack::kITSpid){
692         Double_t p = esdTrack->P();
693
694         // ITS signal
695         //Double_t itsSignal = esdTrack->GetITSsignal();
696         
697         // ITS dEdx
698         Double_t dEdxSamples[4];
699         esdTrack->GetITSdEdxSamples(dEdxSamples);
700         Int_t nSamples = 0;
701         Double_t dEdxSum = 0.;
702         sprintf(hname, "hITS_dEdx_%s", typeName[species]);
703         for(Int_t i=0; i<4; ++i){
704           if(dEdxSamples[i] > 0){
705             nSamples++;
706             fOutput->Fill(hname, i+1, p, dEdxSamples[i]);
707             dEdxSum += dEdxSamples[i];
708           }
709         }
710         if(4 == nSamples)fOutput->Fill(hname, 0, p, dEdxSum);
711         fOutput->Fill("hITS_dEdx_nSamples", nSamples);
712
713         Double_t signal = esdTrack->GetITSsignal();
714         sprintf(hname, "hITS_Signal_%s", typeName[species]);
715         fOutput->Fill(hname, p, signal);
716         data[4] = signal;
717         
718         // ITS number of signas
719         Double_t nsigma = fESDpid->NumberOfSigmasITS(esdTrack,(AliPID::EParticleType)typePID[species]);
720         sprintf(hname, "hITS_nSigma_%s",  typeName[species]);
721         fOutput->Fill(hname, p, nsigma);
722         
723         // ITS PID response
724         Double_t itsPID[5] = {-1, -1, -1, -1, -1};
725         esdTrack->GetITSpid(itsPID);
726         Int_t ix = GetMaxPID(itsPID);
727         sprintf(hname, "hITS_PID_p_%s", typeName[species]);
728         fOutput->Fill(hname, p, ix);
729       }//.. kITSpid
730       
731       //
732       // TPC
733       //
734       if(status & AliESDtrack::kTPCpid){
735         // Make TPC clusters Plot
736         data[5] = esdTrack->GetTPCNcls();
737         FillTPCinfo(esdTrack, species);
738
739         Double_t p = esdTrack->GetInnerParam() ? esdTrack->GetInnerParam()->P() : esdTrack->P();
740         // TPC dEdx
741         Double_t dEdx = esdTrack->GetTPCsignal();
742         sprintf(hname, "hTPC_dEdx_%s", typeName[species]);
743         fOutput->Fill(hname, p, dEdx);
744         data[6] = dEdx;
745         
746         //TPC number of sigmas
747         Double_t nsigma = fESDpid->NumberOfSigmasTPC(esdTrack,(AliPID::EParticleType)typePID[species]);
748         sprintf(hname, "hTPC_nSigma_%s",  typeName[species]);
749         fOutput->Fill(hname, p, nsigma);
750         data[7] = nsigma;
751
752         // TPC PID response
753         sprintf(hname, "hTPC_PID_p_%s", typeName[species]);
754         Double_t tpcPID[5] = {-1, -1, -1, -1, -1};
755         esdTrack->GetTPCpid(tpcPID);
756         Int_t ix = GetMaxPID(tpcPID);
757         fOutput->Fill(hname, p, ix);    
758       }//.. kTPCpid
759       
760       //
761       // TRD
762       //
763       
764       if(status & AliESDtrack::kTRDpid){
765         Double_t p = esdTrack->GetOuterParam() ? esdTrack->GetOuterParam()->P() : esdTrack->P();
766         
767         // TRD number of tracklets
768         Int_t ntrk = esdTrack->GetTRDntrackletsPID();
769         sprintf(hname, "hTRD_trk_%s", typeName[species]);
770         fOutput->Fill(hname, p, ntrk);
771         data[8] = ntrk;
772
773         // TRD PID response
774         sprintf(hname, "hTRD_PID_p_%s", typeName[species]);
775         Double_t trdPID[5] = {-1., -1., -1., -1., -1};
776         esdTrack->GetTRDpid(trdPID);
777         Int_t ix = GetMaxPID(trdPID);
778         fOutput->Fill(hname, p, ix);
779         // TRD n clusters
780         Int_t ncls = esdTrack->GetTRDncls();
781         sprintf(hname, "hTRD_cls_%s", typeName[species]);
782         fOutput->Fill(hname, p, ncls);
783         data[9] = ncls;
784
785         // TRD - per tracklet - dEdx per, likelihood
786         sprintf(hname, "hTRD_Nslc_%s", typeName[species]);
787         sprintf(hname2, "hTRD_slc_%s", typeName[species]);
788         sprintf(hname3, "hTRD_dEdx_%s", typeName[species]);
789         Int_t nSlices = esdTrack->GetNumberOfTRDslices();
790         Double_t sumTot = 0.;
791         Int_t not0Tot = 0;
792         for(Int_t l=0; l< 6; ++l){
793           Double_t trkData[7] = {-1.,-1, -1, -1, -1, -1, -1};
794           trkData[0] = species;
795           trkData[1] = p;
796           trkData[2] = l;
797           trkData[3] = ntrk;
798           trkData[5] = ncls;        
799           Double_t sum = 0.;
800           Int_t not0 = 0;
801           for(Int_t s=0; s<nSlices; ++s){
802             Double_t slice = esdTrack->GetTRDslice(l, s);
803             sum += slice;
804             if(slice > 0){
805               not0 += 1;
806               fOutput->Fill(hname2, p, s);
807             }
808           }//..slices
809
810           trkData[4] = not0;
811           fOutput->Fill(hname, p, not0);
812           fOutput->Fill(hname3, p, sum);
813           if(sum > 0){
814             sumTot += sum;
815             not0Tot += 1;
816           }
817           // lkelihoods per layer
818           if(not0Tot > 0 && fNNref){
819             Double_t likelihoods[5] = {-1., -1., -1., -1., -1};
820             TRDlikeTracklet(l, esdTrack, likelihoods);
821             trkData[6] = likelihoods[0];
822             //printf(" -D: species: %i, P; %f : %f, s: %f\n", species, p, likelihoods[0], s);
823           }
824           if(not0Tot) fOutput->Fill("hTRDtracklets", trkData);
825         }//..layers
826         // average dEx per number of tracklets
827         if(0 < not0Tot)
828           data[10] = sumTot / not0Tot;
829       }//.. kTRDpid
830       
831       
832       //
833       // TOF
834       //
835       if(status & AliESDtrack::kTOFpid){
836         Double_t p = esdTrack->GetOuterParam() ? esdTrack->GetOuterParam()->P() : esdTrack->P();
837         Double_t t0 = fESDpid->GetTOFResponse().GetTimeZero();
838
839         //TOF beta
840         sprintf(hname, "hTOF_beta_%s", typeName[species]);
841         Float_t beta = TOFbeta(esdTrack);       
842         fOutput->Fill(hname, p, beta);
843         fOutput->Fill("hTOF_beta_all", p, beta);
844         // TOF nSigma
845         Double_t nsigma = fESDpid->NumberOfSigmasTOF(esdTrack,(AliPID::EParticleType)typePID[species], t0);
846         sprintf(hname, "hTOF_nSigma_%s",  typeName[species]);
847         fOutput->Fill(hname, p, nsigma);
848         if(beta > 0.97 && beta < 1.03){
849           data[11] = 1;
850         }
851         
852         // TOF PID response
853         sprintf(hname, "hTOF_PID_p_%s", typeName[species]);
854         Double_t tofPID[5] = {-1., -1., -1., -1., -1};
855         esdTrack->GetTOFpid(tofPID);
856         Int_t ix = GetMaxPID(tofPID);
857         fOutput->Fill(hname, p, ix);
858         
859         // time of flight QA
860         // - distribution of (time - t0 - pion_time)
861         Double_t times[5];
862         esdTrack->GetIntegratedTimes(times);
863         Double_t tItrackL = esdTrack->GetIntegratedLength();
864         Double_t tTOFsignal = esdTrack->GetTOFsignal();
865         Double_t dT = tTOFsignal - t0 - times[2];
866         fOutput->Fill("hTOF_qa_TmT0mT", run*1.0, dT);
867         fOutput->Fill("hTOF_qa_length", run*1.0, tItrackL);
868
869         
870       }//.. kTOFpid
871       // temporary - the PIDsparse needs rebuilding
872       //fOutput->Fill("PIDsparse", data);
873     }
874     // AOD - comming soon
875     else{
876       continue;
877     }
878   }// .. tracks in TObjArray
879   
880
881 }
882 //__________________________________________
883 void AliHFEpidQA:: CheckEvent(){
884   //
885   // check some event variables
886   //
887
888   // check the T0 as a function of run number (less than one bin per number
889   Double_t t0 = fESDpid->GetTOFResponse().GetTimeZero();
890   Int_t run = fEvent->GetRunNumber();
891   Double_t data[2] = {run*1.0, t0*1000.};
892   fOutput->Fill("hEvent_T0", data);
893   
894
895 }
896 //__________________________________________
897 TList *AliHFEpidQA::GetOutput(){
898   //
899   // Getter for Output histograms
900   //
901   return fOutput->GetList();
902 }
903
904 //__________________________________________
905 TList *AliHFEpidQA::GetV0pidQA(){
906   //
907   // Getter for V0 PID QA histograms
908   //
909   return fV0pid->GetListOfQAhistograms();
910 }
911
912 //__________________________________________
913 TList *AliHFEpidQA::GetV0pidMC(){
914   //
915   // Getter for V0 PID QA histograms
916   //
917   if(fV0pidMC)
918     return fV0pidMC->GetListOfQAhistograms();
919   return NULL;
920 }
921
922 //__________________________________________
923 void AliHFEpidQA::RecalculateTRDpid(AliESDtrack * /*track*/, Double_t * /*pidProbs*/) const{
924   //  fTRDpidResponse->MakePID(track);
925   //  track->GetTRDpid(pidProbs);
926 }
927
928 //__________________________________________
929 void AliHFEpidQA::RecalculateTRDpid(AliAODTrack * /*track*/, Double_t * /*pidProbs*/) const{
930   //  fTRDpidResponse->MakePID(track, pidProbs);
931 }
932 //___________________________________________________________________
933 Float_t AliHFEpidQA::TOFbeta(AliESDtrack * const track) const {
934   // computes the TOF beta
935   Double_t l = track->GetIntegratedLength();  // cm
936   Double_t t = track->GetTOFsignal();
937   Double_t t0 = fESDpid->GetTOFResponse().GetTimeZero(); // ps
938
939   //printf("-D: l: %f, t: %f, t0: %f\n", l, t, t0);
940
941   if(l < 360. || l > 800.) return 0.;
942   if(t <= 0.) return 0.;
943   if(t0 >999990.0) return 0.;
944   
945
946   t -= t0; // subtract the T0
947
948   l *= 0.01;  // cm ->m
949   t *= 1e-12; //ps -> s
950
951   
952   Double_t v = l / t;
953   Float_t beta = v / TMath::C();
954
955   return beta;
956 }
957 //____________________________________________
958 Int_t AliHFEpidQA::GetMaxPID(const Double_t *pidProbs) const {
959   //
960   // return the index of maximal PID probability
961   //
962   Int_t ix = -1;
963   Double_t tmp = 0.2;
964   for(Int_t i=0; i<5; ++i){
965     if(pidProbs[i] > tmp){
966       ix = i;
967       tmp = pidProbs[i];
968     }
969   }
970   return ix;
971 }
972 //_____________________________________________
973 Int_t AliHFEpidQA::GetPDG(Int_t species){
974   //
975   // return the PDG particle code
976   //
977
978   Int_t pdg = 0;
979
980   switch(species){
981   case  AliHFEV0cuts::kRecoElectron:
982     pdg = TMath::Abs(kElectron);
983     break;
984   case  AliHFEV0cuts::kRecoPionK0:
985     pdg = TMath::Abs(kPiPlus);
986     break;
987   case  AliHFEV0cuts::kRecoPionL:
988     pdg = TMath::Abs(kPiPlus);
989     break;
990   case  AliHFEV0cuts::kRecoProton:
991     pdg = TMath::Abs(kProton);
992     break;
993   default:  // non investigated species
994     AliDebug(3, "Species not recognised");
995     return 0;
996   }  
997
998   return pdg;
999
1000 }
1001
1002 //_____________________________________________
1003 TObjArray * AliHFEpidQA::MakeTrackList(const TObjArray *tracks) const {
1004   //
1005   // convert list of AliHFEV0Info into a list of AliVParticle
1006   //
1007   TObjArray *output = new TObjArray;
1008   TIter trackInfos(tracks);
1009   AliHFEV0info *trackInfo = NULL;
1010   while((trackInfo = dynamic_cast<AliHFEV0info *>(trackInfos())))
1011     output->Add(trackInfo->GetTrack());
1012
1013   return output;
1014 }
1015
1016 //_____________________________________________
1017 TObjArray * AliHFEpidQA::MakeCleanListElectrons(const TObjArray *electrons) const {
1018   //
1019   // Cleanup electron sample using TPC PID
1020   // PID requirement will allways be implemented to the pair
1021   // Strategy
1022   //
1023   TObjArray *tracks = new TObjArray;
1024   TIter candidates(electrons);
1025   AliESDEvent *esd; AliAODEvent *aod;
1026   AliHFEV0info *hfetrack;
1027   // const Double_t kSigmaTight = 1.;
1028   // const Double_t kSigmaLoose = 4.;
1029   const Double_t kSigmaTight = 2.;
1030   const Double_t kSigmaLoose = 2.;
1031   const Double_t shift = 0.57;
1032   if((esd = dynamic_cast<AliESDEvent *>(fEvent))){
1033     AliESDtrack *track = NULL, *partnerTrack = NULL;
1034     while((hfetrack = dynamic_cast<AliHFEV0info *>(candidates()))){
1035       track = dynamic_cast<AliESDtrack *>(hfetrack->GetTrack());
1036       partnerTrack = esd->GetTrack(hfetrack->GetPartnerID());
1037       Double_t nSigmaTrack = TMath::Abs(fESDpid->NumberOfSigmasTPC(track, AliPID::kElectron) - shift);
1038       Double_t nSigmaPartner = TMath::Abs(fESDpid->NumberOfSigmasTPC(partnerTrack, AliPID::kElectron) - shift);
1039       if((nSigmaTrack < kSigmaTight && nSigmaPartner < kSigmaLoose) || (nSigmaTrack < kSigmaLoose && nSigmaPartner < kSigmaTight))
1040         tracks->Add(track);
1041     }
1042   } else {
1043     aod = dynamic_cast<AliAODEvent *>(fEvent);
1044     AliAODTrack *track = NULL, *partnerTrack = NULL;
1045     while((hfetrack = dynamic_cast<AliHFEV0info *>(candidates()))){
1046       track = dynamic_cast<AliAODTrack *>(hfetrack->GetTrack());
1047       partnerTrack = aod->GetTrack(hfetrack->GetPartnerID());
1048       // will be coming soon
1049     }
1050   }
1051   return tracks;
1052 }
1053 //___________________________________________________________
1054 void AliHFEpidQA::CheckTenderV0pid(const TObjArray * const particles, Int_t species){
1055
1056   //
1057   // retrieve the status bits from the TObject used to flag
1058   // the V0 daughter tracks and return the found PID
1059   // 0 - electron, 1 - pion, 2 - proton
1060   // 
1061
1062   const Int_t id[5] = {0, 1, 1, -1, 2}; // convert AliHFEpid to simple 0, 1, 2
1063
1064   AliVParticle *recTrack = NULL;
1065   TIter trackIter(particles); 
1066   while((recTrack = dynamic_cast<AliVParticle *>(trackIter()))){
1067     if(!TString(recTrack->IsA()->GetName()).CompareTo("AliESDtrack")){
1068       // case ESD
1069       AliESDtrack *track = dynamic_cast<AliESDtrack *>(recTrack);
1070       if(!track) continue;
1071       Int_t tPID = GetTenderV0pid(track);
1072       fOutput->Fill("h_tender_check_01", id[species]*1.0, tPID);
1073       fOutput->Fill("h_tender_check_02", id[species], id[species]*1.0, tPID);
1074       
1075     } //.. case ESD      
1076     
1077   }//.. iterate of tracks  
1078 }
1079 //___________________________________________________________
1080 Int_t AliHFEpidQA::GetTenderV0pid(AliESDtrack * const track){
1081   //
1082   // retrieve the PID nformation stored in the status flags by the train tender
1083   //
1084   
1085
1086   Int_t pid = -1;
1087   if(!track){
1088     return pid;
1089   }
1090   
1091   Int_t nTimes = 0;
1092   
1093   if(track->TestBit(2<<14)){
1094     pid = 0;
1095     nTimes++;
1096   }
1097   if(track->TestBit(2<<15)){
1098     pid = 1;
1099     nTimes++;
1100   }
1101   if(track->TestBit(2<<16)){
1102     pid = 2;
1103     nTimes++;
1104   }
1105   
1106   if(nTimes > 1){
1107     AliWarning("V0 track labeled multiple times by the V0 tender");
1108     pid = -1;
1109   }
1110     
1111   return pid;
1112
1113 }
1114 //___________________________________________________________
1115 Double_t AliHFEpidQA::TRDlikeTracklet(Int_t layer, AliESDtrack * const track, Double_t *likelihood){
1116   //
1117   // compute the TRD electron likelihoods for 1 tracklet
1118   // based on teh AliTRDpidRecalculator in train/until/tender
1119   // returns sum of the likelihoods (which should be 1)
1120   //
1121
1122   const Double_t cScaleGain = 1./ 16000.;
1123   const Float_t pBins[11] ={0.6, 0.8, 1.0, 1.5, 2.0, 3.0, 4.0, 5.0, 6.0, 8.0, 10.0}; // momentum bins
1124
1125   if(!track) return kFALSE;
1126   Float_t p = track->GetTRDmomentum(layer); // momentum for a tracklet in the ESDtrack
1127   if(p < 0) return kFALSE;
1128
1129   Int_t mombin = TRDmomBin(p);  // momentum bin
1130   Float_t dEdxTRDsum = 0;              // dEdxTRDsum for checking if tracklet is available
1131   Float_t dEdxTRD[8];           // dEdx for a tracklet in the ESD slices
1132   Double_t ddEdxTRD[8];         // dEdx as Double_t for TMultiLayerPerceptron::Evaluate()
1133   
1134   Double_t prob[AliPID::kSPECIES];  // probabilities for all species in all layers
1135
1136   for(Int_t is = 0; is < AliPID::kSPECIES; is++){
1137     likelihood[is] = 0.2;               // init probabilities
1138     prob[is] = 0.2; 
1139   }
1140
1141   Double_t sum = 0.;
1142
1143   for(Int_t islice = 0; islice<8; islice++){
1144     dEdxTRD[islice]=0.;                                    // init dE/dx
1145     ddEdxTRD[islice]=0.;                                   // init dE/dx
1146     dEdxTRD[islice]=track->GetTRDslice(layer,islice);       // get the dE/dx values
1147     dEdxTRDsum += dEdxTRD[islice];
1148     ddEdxTRD[islice]=(Double_t)dEdxTRD[islice]*cScaleGain;            // rescale dE/dx
1149     
1150   }
1151   for(Int_t is = 0; is < AliPID::kSPECIES; is++){
1152     Double_t probloc1, probloc2; 
1153     if(mombin == 0 && mombin < pBins[0]){          //  calculate PID for p > 0.6 GeV/c
1154       prob[is] = fNet[mombin]->Evaluate(is, ddEdxTRD);
1155     }
1156     else if(mombin == 10 && mombin >= pBins[10]){   //  calculate PID for p >= 10.0 GeV/c
1157       prob[is] = fNet[mombin]->Evaluate(is, ddEdxTRD);
1158     }
1159     else{                                                    //  standard calculation
1160       Int_t mombin1 = 0, mombin2 = 0;             // lower and upper momentum bin      
1161       if(p < pBins[mombin]) {mombin1 = mombin -1; mombin2 = mombin;}
1162       if(p >= pBins[mombin]) {mombin1 = mombin; mombin2 = mombin+1;}
1163       probloc1 = fNet[mombin1]->Evaluate(is, ddEdxTRD);
1164       probloc2 = fNet[mombin2]->Evaluate(is, ddEdxTRD);
1165       // weighting of the probability with regard to the track momentum
1166       prob[is] = probloc1 + (probloc2-probloc1)*(p-pBins[mombin1])/(pBins[mombin2]-pBins[mombin1]);
1167     }
1168     likelihood[is] = prob[is];
1169     sum += likelihood[is];
1170   }     
1171   
1172   return sum;
1173 }
1174 //__________________________________________________________________________
1175 Int_t AliHFEpidQA::TRDmomBin(Double_t p) const {
1176   //
1177   // compute the momentum bin position
1178   // 
1179
1180   const Float_t pBins[11] ={0.6, 0.8, 1.0, 1.5, 2.0, 3.0, 4.0, 5.0, 6.0, 8.0, 10.0}; // momentum bins
1181   
1182   Int_t pBin1 = -1;                                    // check bin1
1183   Int_t pBin2 = -1;                                    // check bin2
1184
1185   if(p < 0) return -1;                                 // return -1 if momentum < 0
1186   if(p < pBins[0]) return 0;                                      // smallest momentum bin
1187   if(p >= pBins[10]) return 10; // largest momentum bin
1188
1189
1190   // calculate momentum bin for non extremal momenta
1191   for(Int_t iMomBin = 1; iMomBin < 11; iMomBin++){
1192     if(p < pBins[iMomBin]){
1193       pBin1 = iMomBin - 1;
1194       pBin2 = iMomBin;
1195     }
1196     else
1197       continue;
1198
1199     if(p - pBins[pBin1] >= pBins[pBin2] - p){
1200       return pBin2;
1201     }
1202     else{
1203       return pBin1;
1204     }
1205   }
1206
1207   return -1;
1208
1209
1210 }
1211 //__________________________________________________________________________
1212
1213