]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/dielectron/AliAnalysisTaskDielectronEfficiency.cxx
Updating the functionality of AliAnalysisHadEtCorrections to accomodate centrality...
[u/mrichter/AliRoot.git] / PWG3 / dielectron / AliAnalysisTaskDielectronEfficiency.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 /* $Id$ */
17
18 //#####################################################
19 //#                                                   #
20 //#  Simple efficiency study for dielectrons          #
21 //#  Author: Jens Wiechula Jens.Wiechula@cern.ch      #
22 //#                                                   #
23 //#####################################################
24
25 #include <TParticle.h>
26 #include <TParticlePDG.h>
27 #include <TDatabasePDG.h>
28 #include <TROOT.h>
29 #include "TChain.h"
30 #include <TCanvas.h>
31 // #include "TTree.h"
32 #include <TH1.h>
33 #include <TH2F.h>
34 #include <THashList.h>
35
36 #include "AliAnalysisManager.h"
37
38 #include "AliESDInputHandler.h"
39 #include "AliMCEventHandler.h"
40 #include "AliMCEvent.h"
41 #include "AliVEvent.h"
42 #include "AliESDEvent.h"
43 #include "AliESDtrack.h"
44 #include "AliStack.h"
45 #include "AliKFParticle.h"
46 #include "AliESDtrackCuts.h"
47 #include "AliKineTrackCuts.h"
48 #include "AliLog.h"
49
50 #include "AliDielectronHistos.h"
51 #include "AliAnalysisTaskDielectronEfficiency.h"
52
53 ClassImp(AliAnalysisTaskDielectronEfficiency)
54
55
56 //=================================================================================
57 AliAnalysisTaskDielectronEfficiency::AliAnalysisTaskDielectronEfficiency() :
58   AliAnalysisTask(),
59   fInputEvent(0),
60   fHist(0),
61   fESDtrackCuts(0),
62   fKineCutsLegs(0),
63   fKineCutsMother(0),
64   fIdMCMother(443),
65   fIdMCDaughterP(-11),
66   fIdMCDaughterN(11),
67   fPDG(TDatabasePDG::Instance())
68 {
69 }
70
71 //=================================================================================
72 AliAnalysisTaskDielectronEfficiency::AliAnalysisTaskDielectronEfficiency(const char *name) :
73   AliAnalysisTask(name,name),
74   fInputEvent(0),
75   fHist(0),
76   fESDtrackCuts(new AliESDtrackCuts),
77   fKineCutsLegs(new AliKineTrackCuts),
78   fKineCutsMother(new AliKineTrackCuts),
79   fIdMCMother(443),
80   fIdMCDaughterP(-11),
81   fIdMCDaughterN(11),
82   fPDG(TDatabasePDG::Instance())
83 {
84   //
85   // named constructor. This is the one that should be used by the user, oterwise the
86   // essential objects are not created!
87   //
88   DefineInput(0,TChain::Class());
89   DefineOutput(0, THashList::Class());
90 }
91
92 //=================================================================================
93 AliAnalysisTaskDielectronEfficiency::~AliAnalysisTaskDielectronEfficiency()
94 {
95   //
96   // dtor
97   //
98   if (fESDtrackCuts)   delete fESDtrackCuts;
99   if (fKineCutsLegs)   delete fKineCutsLegs;
100   if (fKineCutsMother) delete fKineCutsMother;
101 }
102 //=================================================================================
103 void AliAnalysisTaskDielectronEfficiency::CreateOutputObjects() {
104   //
105   // Create histograms
106   // Called once
107   //
108
109   //-------------------
110   // MC truth produced
111   fHist=new AliDielectronHistos;
112   fHist->AddClass("MC;MCcut;DataSameMother;Event;DataCuts;DataTRDCuts");
113
114   fHist->UserHistogram("MC",    "JpsiMCPt"  ,"MC Jpsi;Pt [GeV]"         ,100,0,10);
115   fHist->UserHistogram("MC",    "mass"    ,"MC Jpsi; Inv.Mass [GeV]" ,100,0,4);
116   fHist->UserHistogram("MC",    "e+Pt"    ,"MC e+ from JPsi;Pt [GeV]" ,100,0,10);
117   fHist->UserHistogram("MC",    "e-Pt"    ,"MC e- from JPsi;Pt [GeV]" ,100,0,10);
118   fHist->UserHistogram("MC", "dndyPt"      ,"MC Jpsi procution; Rapidity;Pt",100,-4,4,100,0,10);
119   fHist->UserHistogram("MC","dndy"        ,"MC dNdy Jpsi;Rapidity;Entries/event"    ,100,-4,4);
120   fHist->GetHistogram("MC","dndy")->Sumw2();
121
122   //-------------------
123   //MC truth after cuts
124   fHist->UserHistogram("MCcut",    "JpsiMCPt"  ,"MC Jpsi;Pt [GeV]"         ,100,0,10);
125   fHist->UserHistogram("MCcut",    "mass"    ,"MC Jpsi; Inv.Mass [GeV]" ,100,0,4);
126   fHist->UserHistogram("MCcut",    "e+Pt"    ,"MC e+ from JPsi;Pt [GeV]" ,100,0,10);
127   fHist->UserHistogram("MCcut",    "e-Pt"    ,"MC e- from JPsi;Pt [GeV]" ,100,0,10);
128   fHist->UserHistogram("MCcut", "dndyPt"      ,"MC Jpsi procution; Rapidity;Pt",100,-4,4,100,0,10);
129   fHist->UserHistogram("MCcut","dndy"        ,"MC dNdy Jpsi;Rapidity;Entries/event"    ,100,-4,4);
130   fHist->GetHistogram("MCcut","dndy")->Sumw2();
131   
132   //-----------------
133   //reconstructed data with cuts on MC truth
134   fHist->UserHistogram("DataSameMother","JpsiMCPt","Rec Jpsi; MC Pt [GeV]", 100,0,10);
135   fHist->UserHistogram("DataSameMother","dndyPtMC" ,"Rec Jpsi procution; Rapidity;Pt",100,-4,4,100,0,10);
136   fHist->UserHistogram("DataSameMother", "e+Pt"   ,"Rec e+ from JPsi;MC Pt [GeV]" ,100,0,10);
137   fHist->UserHistogram("DataSameMother", "e-Pt"   ,"Rec e- from JPsi;MC Pt [GeV]" ,100,0,10);
138   fHist->UserHistogram("DataSameMother","dndy"    ,"Rec Jpsi;Rapidity;Entries/event"             ,100,-4,4);
139   fHist->GetHistogram("DataSameMother","dndy")->Sumw2();
140   
141   fHist->UserHistogram("DataSameMother","mass"    ,"Rec Jpsi (KF); Inv.Mass [GeV]" ,100,0,4);
142   fHist->UserHistogram("DataSameMother","JpsiPt"  ,"Rec Jpsi (KF); Pt [GeV]"       ,100,0,10);
143   fHist->UserHistogram("DataSameMother","Chi2"    ,"Rec Jpsi (KF); #Chi^{2}"       ,100,0,50);
144   fHist->UserHistogram("DataSameMother","dndyPt"   ,"Rec Jpsi procution (KF); Rapidity;Pt",100,-4,4,100,0,10);
145   //------------------
146   // reconstructed data after ESD track cuts and cuts on MC truth
147   //------------------
148   fHist->UserHistogram("DataCuts","JpsiMCPt","Rec Jpsi; MC Pt [GeV]", 100,0,10);
149   fHist->UserHistogram("DataCuts","dndyPtMC" ,"Rec Jpsi procution; Rapidity;Pt",100,-4,4,100,0,10);
150   fHist->UserHistogram("DataCuts", "e+Pt"   ,"Rec e+ from JPsi;MC Pt [GeV]" ,100,0,10);
151   fHist->UserHistogram("DataCuts", "e-Pt"   ,"Rec e- from JPsi;MC Pt [GeV]" ,100,0,10);
152   fHist->UserHistogram("DataCuts","dndy"    ,"Rec Jpsi;Rapidity;Entries/event"             ,100,-4,4);
153   fHist->GetHistogram("DataCuts","dndy")->Sumw2();
154   
155   fHist->UserHistogram("DataCuts","mass"    ,"Rec Jpsi (KF); Inv.Mass [GeV]" ,100,0,4);
156   fHist->UserHistogram("DataCuts","JpsiPt"  ,"Rec Jpsi (KF); Pt [GeV]"       ,100,0,10);
157   fHist->UserHistogram("DataCuts","Chi2"    ,"Rec Jpsi (KF); #Chi^{2}"       ,100,0,50);
158   fHist->UserHistogram("DataCuts","dndyPt"   ,"Rec Jpsi procution (KF); Rapidity;Pt",100,-4,4,100,0,10);
159   //------------------
160   // after ESD track cuts + TRD cuts + MC cuts
161   //------------------
162   fHist->UserHistogram("DataTRDCuts","JpsiMCPt","Rec Jpsi; MC Pt [GeV]", 100,0,10);
163   fHist->UserHistogram("DataTRDCuts","dndyPtMC" ,"Rec Jpsi procution; Rapidity;Pt",100,-4,4,100,0,10);
164   fHist->UserHistogram("DataTRDCuts", "e+Pt"   ,"Rec e+ from JPsi;MC Pt [GeV]" ,100,0,10);
165   fHist->UserHistogram("DataTRDCuts", "e-Pt"   ,"Rec e- from JPsi;MC Pt [GeV]" ,100,0,10);
166   fHist->UserHistogram("DataTRDCuts","dndy"    ,"Rec Jpsi;Rapidity;Entries/event"             ,100,-4,4);
167   fHist->GetHistogram("DataTRDCuts","dndy")->Sumw2();
168   
169   fHist->UserHistogram("DataTRDCuts","mass"    ,"Rec Jpsi (KF); Inv.Mass [GeV]" ,100,0,4);
170   fHist->UserHistogram("DataTRDCuts","JpsiPt"  ,"Rec Jpsi (KF); Pt [GeV]"       ,100,0,10);
171   fHist->UserHistogram("DataTRDCuts","Chi2"    ,"Rec Jpsi (KF); #Chi^{2}"       ,100,0,50);
172   fHist->UserHistogram("DataTRDCuts","dndyPt"   ,"Rec Jpsi procution (KF); Rapidity;Pt",100,-4,4,100,0,10);
173
174   //-----------------
175   //Event information
176   //-----------------
177   fHist->UserHistogram("Event","NEvents","Number of events",1,0,1);
178
179 }
180
181 // //____________________________________________________________
182 void AliAnalysisTaskDielectronEfficiency::ConnectInputData(Option_t *) {
183   //
184   // connect the input data
185   //
186   fInputEvent=0;
187   TTree* tree=dynamic_cast<TTree*>(GetInputData(0));
188   if (!tree) {
189     printf("ERROR: Could not read chain from input slot 0\n");
190   } else {
191     AliInputEventHandler *eventH = dynamic_cast<AliInputEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
192     if (!eventH) {
193       AliError("Could not get ESDInputHandler");
194     } else {
195       fInputEvent = eventH->GetEvent();
196       AliInfo("*** CONNECTED NEW EVENT ****");
197     }
198   }
199 }
200
201 //=================================================================================
202 void AliAnalysisTaskDielectronEfficiency::Exec(Option_t *) {
203   //
204   // Main loop. Called for every event
205   // Process the event in FillPlots and post the data afterwards
206   //
207   if (!fInputEvent) {
208     Printf("ERROR: Could not get input event\n");
209     return;
210   }
211
212   FillPlots(fInputEvent);
213   
214   PostData(0, const_cast<THashList*>(fHist->GetHistogramList()));
215 }
216
217
218 //====================================================================================
219 void AliAnalysisTaskDielectronEfficiency::FillPlots(AliVEvent *event)
220 {
221   //
222   // Fill histograms
223   //
224   AliESDEvent *esd=dynamic_cast<AliESDEvent*>(event);
225   if (!esd) return;
226   Int_t ntrack=esd->GetNumberOfTracks();
227
228   // Fetch Stack 
229   AliMCEventHandler *mcH = (AliMCEventHandler*) ((AliAnalysisManager::GetAnalysisManager())->GetMCtruthEventHandler());
230   if(!mcH) {
231     AliError("No MC handler found\n")
232     return;
233   }
234   
235   AliMCEvent *mcev=mcH->MCEvent();
236   if (!mcev){
237     AliError("No MC event found\n")
238     return;
239   }
240   AliStack*  pStack = mcev->Stack();
241   
242   if (!pStack) return;
243
244   //fill event info
245   fHist->Fill("Event","NEvents",0);
246
247   //fill MC histograms
248   FillMCInfo(pStack);
249
250   //
251   Float_t massMother=0;
252   if (fIdMCMother>-1) fPDG->GetParticle(fIdMCMother)->Mass();
253   
254   TLorentzVector v;
255   //loop over all tracks
256   for (Int_t itrack=0; itrack<ntrack; ++itrack){
257     //negative particles only in this loop
258     AliESDtrack *trackN=esd->GetTrack(itrack);
259     if (trackN->Charge()!=-1) continue;
260
261     //MC truth
262     Int_t labelN=TMath::Abs(trackN->GetLabel());
263     if (labelN<0) continue;
264     TParticle *pN=pStack->Particle(labelN);
265     Int_t pdgN=pN->GetPdgCode();
266
267     //MC mother
268     Int_t idMotherN=pN->GetFirstMother();
269     TParticle *motherN=0;
270     Int_t pdgMotherN=-1;
271     if (fIdMCMother>-1&&idMotherN>-1){
272       motherN=pStack->Particle(idMotherN);
273       pdgMotherN=motherN->GetPdgCode();
274     }
275     
276     for (Int_t itrack2=0; itrack2<ntrack; ++itrack2){
277       //positive particles only in this loop
278       AliESDtrack *trackP=esd->GetTrack(itrack2);
279       if (trackP->Charge()!=1) continue;
280       
281       //MC truth
282       Int_t labelP=TMath::Abs(trackP->GetLabel());
283       if (labelP<0) continue;
284       TParticle *pP=pStack->Particle(labelP);
285       //       Int_t pdgP=pP->GetPdgCode();
286
287       //MC mother
288       Int_t idMotherP=pP->GetFirstMother();
289 //       TParticle *motherP=0;
290       //       Int_t pdgMotherP=0;
291       if (idMotherP>-1){
292 //         motherP=pStack->Particle(idMotherP);
293       //       pdgMotherP=motherP->GetPdgCode();
294       }
295       //===============
296       //Fill histograms
297       //===============
298       Bool_t motherOK=kFALSE;
299       if (fIdMCMother==-1) motherOK=kTRUE;
300       else if (pdgMotherN==fIdMCMother) motherOK=kTRUE;
301
302       if (pdgN==fIdMCDaughterN && motherOK){ //electron and mother is fIdMCMother
303         AliKFParticle electron(*trackN,11);
304         AliKFParticle positron(*trackP,-11);
305         AliKFParticle jpsi(electron);
306         jpsi+=positron;
307
308         Bool_t sameMother=kFALSE;
309         Bool_t motherCutOK=kFALSE;
310
311         if (fIdMCMother==-1) {
312           //accept as same mother if we don't requite
313           sameMother=kTRUE;
314           motherCutOK=kTRUE;
315         }
316         else
317         {
318           if (idMotherN==idMotherP) sameMother=kTRUE;
319           if (fKineCutsMother->IsSelected(motherN)) motherCutOK=kTRUE;
320         }
321
322         if ( sameMother  &&  // same mother
323              motherCutOK && //cuts mother MC truth
324             fKineCutsLegs->IsSelected(pN) && fKineCutsLegs->IsSelected(pP)){ //cuts legs MC truth
325
326           //MC only data
327           if (motherN){
328             fHist->Fill("DataSameMother","JpsiMCPt",motherN->Pt());
329             fHist->Fill("DataSameMother","e+Pt",pP->Pt());
330             fHist->Fill("DataSameMother","e-Pt",pN->Pt());
331             v.SetPxPyPzE(motherN->Px(),motherN->Py(),motherN->Pz(),motherN->Energy());
332             fHist->Fill("DataSameMother","dndy",v.Rapidity());
333             fHist->Fill("DataSameMother","dndyPtMC",v.Rapidity(),motherN->Pt());
334           }
335
336           //reconstructed data
337           v.SetPtEtaPhiM(jpsi.GetPt(),jpsi.GetEta(),jpsi.GetPhi(),massMother);
338 //           printf("Jpsi: %f,%f,%f,%f\n",jpsi.GetPt(),jpsi.GetEta(),jpsi.GetPhi(),massMother);
339           fHist->Fill("DataSameMother","JpsiPt",jpsi.GetPt());
340           fHist->Fill("DataSameMother","Chi2",jpsi.GetChi2()/jpsi.GetNDF());
341           fHist->Fill("DataSameMother","mass",jpsi.GetMass());
342           fHist->Fill("DataSameMother","dndyPt",v.Rapidity(),jpsi.GetPt());
343
344           //histograms after ESD cuts
345           if (fESDtrackCuts->IsSelected(trackN)&&fESDtrackCuts->IsSelected(trackP)){
346             //MC only
347             if (motherN) {
348               fHist->Fill("DataCuts","JpsiMCPt",motherN->Pt());
349               fHist->Fill("DataCuts","e+Pt",pP->Pt());
350               fHist->Fill("DataCuts","e-Pt",pN->Pt());
351               v.SetPxPyPzE(motherN->Px(),motherN->Py(),motherN->Pz(),motherN->Energy());
352               fHist->Fill("DataCuts","dndy",v.Rapidity());
353               fHist->Fill("DataCuts","dndyPtMC",motherN->Eta(),motherN->Pt());
354             }
355
356           //reconstructed data
357             v.SetPtEtaPhiM(jpsi.GetPt(),jpsi.GetEta(),jpsi.GetPhi(),massMother);
358             fHist->Fill("DataCuts","JpsiPt",jpsi.GetPt());
359             fHist->Fill("DataCuts","Chi2",jpsi.GetChi2()/jpsi.GetNDF());
360             fHist->Fill("DataCuts","mass",jpsi.GetMass());
361             fHist->Fill("DataCuts","dndyPt",v.Rapidity(),jpsi.GetPt());
362
363             //Additional TRD cuts
364             if ( ((trackN->GetStatus()&AliESDtrack::kTRDrefit)!=0) && trackN->GetTRDntrackletsPID()>4 ){
365               if ( ((trackP->GetStatus()&AliESDtrack::kTRDrefit)!=0) && trackP->GetTRDntrackletsPID()>4 ){
366                 if (motherN){
367                   fHist->Fill("DataTRDCuts","JpsiMCPt",motherN->Pt());
368                   fHist->Fill("DataTRDCuts","e+Pt",pP->Pt());
369                   fHist->Fill("DataTRDCuts","e-Pt",pN->Pt());
370                   v.SetPxPyPzE(motherN->Px(),motherN->Py(),motherN->Pz(),motherN->Energy());
371                   fHist->Fill("DataTRDCuts","dndy",v.Rapidity());
372                   fHist->Fill("DataTRDCuts","dndyPtMC",v.Rapidity(),motherN->Pt());
373                 }
374                 
375                 //reconstructed data
376                 v.SetPtEtaPhiM(jpsi.GetPt(),jpsi.GetEta(),jpsi.GetPhi(),massMother);
377                 fHist->Fill("DataTRDCuts","JpsiPt",jpsi.GetPt());
378                 fHist->Fill("DataTRDCuts","Chi2",jpsi.GetChi2()/jpsi.GetNDF());
379                 fHist->Fill("DataTRDCuts","mass",jpsi.GetMass());
380                 fHist->Fill("DataTRDCuts","dndyPt",v.Rapidity(),jpsi.GetPt());
381               }
382             }
383           }
384         }
385       }
386     }
387   }
388 }
389
390 //____________________________________________________________
391 // Int_t AliAnalysisTaskDielectronEfficiency::Merge(TList *list)
392 // {
393 //   //
394 //   // Merge function
395 //   //
396 //   if (!list) return 0;
397 //   if (list->IsEmpty()) return 1;
398 // 
399 //   TIter next(list);
400 //   while ( (TObject *o=next()) ){
401 //     AliAnalysisTaskDielectronEfficiency *task=dynamic_cast<AliAnalysisTaskDielectronEfficiency*>o;
402 //     if (!o) continue;
403 //     fNev+=task->fNev;
404 //   }
405 // }
406
407
408 void AliAnalysisTaskDielectronEfficiency::FillMCInfo(AliStack * const pStack)
409 {
410   //
411   // fill pure MC histograms
412   //
413   
414   TLorentzVector v;
415   //Fill MC info
416   for (Int_t ipart=0; ipart<pStack->GetNtrack(); ++ipart){
417     TParticle *part=pStack->Particle(ipart);
418 //     printf("Particle %d\n",part->GetPdgCode());
419     if (part->GetPdgCode()!=fIdMCMother || part->GetNDaughters()!=2) continue;
420     TParticle *d1=pStack->Particle(part->GetFirstDaughter());
421     TParticle *d2=pStack->Particle(part->GetLastDaughter());
422     TParticle *dP=0;
423     TParticle *dN=0;
424     if (fPDG->GetParticle(d1->GetPdgCode())->Charge()>0){
425       dP=d1;
426       dN=d2;
427     }else{
428       dP=d2;
429       dN=d1;
430     }
431     if ( dP->GetPdgCode()!=fIdMCDaughterP || dN->GetPdgCode()!=fIdMCDaughterN ) continue;
432     v.SetPxPyPzE(part->Px(),part->Py(),part->Pz(),part->Energy());
433     fHist->Fill("MC","JpsiMCPt",part->Pt());
434     fHist->Fill("MC","dndy",v.Rapidity());
435     fHist->Fill("MC","dndyPt",v.Rapidity(),part->Pt());
436     fHist->Fill("MC","e-Pt",dN->Pt());
437     fHist->Fill("MC","e+Pt",dP->Pt());
438     //e+ e- inv mass
439     TLorentzVector vE;
440     vE.SetPxPyPzE(dN->Px(),dN->Py(),dN->Pz(),dN->Energy());
441     TLorentzVector vP;
442     vP.SetPxPyPzE(dP->Px(),dP->Py(),dP->Pz(),dP->Energy());
443     fHist->Fill("MC","mass",(vE+vP).M());
444
445
446     //cuts
447     if (!fKineCutsMother->IsSelected(part)) continue;
448     if (!fKineCutsLegs->IsSelected(d1) || !fKineCutsLegs->IsSelected(d2) ) continue;
449       
450     v.SetPxPyPzE(part->Px(),part->Py(),part->Pz(),part->Energy());
451     fHist->Fill("MCcut","JpsiMCPt",part->Pt());
452     fHist->Fill("MCcut","dndy",v.Rapidity());
453     fHist->Fill("MCcut","dndyPt",v.Rapidity(),part->Pt());
454     fHist->Fill("MCcut","e-Pt",dN->Pt());
455     fHist->Fill("MCcut","e+Pt",dP->Pt());
456   }
457 }
458 void AliAnalysisTaskDielectronEfficiency::SetupDefaultCuts(Int_t type)
459 {
460   //
461   // setup standard ESD track cuts
462   //
463
464   if (type==0){
465     //ESD cuts
466     fESDtrackCuts->SetMaxDCAToVertexZ(3.0);
467     fESDtrackCuts->SetMaxDCAToVertexXY(3.0);
468     fESDtrackCuts->SetRequireTPCRefit(kTRUE);
469     fESDtrackCuts->SetRequireITSRefit(kTRUE);
470     fESDtrackCuts->SetAcceptKinkDaughters(kFALSE);
471     fESDtrackCuts->SetClusterRequirementITS(AliESDtrackCuts::kSPD,AliESDtrackCuts::kAny);
472     
473     fESDtrackCuts->SetMinNClustersTPC(50);
474     fESDtrackCuts->SetMaxChi2PerClusterTPC(4);
475
476     //MC cuts
477     fKineCutsLegs->SetEtaRange(-0.9,0.9);
478     fKineCutsMother->SetRapRange(-0.9,0.9);
479   } else if (type==1) {
480 //     fESDtrackCuts->SetMaxCovDiagonalElements(2, 2, .5, .5, 2);
481     fESDtrackCuts->SetMaxDCAToVertexZ(3.0);
482     fESDtrackCuts->SetMaxDCAToVertexXY(3.0);
483     fESDtrackCuts->SetRequireTPCRefit(kTRUE);
484     fESDtrackCuts->SetRequireITSRefit(kTRUE);
485     fESDtrackCuts->SetAcceptKinkDaughters(kFALSE);
486     fESDtrackCuts->SetClusterRequirementITS(AliESDtrackCuts::kSPD,AliESDtrackCuts::kAny);
487     
488     fESDtrackCuts->SetMinNClustersTPC(50);
489     fESDtrackCuts->SetMaxChi2PerClusterTPC(4);
490     
491     //MC cuts
492 //     fKineCutsLegs->SetEtaRange(-0.9,0.9);
493 //     fKineCutsMother->SetRapRange(-0.9,0.9);
494     
495   }
496 }
497
498 //===================================================================================
499 void AliAnalysisTaskDielectronEfficiency::Terminate(Option_t *) {
500   //
501   // Called once at the end of the query
502   //
503
504   
505   AliDielectronHistos *hist=new AliDielectronHistos;
506   hist->SetHistogramList(*(THashList*)GetOutputData(0));
507   
508   if (hist->GetHistogram("Event","NEvents")){
509     //get number of events
510     Double_t nev=hist->GetHistogram("Event","NEvents")->GetBinContent(1);
511     //
512     //normalise dndy histograms
513     //
514     hist->GetHistogram("DataSameMother","dndy")->Scale(1./nev);
515     hist->GetHistogram("MC","dndy")->Scale(1./nev);
516     hist->GetHistogram("DataCuts","dndy")->Scale(1./nev);
517     hist->GetHistogram("DataTRDCuts","dndy")->Scale(1./nev);
518
519     //
520     // create the efficiency histograms
521     //
522     // hEffTracking/2D only tracking effects, no esd cuts
523     // hEffESDCuts  tracking plus ESD cuts
524     // hEffTRDCuts  tracking plus ESD plus TRD cuts
525     //
526     TH1F *hEffTracking=(TH1F*)hist->GetHistogram("DataSameMother","JpsiMCPt")->Clone("Efficiency");
527     hEffTracking->Divide(hist->GetHistogram("MCcut","JpsiMCPt"));
528     hEffTracking->SetTitle("Efficiencies");
529     hist->UserHistogram("DataSameMother",hEffTracking);
530     
531     TH1F *hEffESDCuts=(TH1F*)hist->GetHistogram("DataCuts","JpsiMCPt")->Clone("Efficiency");
532     hEffESDCuts->Divide(hist->GetHistogram("MCcut","JpsiMCPt"));
533     hEffTracking->SetTitle("Efficiencies");
534     hist->UserHistogram("DataCuts",hEffESDCuts);
535     
536     TH1F *hEffTRDCuts=(TH1F*)hist->GetHistogram("DataTRDCuts","JpsiMCPt")->Clone("Efficiency");
537     hEffTRDCuts->Divide(hist->GetHistogram("MCcut","JpsiMCPt"));
538     hEffTRDCuts->SetTitle("Efficiencies");
539     hist->UserHistogram("DataTRDCuts",hEffTRDCuts);
540     
541     hist->DrawSame("Efficiency");
542
543     //2D efficiencies
544     TH2F *hEffTracking2D=(TH2F*)hist->GetHistogram("DataSameMother","dndyPtMC")->Clone("2DEfficiency");
545     hEffTracking2D->Divide(hist->GetHistogram("MCcut","dndyPt"));
546     hEffTracking2D->SetTitle("2D Efficiency - tracking");
547     hist->UserHistogram("DataSameMother",hEffTracking2D);
548
549     TH2F *hEffESDCuts2D=(TH2F*)hist->GetHistogram("DataCuts","dndyPtMC")->Clone("2DEfficiency");
550     hEffESDCuts2D->Divide(hist->GetHistogram("MCcut","dndyPt"));
551     hEffESDCuts2D->SetTitle("2D Efficiency - quality cuts");
552     hist->UserHistogram("DataCuts",hEffESDCuts2D);
553     
554     TH2F *hEffTRDCuts2D=(TH2F*)hist->GetHistogram("DataTRDCuts","dndyPtMC")->Clone("2DEfficiency");
555     hEffTRDCuts2D->Divide(hist->GetHistogram("MCcut","dndyPt"));
556     hEffTRDCuts2D->SetTitle("2D Efficiency - quality+TRD cuts");
557     hist->UserHistogram("DataTRDCuts",hEffTRDCuts2D);
558     
559 //
560     // Draw all histograms of all histogram classes
561     // Use the Draw functionality of AliDielectronHistos
562     //
563     hist->Draw();
564     
565     //
566     // Draw all histograms with the same name of all classes into one canvas
567     // Use the Draw functionality of AliDielectronHistos
568     //
569     hist->DrawSame("JpsiMCPt");
570   }
571
572   PostData(0, const_cast<THashList*>(hist->GetHistogramList()));
573 }
574