]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/hfe/AliAnalysisTaskHFEemcQA.cxx
920f137f50bef51ab1473e57f60e9fbda43c4351
[u/mrichter/AliRoot.git] / PWGHF / hfe / AliAnalysisTaskHFEemcQA.cxx
1 #include "TChain.h"
2 #include "TTree.h"
3 #include "TH1F.h"
4 #include "TCanvas.h"
5 #include "THnSparse.h"
6
7 #include "AliAnalysisTask.h"
8 #include "AliAnalysisManager.h"
9
10 #include "AliESDEvent.h"
11 #include "AliESDInputHandler.h"
12 #include "AliESDtrackCuts.h"
13 #include "AliAODEvent.h"
14 #include "AliAODHandler.h"
15
16 #include "AliPID.h"
17 #include "AliESDpid.h"
18 #include "AliAODPid.h"
19 #include "AliPIDResponse.h"
20 #include "AliHFEcontainer.h"
21 #include "AliHFEcuts.h"
22 #include "AliHFEpid.h"
23 #include "AliHFEpidBase.h"
24 #include "AliHFEpidQAmanager.h"
25 #include "AliHFEtools.h"
26 #include "AliCFContainer.h"
27 #include "AliCFManager.h"
28
29 #include "AliAnalysisTaskHFEemcQA.h"
30
31 //QA task for EMCAL electron analysis 
32
33 ClassImp(AliAnalysisTaskHFEemcQA)
34   //________________________________________________________________________
35   AliAnalysisTaskHFEemcQA::AliAnalysisTaskHFEemcQA(const char *name) 
36 : AliAnalysisTaskSE(name),
37   fVevent(0),
38   fESD(0),
39   fAOD(0),
40   fpidResponse(0),
41   fFlagSparse(kFALSE),
42   fOutputList(0),
43   fNevents(0),
44   fVtxZ(0),
45   fVtxX(0),
46   fVtxY(0),
47   fTrigMulti(0),
48   fHistClustE(0),
49   fEMCClsEtaPhi(0),
50   fNegTrkIDPt(0),
51   fTrkPt(0),
52   fTrketa(0),
53   fTrkphi(0),
54   fdEdx(0),
55   fTPCNpts(0),
56   fTPCnsig(0),
57   fHistPtMatch(0),
58   fEMCTrkMatch(0),
59   fEMCTrkPt(0),
60   fEMCTrketa(0),
61   fEMCTrkphi(0),
62   fEMCdEdx(0),
63   fEMCTPCnsig(0),
64   fEMCTPCNpts(0),
65   fHistdEdxEop(0),
66   fHistNsigEop(0),
67   fHistEop(0),
68   fEleCanTPCNpts(0),
69   fEleCanTPCNCls(0),
70   fEleCanITSNCls(0),
71   fEleCanITShit(0),
72   fEleCanSPD1(0),
73   fEleCanSPD2(0),
74   fEleCanSPDBoth(0),
75   fEleCanSPDOr(0),
76   fSparseElectron(0),
77   fvalueElectron(0)
78 {
79   // Constructor
80
81   fvalueElectron = new Double_t[6];
82   // Define input and output slots here
83   // Input slot #0 works with a TChain
84   DefineInput(0, TChain::Class());
85   // Output slot #0 id reserved by the base class for AOD
86   // Output slot #1 writes into a TH1 container
87   DefineOutput(1, TList::Class());
88 }
89 //________________________________________________________________________
90 AliAnalysisTaskHFEemcQA::AliAnalysisTaskHFEemcQA() 
91   : AliAnalysisTaskSE("DefaultTask_HfeEMCQA"),
92   fVevent(0),
93   fESD(0),
94   fAOD(0),
95   fpidResponse(0),
96   fFlagSparse(kFALSE),
97   fOutputList(0),
98   fNevents(0),
99   fVtxZ(0),
100   fVtxX(0),
101   fVtxY(0),
102   fTrigMulti(0),
103   fHistClustE(0),
104   fEMCClsEtaPhi(0),
105   fNegTrkIDPt(0),
106   fTrkPt(0),
107   fTrketa(0),
108   fTrkphi(0),
109   fdEdx(0),
110   fTPCNpts(0),
111   fTPCnsig(0),
112   fHistPtMatch(0),
113   fEMCTrkMatch(0),
114   fEMCTrkPt(0),
115   fEMCTrketa(0),
116   fEMCTrkphi(0),
117   fEMCdEdx(0),
118   fEMCTPCnsig(0),
119   fEMCTPCNpts(0),
120   fHistdEdxEop(0),
121   fHistNsigEop(0),
122   fHistEop(0),
123   fEleCanTPCNpts(0),
124   fEleCanTPCNCls(0),
125   fEleCanITSNCls(0),
126   fEleCanITShit(0),
127   fEleCanSPD1(0),
128   fEleCanSPD2(0),
129   fEleCanSPDBoth(0),
130   fEleCanSPDOr(0),
131   fSparseElectron(0),
132   fvalueElectron(0)
133 {
134   //Default constructor
135
136   fvalueElectron = new Double_t[6];
137   // Define input and output slots here
138   // Input slot #0 works with a TChain
139   DefineInput(0, TChain::Class());
140   // Output slot #0 id reserved by the base class for AOD
141   // Output slot #1 writes into a TH1 container
142   // DefineOutput(1, TH1I::Class());
143   DefineOutput(1, TList::Class());
144   //DefineOutput(3, TTree::Class());
145 }
146 //________________________________________________________________________
147 AliAnalysisTaskHFEemcQA::~AliAnalysisTaskHFEemcQA()
148 {
149   //Destructor 
150   delete fOutputList;
151   delete fSparseElectron;
152   delete []fvalueElectron;
153 }
154 //________________________________________________________________________
155 void AliAnalysisTaskHFEemcQA::UserCreateOutputObjects()
156 {
157   // Create histograms
158   // Called once
159   AliDebug(3, "Creating Output Objects");
160
161   /////////////////////////////////////////////////
162   //Automatic determination of the analysis mode//
163   ////////////////////////////////////////////////
164   AliVEventHandler *inputHandler = dynamic_cast<AliVEventHandler *>(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
165   if(!TString(inputHandler->IsA()->GetName()).CompareTo("AliAODInputHandler")){
166     SetAODAnalysis();
167   } else {
168     SetESDAnalysis();
169   }
170   printf("Analysis Mode: %s Analysis\n", IsAODanalysis() ? "AOD" : "ESD");
171
172   //////////////// 
173   //Output list//
174   ///////////////
175   fOutputList = new TList();
176   fOutputList->SetOwner();  
177
178   fNevents = new TH1F("fNevents","No of events",3,-0.5,2.5);
179   fOutputList->Add(fNevents);
180   fNevents->GetYaxis()->SetTitle("counts");
181   fNevents->GetXaxis()->SetBinLabel(1,"All");
182   fNevents->GetXaxis()->SetBinLabel(2,"With >2 Trks");
183   fNevents->GetXaxis()->SetBinLabel(3,"Vtx_{z}<10cm");
184
185   fVtxZ = new TH1F("fVtxZ","Z vertex position;Vtx_{z};counts",1000,-50,50);
186   fOutputList->Add(fVtxZ);
187
188   fVtxY = new TH1F("fVtxY","Y vertex position;Vtx_{y};counts",1000,-50,50);
189   fOutputList->Add(fVtxY);
190
191   fVtxX = new TH1F("fVtxX","X vertex position;Vtx_{x};counts",1000,-50,50);
192   fOutputList->Add(fVtxX);
193
194   fTrigMulti = new TH2F("fTrigMulti","Multiplicity distribution for different triggers; Trigger type; multiplicity",11,-1,10,2000,0,2000);
195   fOutputList->Add(fTrigMulti);
196
197   fHistClustE = new TH1F("fHistClustE", "EMCAL cluster energy distribution; Cluster E;counts", 500, 0.0, 50.0);
198   fOutputList->Add(fHistClustE);
199
200   fEMCClsEtaPhi = new TH2F("fEMCClsEtaPhi","EMCAL cluster #eta and #phi distribution;#eta;#phi",100,-0.9,0.9,200,0,6.3);
201   fOutputList->Add(fEMCClsEtaPhi);
202
203   fNegTrkIDPt = new TH1F("fNegTrkIDPt", "p_{T} distribution of tracks with negative track id;p_{T} (GeV/c);counts", 500, 0.0, 50.0); 
204   fOutputList->Add(fNegTrkIDPt);
205
206   fTrkPt = new TH1F("fTrkPt","p_{T} distribution of all tracks;p_{T} (GeV/c);counts",1000,0,100);
207   fOutputList->Add(fTrkPt);
208
209   fTrketa = new TH1F("fTrketa","All Track #eta distribution;#eta;counts",100,-1.5,1.5);
210   fOutputList->Add(fTrketa);
211
212   fTrkphi = new TH1F("fTrkphi","All Track #phi distribution;#phi;counts",100,0,6.3);
213   fOutputList->Add(fTrkphi);
214
215   fdEdx = new TH2F("fdEdx","All Track dE/dx distribution;p (GeV/c);dE/dx",200,0,20,500,0,160);
216   fOutputList->Add(fdEdx);
217
218   fTPCNpts = new TH2F("fTPCNpts","All track TPC Npoints used for dE/dx calculation;p (GeV/c);N points",200,0,20,200,0.,200.);
219   fOutputList->Add(fTPCNpts);
220
221   fTPCnsig = new TH2F("fTPCnsig","All Track TPC Nsigma distribution;p (GeV/c);#sigma_{TPC-dE/dx}",1000,0,50,200,-10,10);
222   fOutputList->Add(fTPCnsig);
223
224   fHistPtMatch = new TH1F("fHistPtMatch", "p_{T} distribution of tracks matched to EMCAL;p_{T} (GeV/c);counts",1000, 0.0, 100.0);
225   fOutputList->Add(fHistPtMatch);                                      
226
227   fEMCTrkMatch = new TH2F("fEMCTrkMatch","Distance of EMCAL cluster to its closest track;#phi;z",100,-0.3,0.3,100,-0.3,0.3);
228   fOutputList->Add(fEMCTrkMatch);
229
230   fEMCTrkPt = new TH1F("fEMCTrkPt","p_{T} distribution of tracks with EMCAL cluster;p_{T} (GeV/c);counts",1000,0,100);
231   fOutputList->Add(fEMCTrkPt);
232
233   fEMCTrketa = new TH1F("fEMCTrketa","#eta distribution of tracks matched to EMCAL;#eta;counts",100,-1.5,1.5);
234   fOutputList->Add(fEMCTrketa);
235
236   fEMCTrkphi = new TH1F("fEMCTrkphi","#phi distribution of tracks matched to EMCAL;#phi;counts",100,0,6.3);
237   fOutputList->Add(fEMCTrkphi);
238
239   fEMCdEdx = new TH2F("fEMCdEdx","dE/dx distribution of tracks matched to EMCAL;p (GeV/c);dE/dx",200,0,20,500,0,160);
240   fOutputList->Add(fEMCdEdx);
241     
242   fEMCTPCnsig = new TH2F("fEMCTPCnsig","TPC Nsigma distribution of tracks matched to EMCAL;p (GeV/c);#sigma_{TPC-dE/dx}",1000,0,50,200,-10,10);
243   fOutputList->Add(fEMCTPCnsig);
244
245   fEMCTPCNpts = new TH2F("fEMCTPCNpts","TPC Npoints used for dE/dx for tracks matched to EMCAL;p (GeV/c);N points",200,0,20,200,0.,200.);
246   fOutputList->Add(fEMCTPCNpts);
247
248   fHistEop = new TH2F("fHistEop", "E/p distribution;p_{T} (GeV/c);E/p", 200,0,20,60, 0.0, 3.0);
249   fOutputList->Add(fHistEop);
250
251   fHistdEdxEop = new TH2F("fHistdEdxEop", "E/p vs dE/dx;E/p;dE/dx", 60, 0.0, 3.0, 500,0,160);
252   fOutputList->Add(fHistdEdxEop);
253     
254   fHistNsigEop = new TH2F ("fHistNsigEop", "E/p vs TPC nsig",60, 0.0, 3.0, 200, -10,10);
255   fOutputList->Add(fHistNsigEop);
256
257   fEleCanTPCNpts = new TH2F("fEleCanTPCNpts","TPC Npoints used for dE/dx for electron candidates;p_{T} (GeV/c);N points",200,0,20,200,0,200);
258   fOutputList->Add(fEleCanTPCNpts);
259
260   fEleCanTPCNCls = new TH2F("fEleCanTPCNCls","TPC N clusters for electron candidates;p_{T} (GeV/c);N TPC clusters",200,0,20,171,-0.5,170.5);
261   fOutputList->Add(fEleCanTPCNCls);
262
263   fEleCanITSNCls = new TH2F("fEleCanITSNCls","ITS N clusters for electron candidates;p_{T} (GeV/c);N ITS clusters",200,0,20,8,-0.5,7.5);
264   fOutputList->Add(fEleCanITSNCls);
265
266   fEleCanITShit = new TH1F("fEleCanITShit","ITS hit map;ITS layer;counts",7,-0.5,6.5);
267   fOutputList->Add(fEleCanITShit);
268
269   fEleCanSPD1 = new TH2F("fEleCanSPD1","Hit on SPD layer 1;p_{T} (GeV/c);Hit",200,0,20,1,0,1);
270   fOutputList->Add(fEleCanSPD1);
271
272   fEleCanSPD2 = new TH2F("fEleCanSPD2","Hit on SPD layer 2;p_{T} (GeV/c);Hit",200,0,20,1,0,1);
273   fOutputList->Add(fEleCanSPD2);
274
275   fEleCanSPDBoth = new TH2F("fEleCanSPDBoth","Tracks with hits on both SPD layer;p_{T} (GeV/c);Hit",200,0,20,1,0,1);
276   fOutputList->Add(fEleCanSPDBoth);
277
278   fEleCanSPDOr = new TH2F("fEleCanSPDOr","Tracks with hits on both SPD layer;p_{T} (GeV/c);Hit",200,0,20,1,0,1);
279   fOutputList->Add(fEleCanSPDOr);
280
281   Int_t bins[6]={8,500,200,400,400,400}; //trigger, pt, TPCnsig, E/p, M20, M02
282   Double_t xmin[6]={-0.5,0,-10,0,0,0};
283   Double_t xmax[6]={7.5,25,10,2,2,2};
284   fSparseElectron = new THnSparseD ("Electron","Electron",6,bins,xmin,xmax);
285   fOutputList->Add(fSparseElectron);
286
287   PostData(1,fOutputList);
288 }
289
290 //________________________________________________________________________
291 void AliAnalysisTaskHFEemcQA::UserExec(Option_t *) 
292 {
293   // Main loop
294   // Called for each event
295   // Post output data.
296
297   UInt_t evSelMask=((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected();
298
299   fVevent = dynamic_cast<AliVEvent*>(InputEvent());
300   if (!fVevent) {
301     printf("ERROR: fVEvent not available\n");
302     return;
303   }
304
305   fESD = dynamic_cast<AliESDEvent*>(InputEvent());
306   if (fESD) {
307     //   printf("fESD available\n");
308     //return;
309   }
310
311   ////////////////////
312   //cuts initialised//
313   ///////////////////
314   AliESDtrackCuts* esdTrackCutsH = AliESDtrackCuts::GetStandardITSTPCTrackCuts2011(kFALSE);
315   esdTrackCutsH->SetMaxDCAToVertexXY(2.4);
316   esdTrackCutsH->SetMaxDCAToVertexZ(3.2);
317   esdTrackCutsH->SetDCAToVertex2D(kTRUE);
318
319   fAOD = dynamic_cast<AliAODEvent*>(InputEvent());
320   if (fAOD) {
321     // printf("fAOD available\n");
322     //return;
323   }
324
325   ///////////////////
326   //PID initialised//
327   //////////////////
328   fpidResponse = fInputHandler->GetPIDResponse();
329
330   ////////////////
331   //Event vertex//
332   ///////////////
333   Int_t ntracks;
334   ntracks = fVevent->GetNumberOfTracks();
335   printf("There are %d tracks in this event\n",ntracks);
336
337   fNevents->Fill(0); //all events
338   Double_t Zvertex = -100, Xvertex = -100, Yvertex = -100;
339   const AliVVertex *pVtx = fVevent->GetPrimaryVertex();
340   Double_t NcontV = pVtx->GetNContributors();
341   if(NcontV<2)return;
342   fNevents->Fill(1); //events with 2 tracks
343
344   Zvertex = pVtx->GetZ();  
345   Yvertex = pVtx->GetY();  
346   Xvertex = pVtx->GetX();  
347   fVtxZ->Fill(Zvertex);
348   fVtxX->Fill(Xvertex);
349   fVtxY->Fill(Yvertex);
350
351   /////////////////
352   //trigger check//
353   /////////////////
354   TString firedTrigger;
355   TString TriggerEG1("EG1");
356   TString TriggerEG2("EG2");
357   fVevent->GetFiredTriggerClasses();
358
359   Bool_t EG1tr = kFALSE;
360   Bool_t EG2tr = kFALSE;
361
362   if(firedTrigger.Contains(TriggerEG1))EG1tr = kTRUE;
363   if(firedTrigger.Contains(TriggerEG2))EG2tr = kTRUE;
364
365   Int_t trigger = -1;
366   if (fAOD){
367     Double_t multiplicity=fAOD->GetHeader()->GetRefMultiplicity();
368     fTrigMulti->Fill(-0.5, multiplicity);
369     if(evSelMask & AliVEvent::kAny) fTrigMulti->Fill(0.5, multiplicity);
370     if(evSelMask & AliVEvent::kMB) fTrigMulti->Fill(1.5, multiplicity);
371     if(evSelMask & AliVEvent::kINT7) fTrigMulti->Fill(2.5, multiplicity);
372     if(evSelMask & AliVEvent::kINT8) fTrigMulti->Fill(3.5, multiplicity);
373     if(evSelMask & AliVEvent::kEMC1) fTrigMulti->Fill(4.5, multiplicity);
374     if(evSelMask & AliVEvent::kEMC7) fTrigMulti->Fill(5.5, multiplicity);
375     if(evSelMask & AliVEvent::kEMC8) fTrigMulti->Fill(6.5, multiplicity);
376     if(evSelMask & AliVEvent::kEMCEJE) fTrigMulti->Fill(7.5, multiplicity);
377     if(evSelMask & AliVEvent::kEMCEGA) fTrigMulti->Fill(8.5, multiplicity);
378     if(evSelMask & AliVEvent::kEMCEGA & EG2tr) fTrigMulti->Fill(9.5, multiplicity);
379
380     if(evSelMask & AliVEvent::kMB) trigger =0;
381     if(evSelMask & AliVEvent::kINT7) trigger =1;
382     if(evSelMask & AliVEvent::kINT8) trigger =2;
383     if(evSelMask & AliVEvent::kEMC1) trigger =3;
384     if(evSelMask & AliVEvent::kEMC7) trigger =4;
385     if(evSelMask & AliVEvent::kEMC8) trigger =5;
386     if(evSelMask & AliVEvent::kEMCEJE) trigger =6;
387     if(evSelMask & AliVEvent::kEMCEGA) trigger =7;
388   }
389
390   ////////////////////
391   //event selection//
392   ///////////////////
393   if(fabs(Zvertex>10.0))return; 
394   fNevents->Fill(2); //events after z vtx cut
395
396   /////////////////////////////
397   //EMCAL cluster information//
398   ////////////////////////////
399   Int_t Nclust = 0;
400   Nclust = fVevent->GetNumberOfCaloClusters();
401   for(Int_t icl=0; icl<Nclust; icl++)
402   {
403     AliVCluster *clust = 0x0;
404     clust = fVevent->GetCaloCluster(icl);
405     if(clust && clust->IsEMCAL())
406     {
407       Double_t clustE = clust->E();
408       Float_t  emcx[3]; // cluster pos
409       clust->GetPosition(emcx);
410       TVector3 clustpos(emcx[0],emcx[1],emcx[2]);
411       Double_t emcphi = clustpos.Phi(); 
412       Double_t emceta = clustpos.Eta();
413       fHistClustE->Fill(clustE);
414       fEMCClsEtaPhi->Fill(emceta,emcphi);
415     }
416   }
417
418   /////////////////////////////////
419   //Look for kink mother for AOD//
420   /////////////////////////////////
421   Int_t numberofvertices = 100;
422   if(fAOD) numberofvertices = fAOD->GetNumberOfVertices();
423   Double_t listofmotherkink[numberofvertices];
424   Int_t numberofmotherkink = 0;
425   if(IsAODanalysis())
426   {
427     for(Int_t ivertex=0; ivertex < numberofvertices; ivertex++) {
428       AliAODVertex *aodvertex = fAOD->GetVertex(ivertex);
429       if(!aodvertex) continue;
430       if(aodvertex->GetType()==AliAODVertex::kKink) {
431         AliAODTrack *mother = (AliAODTrack *) aodvertex->GetParent();
432         if(!mother) continue;
433         Int_t idmother = mother->GetID();
434         listofmotherkink[numberofmotherkink] = idmother;
435         numberofmotherkink++;
436       }
437     }
438   } //+++
439
440
441   ///////////////
442   //Track loop///
443   ///////////////
444   for (Int_t iTracks = 0; iTracks < ntracks; iTracks++) {
445
446     AliVParticle* Vtrack = fVevent->GetTrack(iTracks);
447     if (!Vtrack) {
448       printf("ERROR: Could not receive track %d\n", iTracks);
449       continue;
450     }
451     AliVTrack *track = dynamic_cast<AliVTrack*>(Vtrack);
452     AliESDtrack *etrack = dynamic_cast<AliESDtrack*>(Vtrack);
453     AliAODTrack *atrack = dynamic_cast<AliAODTrack*>(Vtrack);
454
455     ////////////////////
456     //Apply track cuts//
457     ////////////////////
458     if(fAOD)
459       if(!atrack->TestFilterMask(AliAODTrack::kTrkGlobalNoDCA)) continue; //mimimum cuts
460
461     if(fESD)
462       if(!esdTrackCutsH->AcceptTrack(etrack))continue;
463
464     //reject kink
465     if(IsAODanalysis()){
466       Bool_t kinkmotherpass = kTRUE;
467       for(Int_t kinkmother = 0; kinkmother < numberofmotherkink; kinkmother++) {
468         if(track->GetID() == listofmotherkink[kinkmother]) {
469           kinkmotherpass = kFALSE;
470           continue;
471         }
472       }
473       if(!kinkmotherpass) continue;
474     }
475     else{
476       if(etrack->GetKinkIndex(0) != 0) continue;
477     }
478
479
480     ////////////////////
481     //Track properties//
482     ///////////////////
483     Double_t dEdx =-999, fTPCnSigma=-999;
484     dEdx = track->GetTPCsignal();
485     fTPCnSigma = fpidResponse->NumberOfSigmasTPC(track, AliPID::kElectron);
486
487     if(track->GetID()<0) fNegTrkIDPt->Fill(track->Pt());
488     fTrkPt->Fill(track->Pt());
489     fTrketa->Fill(track->Eta());
490     fTrkphi->Fill(track->Phi());
491     fdEdx->Fill(track->P(),dEdx);
492     fTPCNpts->Fill(track->P(),track->GetTPCsignalN());
493     fTPCnsig->Fill(track->P(),fTPCnSigma);
494
495     ///////////////////////////
496     //Track matching to EMCAL//
497     //////////////////////////
498     Int_t EMCalIndex = -1;
499     EMCalIndex = track->GetEMCALcluster();
500     if(EMCalIndex < 0) continue;
501     fHistPtMatch->Fill(track->Pt());
502
503     AliVCluster *clustMatch = (AliVCluster*)fVevent->GetCaloCluster(EMCalIndex);
504     if(clustMatch && clustMatch->IsEMCAL())
505     {
506       /////////////////////////////////////////////
507       //Properties of tracks matched to the EMCAL//
508       /////////////////////////////////////////////
509       fEMCTrkMatch->Fill(clustMatch->GetTrackDx(),clustMatch->GetTrackDz());
510       fEMCTrkPt->Fill(track->Pt());
511       fEMCTrketa->Fill(track->Eta());
512       fEMCTrkphi->Fill(track->Phi());
513       fEMCdEdx->Fill(track->P(),dEdx);
514       fEMCTPCnsig->Fill(track->P(),fTPCnSigma);
515       fEMCTPCNpts->Fill(track->P(),track->GetTPCsignalN());
516
517       //E/p distribution
518       Double_t clustMatchE = clustMatch->E();
519       Double_t eop = -1.0;
520       if(track->P()>0)eop = clustMatchE/track->P();
521
522       if(track->Pt()>1.0){
523         fHistdEdxEop->Fill(eop,dEdx);
524         fHistNsigEop->Fill(eop,fTPCnSigma);
525       }
526       fHistEop->Fill(track->Pt(),eop);
527
528       //EID THnsparse
529       fvalueElectron[0] = trigger;
530       fvalueElectron[1] = track->Pt();
531       fvalueElectron[2] = fTPCnSigma;
532       fvalueElectron[3] = eop;
533       fvalueElectron[4] = clustMatch->GetM20();
534       fvalueElectron[5] = clustMatch->GetM02();
535
536       if(fFlagSparse){
537         cout << "filling sparse"<<endl;
538         fSparseElectron->Fill(fvalueElectron);
539       }
540
541       ////////////////////////////////////////////////
542       //Track properties of EMCAL electron cadidates//
543       ///////////////////////////////////////////////
544       if(fTPCnSigma > -1 && fTPCnSigma < 3 && eop>0.8 && eop<1.2){
545         fEleCanTPCNpts->Fill(track->Pt(),track->GetTPCsignalN());
546         fEleCanTPCNCls->Fill(track->Pt(),track->GetTPCNcls());
547
548         Int_t fITSncls=0;
549         for(Int_t l=0;l<6;l++) {
550           if(TESTBIT(track->GetITSClusterMap(),l)) {
551             fEleCanITShit->Fill(l);
552             if(l==0) fEleCanSPD1->Fill(track->Pt(),0.5);
553             if(l==1) fEleCanSPD2->Fill(track->Pt(),0.5);
554             if(l==0 && l==1) fEleCanSPDBoth->Fill(track->Pt(),0.5);
555             if(l==0 || l==1) fEleCanSPDOr->Fill(track->Pt(),0.5);
556             fITSncls++;
557           }
558         }
559         fEleCanITSNCls->Fill(track->Pt(),fITSncls++);
560       }
561     }
562   } //track loop 
563
564   PostData(1, fOutputList);
565 }      
566 //________________________________________________________________________
567 void AliAnalysisTaskHFEemcQA::Terminate(Option_t *) 
568 {
569   // Draw result to the screen
570   // Called once at the end of the query
571
572   fOutputList = dynamic_cast<TList*> (GetOutputData(1));
573   if (!fOutputList) {
574     printf("ERROR: Output list not available\n");
575     return;
576   }
577
578 }