]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/vertexingHF/AliAnalysisTaskSEHFQA.cxx
Coverity
[u/mrichter/AliRoot.git] / PWG3 / vertexingHF / AliAnalysisTaskSEHFQA.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2009, 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 /////////////////////////////////////////////////////////////
17 //
18 // AliAnalysisTaskSE for HF quality assurance
19 //
20 // Author: Chiara Bianchin, chiara.bianchin@pd.infn.it
21 /////////////////////////////////////////////////////////////
22
23 #include <Riostream.h>
24 #include <TClonesArray.h>
25 #include <TCanvas.h>
26 #include <TList.h>
27 #include <TH1F.h>
28 #include <TH2F.h>
29 #include <TDatabasePDG.h>
30
31 #include <AliAnalysisDataSlot.h>
32 #include <AliAnalysisDataContainer.h>
33 #include "AliAnalysisManager.h"
34 #include "AliESDtrack.h"
35 #include "AliVertexerTracks.h"
36 #include "AliPID.h"
37 #include "AliTPCPIDResponse.h"
38 #include "AliAODHandler.h"
39 #include "AliAODEvent.h"
40 #include "AliAODVertex.h"
41 #include "AliAODTrack.h"
42 #include "AliAODMCParticle.h"
43 #include "AliAODMCHeader.h"
44 #include "AliAODRecoDecayHF2Prong.h"
45 #include "AliAODRecoCascadeHF.h"
46 #include "AliAnalysisVertexingHF.h"
47 #include "AliAnalysisTaskSE.h"
48 #include "AliCounterCollection.h"
49 #include "AliRDHFCuts.h"
50 #include "AliRDHFCutsDplustoKpipi.h"
51 #include "AliRDHFCutsD0toKpipipi.h"
52 #include "AliRDHFCutsDstoKKpi.h"
53 #include "AliRDHFCutsDStartoKpipi.h"
54 #include "AliRDHFCutsD0toKpi.h"
55 #include "AliRDHFCutsLctopKpi.h"
56
57 #include "AliAnalysisTaskSEHFQA.h"
58
59
60 ClassImp(AliAnalysisTaskSEHFQA)
61
62 //____________________________________________________________________________
63
64 AliAnalysisTaskSEHFQA::AliAnalysisTaskSEHFQA():AliAnalysisTaskSE(),
65   fNEntries(0x0),
66   fOutputPID(0x0),
67   fOutputTrack(0x0),
68   fOutputCounters(0x0),
69   fOutputCheckCentrality(0x0),
70   fDecayChannel(AliAnalysisTaskSEHFQA::kD0toKpi),
71   fCuts(0x0),
72   fEstimator(AliRDHFCuts::kCentTRK),
73   fReadMC(kFALSE)
74 {
75   //default constructor
76 }
77
78 //____________________________________________________________________________
79 AliAnalysisTaskSEHFQA::AliAnalysisTaskSEHFQA(const char *name, AliAnalysisTaskSEHFQA::DecChannel ch,AliRDHFCuts* cuts):
80   AliAnalysisTaskSE(name),
81   fNEntries(0x0),
82   fOutputPID(0x0),
83   fOutputTrack(0x0),
84   fOutputCounters(0x0),
85   fOutputCheckCentrality(0x0),
86   fDecayChannel(ch),
87   fCuts(0x0),
88   fEstimator(AliRDHFCuts::kCentTRK),
89   fReadMC(kFALSE)
90 {
91   //constructor
92
93   //SetCutObject(cuts);
94   fCuts=cuts;
95
96   // Output slot #1 writes into a TH1F container (number of events)
97   DefineOutput(1,TH1F::Class());  //My private output
98   // Output slot #2 writes into a TList container (PID)
99   DefineOutput(2,TList::Class());  //My private output
100   // Output slot #3 writes into a TList container (Tracks)
101   DefineOutput(3,TList::Class());  //My private output
102   // Output slot #4 writes into a AliRDHFCuts container (cuts)
103   switch(fDecayChannel){
104   case 0:
105     DefineOutput(4,AliRDHFCutsDplustoKpipi::Class());  //My private output
106     break;
107   case 1:
108     DefineOutput(4,AliRDHFCutsD0toKpi::Class());  //My private output
109     break;
110   case 2:
111     DefineOutput(4,AliRDHFCutsDStartoKpipi::Class());  //My private output
112     break;
113   case 3:
114     DefineOutput(4,AliRDHFCutsDstoKKpi::Class());  //My private output
115     break;
116   case 4:
117     DefineOutput(4,AliRDHFCutsD0toKpipipi::Class());  //My private output
118     break;
119   case 5:
120     DefineOutput(4,AliRDHFCutsLctopKpi::Class());  //My private output
121     break;
122   }
123   // Output slot #5 writes into a TList container (AliCounterCollection)
124   DefineOutput(5,TList::Class());  //My private output
125   // Output slot #6 writes into a TList container (TH1F)
126   DefineOutput(6,TList::Class());  //My private output
127 }
128
129 //___________________________________________________________________________
130 AliAnalysisTaskSEHFQA::~AliAnalysisTaskSEHFQA()
131 {
132   //destructor
133   delete fNEntries;
134
135   delete fOutputPID;
136
137   delete fOutputTrack;
138
139   delete fOutputCounters;
140
141   delete fOutputCheckCentrality;
142
143 }
144
145 //___________________________________________________________________________
146 void AliAnalysisTaskSEHFQA::Init(){
147
148   //initialization
149   if(fDebug > 1) printf("AnalysisTaskSEHFQA::Init() \n");
150
151   switch(fDecayChannel){
152   case 0:
153     {
154       AliRDHFCutsDplustoKpipi* copycut=new AliRDHFCutsDplustoKpipi(*(static_cast<AliRDHFCutsDplustoKpipi*>(fCuts)));
155       // Post the data
156       PostData(4,copycut);
157     }
158     break;
159   case 1:
160     {
161       AliRDHFCutsD0toKpi* copycut=new AliRDHFCutsD0toKpi(*(static_cast<AliRDHFCutsD0toKpi*>(fCuts)));
162       // Post the data
163       PostData(4,copycut);
164     }
165     break;
166   case 2:
167     {
168       AliRDHFCutsDStartoKpipi* copycut=new AliRDHFCutsDStartoKpipi(*(static_cast<AliRDHFCutsDStartoKpipi*>(fCuts)));
169       // Post the data
170       PostData(4,copycut);
171     }
172     break;
173   case 3:
174     {
175       AliRDHFCutsDstoKKpi* copycut=new AliRDHFCutsDstoKKpi(*(static_cast<AliRDHFCutsDstoKKpi*>(fCuts)));
176       // Post the data
177       PostData(4,copycut);
178     }
179     break;
180   case 4:
181     {
182       AliRDHFCutsD0toKpipipi* copycut=new AliRDHFCutsD0toKpipipi(*(static_cast<AliRDHFCutsD0toKpipipi*>(fCuts)));
183       // Post the data
184       PostData(4,copycut);
185     }
186     break;
187   case 5:
188     {
189       AliRDHFCutsLctopKpi* copycut=new AliRDHFCutsLctopKpi(*(static_cast<AliRDHFCutsLctopKpi*>(fCuts)));
190       // Post the data
191       PostData(4,copycut);
192     }
193     break;
194
195   default:
196     return;
197   }
198
199
200
201 }
202
203 //___________________________________________________________________________
204 void AliAnalysisTaskSEHFQA::UserCreateOutputObjects()
205 {
206
207   //create the output container
208   if(fDebug > 1) printf("AnalysisTaskSEHFQA::UserCreateOutputObjects() \n");
209
210   //count events
211
212   fNEntries=new TH1F(GetOutputSlot(1)->GetContainer()->GetName(), "Counts the number of events", 9,-0.5,8.5);
213   fNEntries->GetXaxis()->SetBinLabel(1,"nEventsAnal");
214   fNEntries->GetXaxis()->SetBinLabel(2,"Pile-up Rej");
215   fNEntries->GetXaxis()->SetBinLabel(3,"No VertexingHF");
216   fNEntries->GetXaxis()->SetBinLabel(4,"nCandidates(AnCuts)");
217   fNEntries->GetXaxis()->SetBinLabel(5,"EventsWithGoodVtx");
218   fNEntries->GetXaxis()->SetBinLabel(6,"N. of 0SMH");
219   if(fReadMC){
220     fNEntries->GetXaxis()->SetBinLabel(7,"MC Cand from c");
221     fNEntries->GetXaxis()->SetBinLabel(8,"MC Cand from b");
222     fNEntries->GetXaxis()->SetBinLabel(9,"N fakes Trks");
223   } else{
224     fNEntries->GetXaxis()->SetBinLabel(7,"N candidates");
225   }
226
227   fNEntries->GetXaxis()->SetNdivisions(1,kFALSE);
228
229   //PID
230
231   fOutputPID=new TList();
232   fOutputPID->SetOwner();
233   fOutputPID->SetName(GetOutputSlot(2)->GetContainer()->GetName());
234
235   //TOF pid
236   TString hname="hTOFsig";
237   TH1F* hTOFsig=new TH1F(hname.Data(),"Distribution of TOF signal;TOF time [ps];Entries", 100, -2.e3,40.e3);
238
239   hname="hTOFtime";
240   TH1F* hTOFtime=new TH1F(hname.Data(),"Distribution of TOF time Kaon;TOF time(Kaon) [ps];Entries", 1000, 0.,50000.);
241
242   hname="hTOFtimeKaonHyptime";
243   TH2F* hTOFtimeKaonHyptime=new TH2F(hname.Data(),"TOFtime - timeHypothesisForKaon;p[GeV/c];TOFtime - timeHypothesisForKaon [ps]",200,0.,4.,1000,-20000.,20000.);
244
245   hname="hTOFtimeKaonHyptimeAC";
246   TH2F* hTOFtimeKaonHyptimeAC=new TH2F(hname.Data(),"TOFtime - timeHypothesisForKaon;p[GeV/c];TOFtime - timeHypothesisForKaon [ps]",200,0.,4.,1000,-20000.,20000.);
247
248   hname="hTOFsigmaK160";
249   TH2F* hTOFsigmaK160=new TH2F(hname.Data(),"(TOFsignal-timeK)/160 ps;p[GeV/c];(TOFsignal-timeK)/160 ps",200,0.,4.,100,-5,5);
250
251   hname="hTOFsigmaPion160";
252   TH2F* hTOFsigmaPion160=new TH2F(hname.Data(),"(TOFsignal-time#pi)/160 ps;p[GeV/c];(TOFsignal-time#pi)/160 ps",200,0.,4.,100,-5,5);
253
254   hname="hTOFsigmaProton160";
255   TH2F* hTOFsigmaProton160=new TH2F(hname.Data(),"(TOFsignal-timep)/160 ps;p[GeV/c];(TOFsignal-time p)/160 ps",200,0.,4.,100,-5,5);
256
257   hname="hTOFsigmaKSigPid";
258   TH2F* hTOFsigmaKSigPid=new TH2F(hname.Data(),"(TOFsignal-timeK)/tofSigPid;p[GeV/c];(TOFsignal-timeK)/tofSigPid",200,0.,4.,100,-5,5);
259
260   hname="hTOFsigmaPionSigPid";
261   TH2F* hTOFsigmaPionSigPid=new TH2F(hname.Data(),"(TOFsignal-time#pi)/tofSigPid;p[GeV/c];(TOFsignal-time#pi)/tofSigPid",200,0.,4.,100,-5,5);
262
263   hname="hTOFsigmaProtonSigPid";
264   TH2F* hTOFsigmaProtonSigPid=new TH2F(hname.Data(),"(TOFsignal-timep)/tofSigPid;p[GeV/c];(TOFsignal-time p)/tofSigPid",200,0.,4.,100,-5,5);
265
266   hname="hTOFsigPid3sigPion";
267   TH1F* hTOFsigPid3sigPion=new TH1F(hname.Data(),"TOF PID resolution (#pi) [ps]",500,0.,1000.);
268
269   hname="hTOFsigPid3sigKaon";
270   TH1F* hTOFsigPid3sigKaon=new TH1F(hname.Data(),"TOF PID resolution (K) [ps]",500,0.,1000.);
271
272   hname="hTOFsigPid3sigProton";
273   TH1F* hTOFsigPid3sigProton=new TH1F(hname.Data(),"TOF PID resolution (p) [ps]",500,0.,1000.);
274
275
276   //TPC pid
277   hname="hTPCsig";
278   TH1F* hTPCsig=new TH1F(hname.Data(),"Distribution of TPC signal;TPC sig;Entries", 100, 35.,100.);
279
280   hname="hTPCsigvsp";
281   TH2F* hTPCsigvsp=new TH2F(hname.Data(),"TPCsig vs p;TPC p[GeV/c];TPCsig",200,0.,4.,1000,35.,100.);
282  
283   hname="hTPCsigvspAC";
284   TH2F* hTPCsigvspAC=new TH2F(hname.Data(),"TPCsig vs p;TPCp[GeV/c];TPCsig",200,0.,4.,1000,35.,100.);
285
286   hname="hTPCsigmaK";
287   TH2F* hTPCsigmaK=new TH2F(hname.Data(),"TPC Sigma for K as a function of momentum;p[GeV/c];Sigma Kaon",200,0.,4.,200,-5,5);
288
289   hname="hTPCsigmaPion";
290   TH2F* hTPCsigmaPion=new TH2F(hname.Data(),"TPC Sigma for #pi as a function of momentum;p[GeV/c];Sigma #pi",200,0.,4.,200,-5,5);
291
292   hname="hTPCsigmaProton";
293   TH2F* hTPCsigmaProton=new TH2F(hname.Data(),"TPC Sigma for proton as a function of momentum;p[GeV/c];Sigma Proton",200,0.,4.,200,-5,5);
294
295   fOutputPID->Add(hTOFsig);
296   fOutputPID->Add(hTPCsig);
297   fOutputPID->Add(hTOFtime);
298   fOutputPID->Add(hTOFtimeKaonHyptime);
299   fOutputPID->Add(hTOFtimeKaonHyptimeAC);
300   fOutputPID->Add(hTOFsigmaK160);
301   fOutputPID->Add(hTOFsigmaPion160);
302   fOutputPID->Add(hTOFsigmaProton160);
303   fOutputPID->Add(hTOFsigmaKSigPid);
304   fOutputPID->Add(hTOFsigmaPionSigPid);
305   fOutputPID->Add(hTOFsigmaProtonSigPid);
306   fOutputPID->Add(hTOFsigPid3sigPion);
307   fOutputPID->Add(hTOFsigPid3sigKaon);
308   fOutputPID->Add(hTOFsigPid3sigProton);
309  fOutputPID->Add(hTPCsigvsp);
310   fOutputPID->Add(hTPCsigvspAC);
311   fOutputPID->Add(hTPCsigmaK);
312   fOutputPID->Add(hTPCsigmaPion);
313   fOutputPID->Add(hTPCsigmaProton);
314
315   //quality of the tracks
316
317   fOutputTrack=new TList();
318   fOutputTrack->SetOwner();
319   fOutputTrack->SetName(GetOutputSlot(3)->GetContainer()->GetName());
320
321   hname="hnClsITS";
322   TH1F* hnClsITS=new TH1F(hname.Data(),"Distribution of number of ITS clusters;nITScls;Entries",7,-0.5,6.5);
323
324   hname="hnClsITS-SA";
325   TH1F* hnClsITSSA=new TH1F(hname.Data(),"Distribution of number of ITS clusters(ITS-SA);nITScls;Entries",7,-0.5,6.5);
326
327   hname="hnClsSPD";
328   TH1F* hnClsSPD=new TH1F(hname.Data(),"Distribution of number of SPD clusters;nSPDcls;Entries",3,-0.5,2.5);
329
330   hname="hptGoodTr";
331   TH1F* hptGoodTr=new TH1F(hname.Data(),"Pt distribution of 'good' tracks;p_{t}[GeV];Entries/0.05 GeV/c",400,0.,20.);
332   hptGoodTr->SetTitleOffset(1.3,"Y");
333
334   hname="hdistrGoodTr";
335   TH1F* hdistrGoodTr=new TH1F(hname.Data(),"Distribution of number of 'good' tracks per event;no.good-tracks/ev;Entries",4000,-0.5,3999.5);
336   hdistrGoodTr->SetTitleOffset(1.3,"Y");
337
338   hname="hd0";
339   TH1F* hd0=new TH1F(hname.Data(),"Impact parameter distribution of 'good' tracks;d_{0}[cm];Entries/10^{3} cm",200,-0.1,0.1);
340
341   fOutputTrack->Add(hnClsITS);
342   fOutputTrack->Add(hnClsITSSA);
343   fOutputTrack->Add(hnClsSPD);
344   fOutputTrack->Add(hptGoodTr);
345   fOutputTrack->Add(hdistrGoodTr);
346   fOutputTrack->Add(hd0);
347   
348   if(fCuts->GetUseCentrality()){
349
350     //Centrality (Counters)
351     fOutputCounters=new TList();
352     fOutputCounters->SetOwner();
353     fOutputCounters->SetName(GetOutputSlot(5)->GetContainer()->GetName());
354
355     AliCounterCollection *stdEstimator=new AliCounterCollection("stdEstimator");
356     stdEstimator->AddRubric("run",500000);
357     stdEstimator->AddRubric("centralityclass","0_10/10_20/20_30/30_40/40_50/50_60/60_70/70_80/80_90/90_100");
358     stdEstimator->Init();
359     AliCounterCollection *secondEstimator=new AliCounterCollection("secondEstimator");
360     secondEstimator->AddRubric("run",500000);
361     secondEstimator->AddRubric("centralityclass","0_10/10_20/20_30/30_40/40_50/50_60/60_70/70_80/80_90/90_100/-990_-980");
362     secondEstimator->Init();
363
364     fOutputCounters->Add(stdEstimator);
365     fOutputCounters->Add(secondEstimator);
366
367     //Centrality (Checks)
368     fOutputCheckCentrality=new TList();
369     fOutputCheckCentrality->SetOwner();
370     fOutputCheckCentrality->SetName(GetOutputSlot(6)->GetContainer()->GetName());
371
372     hname="hNtrackletsIn";
373     TH1F* hNtrackletsIn=new TH1F(hname.Data(),"Number of tracklets in Centrality range;ntracklets;Entries",5000,-0.5,4999.5);
374
375     hname="hMultIn";
376     TH1F* hMultIn=new TH1F(hname.Data(),"Multiplicity;multiplicity in Centrality range;Entries",10000,-0.5,9999.5);
377
378     hname="hNtrackletsOut";
379     TH1F* hNtrackletsOut=new TH1F(hname.Data(),"Number of tracklets out of Centrality range;ntracklets;Entries",5000,-0.5,4999.5);
380
381     hname="hMultOut";
382     TH1F* hMultOut=new TH1F(hname.Data(),"Multiplicity out of Centrality range;multiplicity;Entries",10000,-0.5,9999.5);
383
384     fOutputCheckCentrality->Add(hNtrackletsIn);
385     fOutputCheckCentrality->Add(hNtrackletsOut);
386     fOutputCheckCentrality->Add(hMultIn);
387     fOutputCheckCentrality->Add(hMultOut);
388   } else{
389     hname="hNtracklets";
390     TH1F* hNtracklets=new TH1F(hname.Data(),"Number of tracklets;ntracklets;Entries",5000,-0.5,4999.5);
391
392     hname="hMult";
393     TH1F* hMult=new TH1F(hname.Data(),"Multiplicity;multiplicity;Entries",10000,-0.5,9999.5);
394     fOutputTrack->Add(hNtracklets);
395     fOutputTrack->Add(hMult);
396   }
397
398   // Post the data
399   PostData(1,fNEntries);
400   PostData(2,fOutputPID);
401   PostData(3,fOutputTrack);
402   PostData(4,fCuts);
403   PostData(5,fOutputCounters);
404   PostData(6,fOutputCheckCentrality);
405 }
406
407 //___________________________________________________________________________
408 void AliAnalysisTaskSEHFQA::UserExec(Option_t */*option*/)
409 {
410   // Execute analysis for current event
411
412   AliAODEvent *aod = dynamic_cast<AliAODEvent*> (InputEvent());
413   if(fDebug>2) printf("Analysing decay %d\n",fDecayChannel);
414   // Post the data already here
415   PostData(1,fNEntries);
416   PostData(2,fOutputPID);
417   PostData(3,fOutputTrack);
418   PostData(4,fCuts);
419   PostData(5,fOutputCounters);
420   PostData(6,fOutputCheckCentrality);
421
422   TClonesArray *arrayProng =0;
423   Int_t pdg=0;
424   Int_t *pdgdaughters=0x0;
425
426   if(!aod && AODEvent() && IsStandardAOD()) { 
427     // In case there is an AOD handler writing a standard AOD, use the AOD 
428     // event in memory rather than the input (ESD) event.    
429     aod = dynamic_cast<AliAODEvent*> (AODEvent());
430     // in this case the braches in the deltaAOD (AliAOD.VertexingHF.root)
431     // have to taken from the AOD event hold by the AliAODExtension
432     AliAODHandler* aodHandler = (AliAODHandler*) 
433       ((AliAnalysisManager::GetAnalysisManager())->GetOutputEventHandler());
434     if(aodHandler->GetExtensions()) {
435       
436       AliAODExtension *ext = (AliAODExtension*)aodHandler->GetExtensions()->FindObject("AliAOD.VertexingHF.root");
437       AliAODEvent *aodFromExt = ext->GetAOD();
438    
439    
440       
441       switch(fDecayChannel){
442       case 0:
443         arrayProng=(TClonesArray*)aodFromExt->GetList()->FindObject("Charm3Prong");
444         pdg=411;
445         if(fReadMC){
446           pdgdaughters =new Int_t[3];
447           pdgdaughters[0]=211;//pi
448           pdgdaughters[1]=321;//K
449           pdgdaughters[2]=211;//pi
450         }
451         break; 
452       case 1:
453         arrayProng=(TClonesArray*)aodFromExt->GetList()->FindObject("D0toKpi");
454         pdg=421;
455         if(fReadMC){
456           pdgdaughters =new Int_t[2];
457           pdgdaughters[0]=211;//pi 
458           pdgdaughters[1]=321;//K
459         }
460         break; 
461       case 2:
462         arrayProng=(TClonesArray*)aodFromExt->GetList()->FindObject("Dstar");
463         pdg=413;
464         if(fReadMC){
465           pdgdaughters =new Int_t[3];
466           pdgdaughters[1]=211;//pi
467           pdgdaughters[0]=321;//K
468           pdgdaughters[2]=211;//pi (soft?)
469         }
470         break; 
471       case 3:
472         arrayProng=(TClonesArray*)aodFromExt->GetList()->FindObject("Charm3Prong");
473         pdg=431;
474         if(fReadMC){
475           pdgdaughters =new Int_t[3];
476           pdgdaughters[0]=321;//K
477           pdgdaughters[1]=321;//K
478           pdgdaughters[2]=211;//pi
479         }
480         break; 
481       case 4:
482         arrayProng=(TClonesArray*)aodFromExt->GetList()->FindObject("Charm4Prong");
483         pdg=421;
484         if(fReadMC){
485           pdgdaughters =new Int_t[4];
486           pdgdaughters[0]=321;
487           pdgdaughters[1]=211;
488           pdgdaughters[2]=211;
489           pdgdaughters[3]=211;
490         }
491         break; 
492       case 5:
493         arrayProng=(TClonesArray*)aodFromExt->GetList()->FindObject("Charm3Prong");
494         pdg=4122;
495         if(fReadMC){
496           pdgdaughters =new Int_t[3];
497           pdgdaughters[0]=2212;//p
498           pdgdaughters[1]=321;//K
499           pdgdaughters[2]=211;//pi
500         }
501         break; 
502       }
503     }
504   } else if(aod) {
505     switch(fDecayChannel){
506     case 0:
507       arrayProng=(TClonesArray*)aod->GetList()->FindObject("Charm3Prong");
508       pdg=411;
509       if(fReadMC){
510         pdgdaughters =new Int_t[3];
511         pdgdaughters[0]=211;//pi
512         pdgdaughters[1]=321;//K
513         pdgdaughters[2]=211;//pi
514       }
515       break; 
516     case 1:
517       arrayProng=(TClonesArray*)aod->GetList()->FindObject("D0toKpi");
518       pdg=421;
519       if(fReadMC){
520         pdgdaughters =new Int_t[2];
521         pdgdaughters[0]=211;//pi 
522         pdgdaughters[1]=321;//K
523       }
524       break; 
525     case 2:
526       arrayProng=(TClonesArray*)aod->GetList()->FindObject("Dstar");
527       pdg=413;
528       if(fReadMC){
529         pdgdaughters[1]=211;//pi
530         pdgdaughters[0]=321;//K
531         pdgdaughters[2]=211;//pi (soft?)
532       }
533       break; 
534     case 3:
535       arrayProng=(TClonesArray*)aod->GetList()->FindObject("Charm3Prong");
536       pdg=431;
537       if(fReadMC){
538         pdgdaughters =new Int_t[3];
539         pdgdaughters[0]=321;//K
540         pdgdaughters[1]=321;//K
541         pdgdaughters[2]=211;//pi
542       }
543       break; 
544     case 4:
545       arrayProng=(TClonesArray*)aod->GetList()->FindObject("Charm4Prong");
546       pdg=421;
547       if(fReadMC){
548         pdgdaughters =new Int_t[4];
549         pdgdaughters[0]=321;
550         pdgdaughters[1]=211;
551         pdgdaughters[2]=211;
552         pdgdaughters[3]=211;
553       }
554       break; 
555     case 5:
556       arrayProng=(TClonesArray*)aod->GetList()->FindObject("Charm3Prong");
557       pdg=4122;
558       if(fReadMC){
559         pdgdaughters =new Int_t[3];
560         pdgdaughters[0]=2212;//p
561         pdgdaughters[1]=321;//K
562         pdgdaughters[2]=211;//pi
563       }
564       break; 
565     }
566   }
567   Bool_t isSimpleMode=kFALSE;
568   if(!arrayProng) {
569     AliInfo("Branch not found! The output will contain only trak related histograms\n");
570     isSimpleMode=kTRUE;
571     fNEntries->Fill(2);
572   }
573   
574   TClonesArray *mcArray = 0;
575   AliAODMCHeader *mcHeader = 0;
576
577   //check if MC
578   if(fReadMC) {
579     // load MC particles
580     mcArray = (TClonesArray*)aod->GetList()->FindObject(AliAODMCParticle::StdBranchName());
581     if(!mcArray) {
582       printf("AliAnalysisTaskSEHFQA::UserExec: MC particles branch not found!\n");
583       delete [] pdgdaughters;
584       return;
585     }
586     
587     // load MC header
588     mcHeader = (AliAODMCHeader*)aod->GetList()->FindObject(AliAODMCHeader::StdBranchName());
589     if(!mcHeader) {
590       printf("AliAnalysisTaskSEHFQA::UserExec: MC header branch not found!\n");
591       delete [] pdgdaughters;
592       return;
593     }
594   }
595   if(!aod) {delete [] pdgdaughters;return;}
596   // fix for temporary bug in ESDfilter 
597   // the AODs with null vertex pointer didn't pass the PhysSel
598   if(!aod->GetPrimaryVertex() || TMath::Abs(aod->GetMagneticField())<0.001) {    delete [] pdgdaughters;return;}
599
600   // count event
601   fNEntries->Fill(0); 
602
603   //count events with good vertex
604   // AOD primary vertex
605   AliAODVertex *vtx1 = (AliAODVertex*)aod->GetPrimaryVertex();
606   TString primTitle = vtx1->GetTitle();
607   if(primTitle.Contains("VertexerTracks") && vtx1->GetNContributors()>0) fNEntries->Fill(4);
608
609   // trigger class for PbPb C0SMH-B-NOPF-ALLNOTRD, C0SMH-B-NOPF-ALL
610   TString trigclass=aod->GetFiredTriggerClasses();
611   if(trigclass.Contains("C0SMH-B-NOPF-ALLNOTRD") || trigclass.Contains("C0SMH-B-NOPF-ALL")) fNEntries->Fill(5);
612
613
614   //select event
615   if(!fCuts->IsEventSelected(aod)) {
616     // rejected for pileup
617     if(fCuts->GetWhyRejection()==1) fNEntries->Fill(1);
618     delete [] pdgdaughters;
619     return;
620   }
621
622   if(fCuts->GetUseCentrality()){
623     Int_t runNumber = aod->GetRunNumber();
624     Int_t stdCent = (Int_t)(fCuts->GetCentrality(aod)+0.5);
625     Int_t secondCent = (Int_t)(fCuts->GetCentrality(aod,fEstimator)+0.5);
626     Int_t mincent=stdCent-stdCent%10;
627     ((AliCounterCollection*)fOutputCounters->FindObject("stdEstimator"))->Count(Form("centralityclass:%d_%d/Run:%d",mincent,mincent+10,runNumber));
628     mincent=secondCent-secondCent%10;
629     ((AliCounterCollection*)fOutputCounters->FindObject("secondEstimator"))->Count(Form("centralityclass:%d_%d/Run:%d",mincent,mincent+10,runNumber));
630
631     if(stdCent<fCuts->GetMinCentrality() || stdCent>fCuts->GetMaxCentrality()){
632       ((TH1F*)fOutputCheckCentrality->FindObject("hNtrackletsOut"))->Fill(aod->GetTracklets()->GetNumberOfTracklets());
633       ((TH1F*)fOutputCheckCentrality->FindObject("hMultOut"))->Fill(aod->GetHeader()->GetRefMultiplicity());
634     }else{
635       ((TH1F*)fOutputCheckCentrality->FindObject("hNtrackletsIn"))->Fill(aod->GetTracklets()->GetNumberOfTracklets());
636       ((TH1F*)fOutputCheckCentrality->FindObject("hMultIn"))->Fill(aod->GetHeader()->GetRefMultiplicity());
637     }
638   } else{
639     ((TH1F*)fOutputTrack->FindObject("hNtracklets"))->Fill(aod->GetTracklets()->GetNumberOfTracklets());
640     ((TH1F*)fOutputTrack->FindObject("hMult"))->Fill(aod->GetHeader()->GetRefMultiplicity());
641   }
642
643   Int_t ntracks=0;
644   Int_t isGoodTrack=0;
645
646   if(aod) ntracks=aod->GetNTracks();
647
648   //loop on tracks in the event
649   for (Int_t k=0;k<ntracks;k++){
650     AliAODTrack* track=aod->GetTrack(k);
651     AliAODPidHF* pidHF=fCuts->GetPidHF();
652     AliAODPid *pid = track->GetDetPid();
653     if(!pid)  {if (fDebug>1)cout<<"No AliAODPid found"<<endl; continue;}
654
655     Double_t times[AliPID::kSPECIES];
656     pid->GetIntegratedTimes(times);
657     
658     Double_t tofRes[AliPID::kSPECIES];
659     pid->GetTOFpidResolution(tofRes);
660
661     //check TOF
662     if(pidHF && pidHF->CheckStatus(track,"TOF")){
663       ((TH1F*)fOutputPID->FindObject("hTOFtime"))->Fill(times[AliPID::kProton]);
664       ((TH2F*)fOutputPID->FindObject("hTOFtimeKaonHyptime"))->Fill(track->P(),pid->GetTOFsignal()-times[3]); //3 is kaon
665       ((TH1F*)fOutputPID->FindObject("hTOFsig"))->Fill(pid->GetTOFsignal());
666       if (pid->GetTOFsignal()< 0) ((TH1F*)fOutputPID->FindObject("hTOFsig"))->Fill(-1);
667
668       // test a "simple" 160 ps TOF sigma PID
669       ((TH2F*)fOutputPID->FindObject("hTOFsigmaK160"))->Fill(track->P(),(pid->GetTOFsignal()-times[AliPID::kKaon])/160);
670       ((TH2F*)fOutputPID->FindObject("hTOFsigmaPion160"))->Fill(track->P(),(pid->GetTOFsignal()-times[AliPID::kPion])/160);
671       ((TH2F*)fOutputPID->FindObject("hTOFsigmaProton160"))->Fill(track->P(),(pid->GetTOFsignal()-times[AliPID::kProton])/160);
672       
673       // test TOF sigma PID
674       if (tofRes[2] != 0.) {   // protection against 'old' AODs...
675         ((TH2F*)fOutputPID->FindObject("hTOFsigmaKSigPid"))->Fill(track->P(),(pid->GetTOFsignal()-times[AliPID::kKaon])/tofRes[3]);
676         ((TH2F*)fOutputPID->FindObject("hTOFsigmaPionSigPid"))->Fill(track->P(),(pid->GetTOFsignal()-times[AliPID::kPion])/tofRes[2]);
677         ((TH2F*)fOutputPID->FindObject("hTOFsigmaProtonSigPid"))->Fill(track->P(),(pid->GetTOFsignal()-times[AliPID::kProton])/tofRes[4]);
678         for (Int_t iS=2; iS<5; iS++){ //we plot TOF Pid resolution for 3-sigma identified particles
679           if ( (TMath::Abs(times[iS]-pid->GetTOFsignal())/tofRes[iS])<3.){
680             switch (iS) {
681             case AliPID::kPion:
682               ((TH1F*)fOutputPID->FindObject("hTOFsigPid3sigPion"))->Fill(tofRes[iS]);
683               break;
684             case AliPID::kKaon:
685               ((TH1F*)fOutputPID->FindObject("hTOFsigPid3sigKaon"))->Fill(tofRes[iS]);
686               break;
687             case AliPID::kProton:
688               ((TH1F*)fOutputPID->FindObject("hTOFsigPid3sigProton"))->Fill(tofRes[iS]);
689               break;
690             default:
691               break;
692             }
693           }
694         }
695       }
696     }//if TOF status
697
698     if(pidHF && pidHF->CheckStatus(track,"TPC")){ 
699       Double_t alephParameters[5];
700       //this is recommended for LHC10d
701       alephParameters[0] = 1.34490e+00/50.;
702       alephParameters[1] =  2.69455e+01;
703       alephParameters[2] =  TMath::Exp(-2.97552e+01);
704       alephParameters[3] = 2.35339e+00;
705       alephParameters[4] = 5.98079e+00;
706       /*
707       //this is recommended for LHC10bc
708       alephParameters[0] = 0.0283086/0.97;
709       alephParameters[1] = 2.63394e+01;
710       alephParameters[2] = 5.04114e-11;
711       alephParameters[3] = 2.12543e+00;
712       alephParameters[4] = 4.88663e+00;
713       */
714       AliTPCPIDResponse* tpcres=new AliTPCPIDResponse();
715       tpcres->SetBetheBlochParameters(alephParameters[0],alephParameters[1],alephParameters[2],alephParameters[3],alephParameters[4]);
716       Double_t TPCp=pid->GetTPCmomentum();
717       Double_t TPCsignal=pid->GetTPCsignal();
718       ((TH1F*)fOutputPID->FindObject("hTPCsig"))->Fill(TPCsignal);
719       ((TH1F*)fOutputPID->FindObject("hTPCsigvsp"))->Fill(TPCp,TPCsignal);
720       //if (pidHF->IsKaonRaw(track, "TOF"))
721       ((TH2F*)fOutputPID->FindObject("hTPCsigmaK"))->Fill(TPCp,tpcres->GetNumberOfSigmas(TPCp,TPCsignal,track->GetTPCNcls(),AliPID::kKaon));
722       //if (pidHF->IsPionRaw(track, "TOF"))
723           ((TH2F*)fOutputPID->FindObject("hTPCsigmaPion"))->Fill(TPCp,tpcres->GetNumberOfSigmas(TPCp,TPCsignal,track->GetTPCNcls(),AliPID::kPion));
724       //if (pidHF->IsProtonRaw(track,"TOF"))
725       ((TH2F*)fOutputPID->FindObject("hTPCsigmaProton"))->Fill(TPCp,tpcres->GetNumberOfSigmas(TPCp,TPCsignal,track->GetTPCNcls(),AliPID::kProton));
726       delete tpcres;
727     }//if TPC status
728
729     //check clusters of the tracks
730     Int_t nclsTot=0,nclsSPD=0;
731     
732     for(Int_t l=0;l<6;l++) {
733       if(TESTBIT(track->GetITSClusterMap(),l)) {
734         nclsTot++; if(l<2) nclsSPD++;
735       }
736     }
737     ((TH1F*)fOutputTrack->FindObject("hnClsITS"))->Fill(nclsTot);
738     ((TH1F*)fOutputTrack->FindObject("hnClsSPD"))->Fill(nclsSPD);
739
740     if(!(track->GetStatus()&AliESDtrack::kTPCin) && track->GetStatus()&AliESDtrack::kITSrefit && !(track->GetStatus()&AliESDtrack::kITSpureSA)){//tracks retrieved in the ITS and not reconstructed in the TPC
741       ((TH1F*)fOutputTrack->FindObject("hnClsITS-SA"))->Fill(nclsTot);
742     }
743
744     if(isSimpleMode){
745
746       if (track->Pt()>0.3 &&
747           track->GetStatus()&AliESDtrack::kTPCrefit &&
748           track->GetStatus()&AliESDtrack::kITSrefit &&
749           /*nclsTot>3 &&*/
750           nclsSPD>0) {//fill hist good tracks
751
752         ((TH1F*)fOutputTrack->FindObject("hptGoodTr"))->Fill(track->Pt());
753         
754         isGoodTrack++;
755       
756         ((TH1F*)fOutputTrack->FindObject("hdistrGoodTr"))->Fill(isGoodTrack);
757       }
758     }//simple mode: no IsSelected on tracks: use "manual" cuts
759       
760   } //end loop on tracks
761
762   if(!isSimpleMode){
763     // loop over candidates
764     Int_t nCand = arrayProng->GetEntriesFast();
765     Int_t ndaugh=3;
766     if(fDecayChannel==AliAnalysisTaskSEHFQA::kD0toKpi) ndaugh=2;
767     if(fDecayChannel==AliAnalysisTaskSEHFQA::kD0toKpipipi) ndaugh=4;
768
769     for (Int_t iCand = 0; iCand < nCand; iCand++) {
770       AliAODRecoDecayHF *d = (AliAODRecoDecayHF*)arrayProng->UncheckedAt(iCand);
771
772       if(fReadMC){ 
773         Int_t labD = d->MatchToMC(pdg,mcArray,ndaugh,pdgdaughters);
774         if(labD>=0){
775           AliAODMCParticle *partD = (AliAODMCParticle*)mcArray->At(labD);
776           Int_t label=partD->GetMother();
777           AliAODMCParticle *mot = (AliAODMCParticle*)mcArray->At(label);
778           while(label>=0){//get first mother
779             mot = (AliAODMCParticle*)mcArray->At(label);
780             label=mot->GetMother();
781           }
782           Int_t pdgMotCode = mot->GetPdgCode();
783         
784           if(TMath::Abs(pdgMotCode)==4) fNEntries->Fill(6); //from primary charm
785           if(TMath::Abs(pdgMotCode)==5) fNEntries->Fill(7); //from beauty
786
787         }
788       }//end MC
789       else fNEntries->Fill(6); //count the candidates (data)
790
791       for(Int_t id=0;id<ndaugh;id++){
792
793         //other histograms to be filled when the cut object is given
794         AliAODTrack* track=(AliAODTrack*)d->GetDaughter(id);
795         if(fReadMC){
796           Int_t label=track->GetLabel();
797           if (label<0)fNEntries->Fill(8);
798         }
799         //track quality
800
801         if (fCuts->IsInFiducialAcceptance(d->Pt(),d->Y(pdg)) && fCuts->IsSelected(d,AliRDHFCuts::kTracks,aod)) {
802         
803           ((TH1F*)fOutputTrack->FindObject("hptGoodTr"))->Fill(track->Pt());
804           isGoodTrack++;
805       
806           ((TH1F*)fOutputTrack->FindObject("hdistrGoodTr"))->Fill(isGoodTrack);
807       
808           ((TH1F*)fOutputTrack->FindObject("hd0"))->Fill(d->Getd0Prong(id));
809   
810           if (fCuts->IsSelected(d,AliRDHFCuts::kAll,aod)){
811           
812             AliAODPid *pid = track->GetDetPid();
813             AliAODPidHF* pidHF=fCuts->GetPidHF();
814             Double_t times[5];
815             pid->GetIntegratedTimes(times);
816             if(pidHF && pidHF->CheckStatus(track,"TOF")) ((TH2F*)fOutputPID->FindObject("hTOFtimeKaonHyptimeAC"))->Fill(track->P(),pid->GetTOFsignal()-times[AliPID::kKaon]);
817             if(pidHF && pidHF->CheckStatus(track,"TPC")) ((TH2F*)fOutputPID->FindObject("hTPCsigvspAC"))->Fill(pid->GetTPCmomentum(),pid->GetTPCsignal());
818
819             fNEntries->Fill(3);
820           } //end analysis cuts
821         } //end acceptance and track cuts
822       } //end loop on tracks in the candidate
823     } //end loop on candidates  
824   }
825
826   delete [] pdgdaughters;
827   PostData(1,fNEntries);
828   PostData(2,fOutputPID);
829   PostData(3,fOutputTrack);
830   PostData(4,fCuts);
831   PostData(5,fOutputCounters);
832   PostData(6,fOutputCheckCentrality);
833 }
834
835 //____________________________________________________________________________
836 void AliAnalysisTaskSEHFQA::Terminate(Option_t */*option*/){
837   //terminate analysis
838
839   fNEntries = dynamic_cast<TH1F*>(GetOutputData(1));
840   if(!fNEntries){
841     printf("ERROR: %s not available\n",GetOutputSlot(1)->GetContainer()->GetName());
842     return;
843   }
844
845   fOutputPID = dynamic_cast<TList*> (GetOutputData(2));
846   if (!fOutputPID) {     
847     printf("ERROR: %s not available\n",GetOutputSlot(2)->GetContainer()->GetName());
848     return;
849   }
850
851   fOutputTrack = dynamic_cast<TList*> (GetOutputData(3));
852   if (!fOutputTrack) {     
853     printf("ERROR: %s not available\n",GetOutputSlot(3)->GetContainer()->GetName());
854     return;
855   }
856
857 }
858