]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/hfe/AliHFEefficiency.cxx
Corrections for Coverity warnings
[u/mrichter/AliRoot.git] / PWG3 / hfe / AliHFEefficiency.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 // Task for Efficiency studies
17 // Used for testing classes AliHFEcontainer and AliHFEfilter
18 // Creates Efficiency Histograms
19 //
20 // Author:
21 //   Markus Fasel <M.Fasel@gsi.de>
22 //
23 #include <TAxis.h>
24 #include <TCanvas.h>
25 #include <TFile.h>
26 #include <TGraphErrors.h>
27 #include <TH1F.h>
28 #include <TH2.h>
29 #include <TIterator.h>
30 #include <TLegend.h>
31 #include <TObjArray.h>
32 #include <TPad.h>
33
34 #include "AliAnalysisManager.h"
35 #include "AliCFAcceptanceCuts.h"
36 #include "AliCFTrackIsPrimaryCuts.h"
37 #include "AliCFContainer.h"
38 #include "AliCFEffGrid.h"
39 #include "AliESDEvent.h"
40 #include "AliHFEcollection.h"
41 #include "AliHFEcontainer.h"
42 #include "AliHFEcutStep.h"
43 #include "AliHFEefficiency.h"
44 #include "AliHFEextraCuts.h"
45 #include "AliHFEtrackFilter.h"
46 #include "AliHFEtools.h"
47 #include "AliMCEvent.h"
48 #include "AliMCEventHandler.h"
49
50 ClassImp(AliHFEefficiency)
51
52 AliHFEefficiency::AliHFEefficiency():
53   AliAnalysisTaskSE(),
54   fFilter(NULL),
55   fMCcut(NULL),
56   fAcceptanceCuts(NULL),
57   fEfficiency(NULL),
58   fOutput(NULL),
59   fCutTRD(kFALSE)
60 {
61   //
62   // Default constructor
63   //
64 }
65
66 AliHFEefficiency::AliHFEefficiency(const Char_t *name):
67   AliAnalysisTaskSE(name),
68   fFilter(NULL),
69   fMCcut(NULL),
70   fAcceptanceCuts(NULL),
71   fEfficiency(NULL),
72   fOutput(NULL),
73   fCutTRD(kFALSE)
74 {
75   //
76   // Main Constructor
77   //
78   DefineOutput(1, AliHFEcontainer::Class());
79   DefineOutput(2, TList::Class());
80
81   SetRunTerminate();
82 }
83
84 AliHFEefficiency::~AliHFEefficiency(){
85   //
86   // Destructor
87   //
88   if(fFilter) delete fFilter;
89   if(fEfficiency) delete fEfficiency;
90   if(fAcceptanceCuts) delete fAcceptanceCuts;
91   if(fOutput) delete fOutput;
92 }
93
94 void AliHFEefficiency::UserCreateOutputObjects(){
95   //
96   // Create the output objects
97   //
98   fEfficiency = new AliHFEcontainer("Efficiency");
99   fEfficiency->SetNumberOfVariables(4);
100   // InitializeVariables
101   fEfficiency->SetVariableName(0, "pt");
102   fEfficiency->MakeLogarithmicBinning(0, 40, 0.1, 10);
103   fEfficiency->SetVariableName(1, "Eta");
104   fEfficiency->MakeLinearBinning(1, 20, -0.9, 0.9);
105   fEfficiency->SetVariableName(2, "phi");
106   fEfficiency->MakeLinearBinning(2 , 18, 0, 2 * TMath::Pi());
107   fEfficiency->SetVariableName(3, "Charge");
108   fEfficiency->MakeLinearBinning(3, 2, -1.1, 1.1);
109   // New container
110   fEfficiency->CreateContainer("MCFilter", "Efficiency for MC Studies", 2);
111   AliCFContainer *cont = fEfficiency->GetCFContainer("MCFilter");
112   cont->SetStepTitle(0, "MC Signal");
113   cont->SetStepTitle(1, "MC Signal in acceptance");
114
115   fFilter = new AliHFEtrackFilter("testfilter");
116   fFilter->GenerateCutSteps();
117   fFilter->InitCF(fEfficiency);
118   fMCcut = fFilter->GetMCSignalCuts();
119   if(fCutTRD){
120     AliHFEcutStep *hfeTRD = fFilter->GetCutStep("HFETRD");
121     AliHFEextraCuts *hfecuts = dynamic_cast<AliHFEextraCuts *>(hfeTRD->GetCut("HFETRDCuts"));
122     if(hfecuts) hfecuts->SetMinTrackletsTRD(4);
123   }
124   AliHFEcutStep *hfeITS = fFilter->GetCutStep("HFEITS");
125   AliHFEextraCuts *hfeitscuts = dynamic_cast<AliHFEextraCuts *>(hfeITS->GetCut("HFEPixelsCuts"));
126   if(hfeitscuts)hfeitscuts->SetRequireITSpixel(AliHFEextraCuts::kFirst);
127
128   AliHFEcutStep *primary = fFilter->GetCutStep("Primary");
129   AliCFTrackIsPrimaryCuts *primcuts = dynamic_cast<AliCFTrackIsPrimaryCuts *>(primary->GetCut("PrimaryCuts"));
130   if(primcuts){ 
131     primcuts->SetMaxDCAToVertexXY(3);
132     primcuts->SetMaxDCAToVertexZ(5);
133   }
134
135   fAcceptanceCuts = new AliCFAcceptanceCuts("Acceptance", "MC Acceptance Cuts");
136   fAcceptanceCuts->SetMinNHitITS(3);
137   fAcceptanceCuts->SetMinNHitTPC(2);
138   fAcceptanceCuts->SetMinNHitTRD(0);
139
140   // create additional histos for testing purpose
141   fOutput = new AliHFEcollection("histos", "QA histograms");
142   fOutput->CreateTH1F("hNtracks","Number of Tracks per Event", 100, 0, 100);
143   fOutput->CreateTH1F("hPt","Pt of the tracks", 40, 0.1, 10, 0);
144   fOutput->CreateTH1F("kinkIndex", "Kink Index", 400, -200, 200);
145   fOutput->CreateTH1F("itspixel", "ITS PIXEL", 2, 0, 2);
146   fOutput->CreateTH1F("mcmother", "Mother PDG", 1000, -500, 500);
147   fOutput->CreateTH2F("ptres", "Pt Resolution", 40, 0.1, 10, 200, -1.5, 1.5, 0);
148 }
149
150 void AliHFEefficiency::UserExec(Option_t *){
151   //
152   // Event Loop
153   // Filter track, fill Efficiency container
154   //
155   AliESDEvent *esdevent = dynamic_cast<AliESDEvent *>(fInputEvent);
156   if(!esdevent){
157     AliError("ESD Event required");
158     return;
159   }
160   fEfficiency->NewEvent();
161   fFilter->SetRecEvent(fInputEvent);
162   if(fMCEvent){
163     AliMCEventHandler *mcH = dynamic_cast<AliMCEventHandler *>(AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
164     if(mcH &&(!mcH->InitOk())) return;
165     if(mcH &&(!mcH->TreeK())) return;
166     if(mcH &&(!mcH->TreeTR())) return;
167     fFilter->SetMC(fMCEvent);
168     FilterMC();
169   }
170   fFilter->FilterTracks(esdevent);
171   TObjArray *tracks = fFilter->GetFilteredTracks();
172   TIterator *iter = tracks->MakeIterator();
173   fOutput->Fill("hNtracks", tracks->GetEntriesFast());
174   AliESDtrack *track = NULL;
175   while((track = dynamic_cast<AliESDtrack *>(iter->Next()))){
176     fOutput->Fill("hPt", track->Pt());  
177     fOutput->Fill("kinkIndex", track->GetKinkIndex(0));
178     if(TESTBIT(track->GetITSClusterMap(), 0))
179       fOutput->Fill("itspixel",1);
180     else
181       fOutput->Fill("itspixel",0);
182     if(fMCEvent){
183       AliMCParticle *mctrack = dynamic_cast<AliMCParticle *>(fMCEvent->GetTrack(TMath::Abs(track->GetLabel())));
184       if(mctrack){
185         Int_t motherLabel = mctrack->Particle()->GetFirstMother();
186         if(motherLabel){
187           AliMCParticle *mother = dynamic_cast<AliMCParticle *>(fMCEvent->GetTrack(motherLabel));
188           if(mother)fOutput->Fill("mcmother", mother->Particle()->GetPdgCode());
189         }
190         fOutput->Fill("ptres", mctrack->Pt(), (track->Pt() - mctrack->Pt())/mctrack->Pt());
191       }
192     }
193   }
194   delete iter;
195   fFilter->Flush();
196   PostData(1, fEfficiency);
197   PostData(2, fOutput->GetList());
198 }
199
200 void AliHFEefficiency::Terminate(Option_t *){
201   //
202   // Evaluate Results
203   //
204   fEfficiency = dynamic_cast<AliHFEcontainer *>(GetOutputData(1));
205   if(!fEfficiency){
206     AliError("No Output data available");
207     return;
208   }
209
210   if(!IsRunningTerminate()) return;
211   PostProcess();
212
213   TList *l = dynamic_cast<TList *>(GetOutputData(2));
214   if(l) DrawPtResolution(l);
215 }
216
217 void AliHFEefficiency::FilterMC(){
218   //
219   // Cut MC tracks
220   //
221   Double_t cont[4] = {0., 0., 0., 0.};
222   AliMCParticle *track = NULL;
223   for(Int_t itrack = 0; itrack < fMCEvent->GetNumberOfTracks(); itrack++){
224     track = dynamic_cast<AliMCParticle *>(fMCEvent->GetTrack(itrack));
225     if(!track) continue;
226     if(!fMCcut->IsSelected(track)) continue;
227     cont[0] = track->Pt();
228     cont[1] = track->Eta();
229     cont[2] = track->Phi();
230     cont[3] = track->Charge()/3;
231     //AliMCParticle *mother = dynamic_cast<AliMCParticle *>(fMCEvent->GetTrack(track->Particle()->GetFirstMother()));
232     //if(TMath::Abs(mother->Particle()->GetPdgCode()) != 443) continue;
233     fEfficiency->FillCFContainer("MCFilter", 0, cont);
234     if(!fAcceptanceCuts->IsSelected(track)) continue;
235     fEfficiency->FillCFContainer("MCFilter", 1, cont);
236   }
237 }
238
239 void AliHFEefficiency::Load(const char* filename){
240   //
241   // Load results for post processing
242   //
243   TFile *input = TFile::Open(filename);
244   AliHFEcontainer *cin = dynamic_cast<AliHFEcontainer *>(input->Get("Efficiency"));
245   if(cin) fEfficiency = dynamic_cast<AliHFEcontainer *>(cin->Clone());
246   input->Close();
247   delete input;
248 }
249
250 void AliHFEefficiency::PostProcess(){
251   //
252   // Calculate Efficiencies
253   //
254   fEfficiency->Print();
255
256   // Prepare containers
257   AliCFContainer *crec  = fEfficiency->MakeMergedCFContainer("merged", "Merged Container", "MCFilter:testfilter_container_signal");
258   AliCFContainer *cmc   = fEfficiency->MakeMergedCFContainer("merged", "Merged Container", "MCFilter:testfilter_container_signalMC");
259   AliCFEffGrid *gridRec = new AliCFEffGrid("effrec", "Efficiency Calculation using rec vars", *crec);
260   AliCFEffGrid *gridMC  = new AliCFEffGrid("effmc", "Efficiency Calculation using mc vars", *cmc);
261
262   TCanvas *ctmp = NULL;
263   for(Int_t ivar = 0; ivar < crec->GetNVar(); ivar++){
264     ctmp = new TCanvas(Form("cEffSignal%s", crec->GetVarTitle(ivar)), Form("Signal Efficiency (%s)", crec->GetVarTitle(ivar)));
265     ctmp->cd();
266     DrawSignalEfficiency(gridRec, crec, ivar);
267     
268     ctmp = new TCanvas(Form("cCutEff%s", crec->GetVarTitle(ivar)), Form("Cut Step Efficiency (%s)", crec->GetVarTitle(ivar)));
269     ctmp->cd();
270     DrawCutEfficiency(gridRec, crec, ivar);
271   
272     ctmp = new TCanvas(Form("cEffSignalMC%s", crec->GetVarTitle(ivar)), Form("Signal Efficiency (%s, MC)", crec->GetVarTitle(ivar)));
273     ctmp->cd();
274     DrawSignalEfficiency(gridMC, cmc, ivar);
275
276     ctmp = new TCanvas(Form("cCutEffMC%s", crec->GetVarTitle(ivar)), Form("Cut Step Efficiency (%s, MC)", crec->GetVarTitle(ivar)));
277     ctmp->cd();
278     DrawCutEfficiency(gridMC, cmc, ivar);
279   } 
280
281   CalculatePTsmearing();
282
283   delete gridRec;
284   delete gridMC;
285   delete crec;
286   delete cmc;
287 }
288
289 void AliHFEefficiency::DrawSignalEfficiency(AliCFEffGrid *eff, AliCFContainer *cont, Int_t var){
290   //
291   // Efficiency of a cut step with respect to the signal
292   //
293   TH1 *hEff = NULL;
294   Bool_t kFirst = kTRUE;
295   TLegend *leg = new TLegend(0.5, 0.7, 0.89, 0.89);
296   for(Int_t istep = 1; istep < cont->GetNStep(); istep++){
297     eff->CalculateEfficiency(istep, 0);
298     hEff = eff->Project(var);
299     hEff->SetTitle(Form("Signal Efficiency (%s)", cont->GetVarTitle(var)));
300     hEff->GetYaxis()->SetRangeUser(0., 1.6);
301     hEff->GetYaxis()->SetTitle("Efficiency");
302     hEff->SetStats(kFALSE);
303     hEff->SetMarkerStyle(20+istep);
304     hEff->SetMarkerColor(kBlue - 5 + istep);
305     hEff->Draw(kFirst ? "ep" : "epsame");
306     kFirst = kFALSE;
307     leg->AddEntry(hEff, cont->GetStepTitle(istep), "p");
308   }
309   leg->SetBorderSize(0);
310   leg->SetFillStyle(0);
311   leg->Draw();
312   gPad->Update();
313 }
314
315 void AliHFEefficiency::DrawCutEfficiency(AliCFEffGrid *eff, AliCFContainer *cont, Int_t var){
316   //
317   // Efficiency of a cut step with respect to the previous one
318   //
319   Bool_t kFirst = kTRUE;
320   TH1 *hEff = NULL;
321   TLegend *leg = new TLegend(0.5, 0.7, 0.89, 0.89);
322   TString effTitle;
323   Int_t start = 0, length = 0;
324   for(Int_t istep = 1; istep < cont->GetNStep(); istep++){
325     eff->CalculateEfficiency(istep, istep -1);
326     hEff = eff->Project(var);
327     effTitle = hEff->GetTitle();
328     start = effTitle.Index("projected"); length = effTitle.Length() - start;
329     effTitle.Remove(start - 1, length + 1);
330     effTitle.ReplaceAll("Efficiency: ", "");
331     hEff->SetTitle(Form("Cut Efficiency (%s)", cont->GetVarTitle(var)));
332     hEff->GetYaxis()->SetRangeUser(0., 1.6);
333     hEff->GetYaxis()->SetTitle("Efficiency");
334     hEff->SetStats(kFALSE);
335     hEff->SetMarkerStyle(20+istep);
336     hEff->SetMarkerColor(kBlue - 5 + istep);
337     hEff->Draw(kFirst ? "ep" : "epsame");
338     kFirst = kFALSE;
339     leg->AddEntry(hEff, effTitle.Data(), "p");
340   }
341   leg->SetBorderSize(0);
342   leg->SetFillStyle(0);
343   leg->Draw();
344   gPad->Update();
345 }
346
347 void AliHFEefficiency::CalculatePTsmearing(){
348   //
349   // Efficiency of a cut step with respect to the same cut step filled with MC 
350   // Information
351   //
352   AliCFContainer *cont = fEfficiency->MakeMergedCFContainer("merged", "Merged Container", "testfilter_container_signal:testfilter_container_signalMC");
353   AliCFEffGrid *grid = new AliCFEffGrid("effc", "Efficiency Calculation", *cont);
354   Bool_t kFirst = kTRUE;
355   TH1 *hEff = NULL;
356   TCanvas *c2 = new TCanvas("cSmear", "pT modification");
357   c2->cd();
358   TLegend *leg1 = new TLegend(0.5, 0.7, 0.89, 0.89);
359   Int_t nStep = cont->GetNStep()/2;
360   for(Int_t istep = 0; istep < nStep; istep++){
361     grid->CalculateEfficiency(istep, istep + nStep);
362     hEff = grid->Project(0);
363     hEff->SetStats(kFALSE);
364     hEff->SetMarkerStyle(20+istep);
365     hEff->SetMarkerColor(kBlue - 5 + istep);
366     hEff->Draw(kFirst ? "ep" : "epsame");
367     kFirst = kFALSE;
368     leg1->AddEntry(hEff, cont->GetStepTitle(istep), "p");
369   }
370   leg1->Draw();
371   gPad->Update();
372
373   delete cont;
374   delete grid;
375 }
376
377 void AliHFEefficiency::DrawPtResolution(const TList * const l){
378   //
379   // Draw pt resolution
380   //
381   TH2 *h2 = dynamic_cast<TH2 *>(l->FindObject("ptres"));
382   if(!h2) return;
383   TGraphErrors *mean = new TGraphErrors(h2->GetNbinsX());
384   TGraphErrors *rms = new TGraphErrors(h2->GetNbinsX());
385
386   Double_t pt = 0., pterr = 0., mymean = 0., myrms = 0., mymeanerr = 0., myrmserr = 0.;
387   TH1 *htmp = NULL;
388   for(Int_t imom = 1; imom <= h2->GetXaxis()->GetLast(); imom++){
389     pt = h2->GetXaxis()->GetBinCenter(imom);
390     pterr = h2->GetXaxis()->GetBinWidth(imom)/2.;
391     htmp = h2->ProjectionY("py", imom, imom);
392     mymean = htmp->GetMean();
393     myrms = htmp->GetRMS();
394     mymeanerr = htmp->GetMeanError();
395     myrmserr = htmp->GetRMSError();
396     delete htmp;
397     mean->SetPoint(imom -1, pt, mymean);
398     rms->SetPoint(imom -1, pt, myrms);
399     mean->SetPointError(imom-1,pterr,mymeanerr);
400     rms->SetPointError(imom-1,pterr,myrmserr);
401   }
402
403   TCanvas *c1 = new TCanvas("cPtShift", "pT Shift");
404   c1->cd();
405   mean->SetMarkerStyle(22);
406   mean->SetMarkerColor(kBlue);
407   mean->Draw("ape");
408   rms->SetMarkerStyle(22);
409   rms->SetMarkerColor(kBlack);
410   rms->Draw("pesame");
411   TLegend *leg = new TLegend(0.5, 0.7, 0.89, 0.89);
412   leg->AddEntry(mean, "Mean", "lp");
413   leg->AddEntry(rms, "RMS", "lp");
414   leg->Draw();
415   gPad->Update();
416
417