3 #include <TCollection.h>
4 #include <TDirectory.h>
12 #include "AliESDtrack.h"
14 #include "AliTPCpidESD.h"
15 #include "AliTOFpidESD.h"
17 #include "AliESDpidCuts.h"
19 ClassImp(AliESDpidCuts)
21 const Int_t AliESDpidCuts::kNcuts = 3;
23 //_____________________________________________________________________
24 AliESDpidCuts::AliESDpidCuts(const Char_t *name, const Char_t *title):
25 AliAnalysisCuts(name, title)
28 , fTPCsigmaCutRequired(0)
29 , fTOFsigmaCutRequired(0)
30 , fCutTPCclusterRatio(0.)
31 , fMinMomentumTOF(0.5)
32 , fHcutStatistics(NULL)
33 , fHcutCorrelation(NULL)
36 // Default constructor
39 fTPCpid = new AliTPCpidESD;
40 fTOFpid = new AliTOFpidESD;
41 memset(fCutTPCnSigma, 0, sizeof(Float_t) * 2);
42 memset(fCutTOFnSigma, 0, sizeof(Float_t) * 2);
44 memset(fHclusterRatio, 0, sizeof(TH1F *) * 2);
45 memset(fHnSigmaTPC, 0, sizeof(TH1F *) * AliPID::kSPECIES * 2);
46 memset(fHnSigmaTOF, 0, sizeof(TH1F *) * AliPID::kSPECIES * 2);
49 //_____________________________________________________________________
50 AliESDpidCuts::AliESDpidCuts(const AliESDpidCuts &ref):
54 , fTPCsigmaCutRequired(ref.fTPCsigmaCutRequired)
55 , fTOFsigmaCutRequired(ref.fTOFsigmaCutRequired)
56 , fCutTPCclusterRatio(ref.fCutTPCclusterRatio)
57 , fMinMomentumTOF(ref.fMinMomentumTOF)
58 , fHcutStatistics(NULL)
59 , fHcutCorrelation(NULL)
64 fTPCpid = new AliTPCpidESD(*ref.fTPCpid);
65 fTOFpid = new AliTOFpidESD(*ref.fTOFpid);
66 memcpy(fCutTPCnSigma, ref.fCutTPCnSigma, sizeof(Float_t) * AliPID::kSPECIES * 2);
67 memcpy(fCutTOFnSigma, ref.fCutTOFnSigma, sizeof(Float_t) * AliPID::kSPECIES * 2);
69 if(ref.fHcutStatistics) fHcutStatistics = dynamic_cast<TH1I *>(ref.fHcutStatistics->Clone());
70 if(ref.fHcutCorrelation) fHcutCorrelation = dynamic_cast<TH2I *>(ref.fHcutCorrelation->Clone());
71 for(Int_t imode = 0; imode < 2; imode++){
72 if(ref.fHclusterRatio[imode]) fHclusterRatio[imode] = dynamic_cast<TH1F *>(ref.fHclusterRatio[imode]->Clone());
73 for(Int_t ispec = 0; ispec < AliPID::kSPECIES; ispec++){
74 if(fHnSigmaTPC[ispec][imode]) fHnSigmaTPC[ispec][imode] = dynamic_cast<TH1F *>(fHnSigmaTPC[ispec][imode]->Clone());
75 if(fHnSigmaTOF[ispec][imode]) fHnSigmaTOF[ispec][imode] = dynamic_cast<TH1F *>(fHnSigmaTPC[ispec][imode]->Clone());
80 //_____________________________________________________________________
81 AliESDpidCuts &AliESDpidCuts::operator=(const AliESDpidCuts &ref){
83 // Assignment operator
90 //_____________________________________________________________________
91 AliESDpidCuts::~AliESDpidCuts(){
98 delete fHcutStatistics;
99 delete fHcutCorrelation;
100 for(Int_t imode = 0; imode < 2; imode++){
101 delete fHclusterRatio[imode];
102 for(Int_t ispec = 0; ispec < AliPID::kSPECIES; ispec++){
103 delete fHnSigmaTPC[ispec][imode];
104 delete fHnSigmaTOF[ispec][imode];
109 //_____________________________________________________________________
110 Bool_t AliESDpidCuts::IsSelected(TObject *o){
114 if(TString(o->IsA()->GetName()).CompareTo("AliESDtrack")){
115 Char_t errormessage[256];
116 sprintf(errormessage, "Provided object not an AliESDtrack: Type %s", o->IsA()->GetName());
117 AliError(errormessage);
120 return AcceptTrack(const_cast<const AliESDtrack *>(dynamic_cast<AliESDtrack *>(o)));
123 //_____________________________________________________________________
124 void AliESDpidCuts::Copy(TObject &c) const {
128 AliESDpidCuts &target = dynamic_cast<AliESDpidCuts &>(c);
130 target.fTPCpid = new AliTPCpidESD;
131 target.fTOFpid = new AliTOFpidESD;
133 target.fCutTPCclusterRatio = fCutTPCclusterRatio;
134 target.fMinMomentumTOF = fMinMomentumTOF;
136 target.fTPCsigmaCutRequired = fTPCsigmaCutRequired;
137 target.fTOFsigmaCutRequired = fTOFsigmaCutRequired;
139 if(fHcutStatistics) target.fHcutStatistics = dynamic_cast<TH1I *>(fHcutStatistics->Clone());
140 if(fHcutCorrelation) target.fHcutCorrelation = dynamic_cast<TH2I *>(fHcutCorrelation->Clone());
141 for(Int_t imode = 0; imode < 2; imode++){
142 if(fHclusterRatio[imode]) target.fHclusterRatio[imode] = dynamic_cast<TH1F *>(fHclusterRatio[imode]->Clone());
143 for(Int_t ispec = 0; ispec < AliPID::kSPECIES; ispec++){
144 if(fHnSigmaTPC[ispec][imode]) target.fHnSigmaTPC[ispec][imode] = dynamic_cast<TH1F *>(fHnSigmaTPC[ispec][imode]->Clone());
145 if(fHnSigmaTOF[ispec][imode]) target.fHnSigmaTOF[ispec][imode] = dynamic_cast<TH1F *>(fHnSigmaTOF[ispec][imode]->Clone());
149 memcpy(target.fCutTPCnSigma, fCutTPCnSigma, sizeof(Float_t) * AliPID::kSPECIES * 2);
150 memcpy(target.fCutTOFnSigma, fCutTOFnSigma, sizeof(Float_t) * AliPID::kSPECIES);
152 AliESDpidCuts::Copy(c);
155 //_____________________________________________________________________
156 Long64_t AliESDpidCuts::Merge(TCollection *coll){
161 if(coll->IsEmpty()) return 1;
162 if(!HasHistograms()) return 0;
164 TIterator *iter = coll->MakeIterator();
166 AliESDpidCuts *ref = NULL;
168 while((o = iter->Next())){
169 ref = dynamic_cast<AliESDpidCuts *>(o);
171 if(!ref->HasHistograms()) continue;
173 fHcutStatistics->Add(ref->fHcutStatistics);
174 fHcutCorrelation->Add(ref->fHcutCorrelation);
175 for(Int_t imode = 0; imode < 2; imode++){
176 fHclusterRatio[imode]->Add(ref->fHclusterRatio[imode]);
177 for(Int_t ispec = 0; ispec < AliPID::kSPECIES; ispec++){
178 fHnSigmaTPC[ispec][imode]->Add(ref->fHnSigmaTPC[ispec][imode]);
179 fHnSigmaTOF[ispec][imode]->Add(ref->fHnSigmaTOF[ispec][imode]);
187 //_____________________________________________________________________
188 void AliESDpidCuts::DefineHistograms(Color_t color){
190 // Swich on QA and create the histograms
192 SetBit(kHasHistograms, kTRUE);
193 fHcutStatistics = new TH1I("fHcutStatistics", "Cut Statistics", kNcuts, 0, kNcuts);
194 fHcutStatistics->SetLineColor(color);
195 fHcutCorrelation = new TH2I("fHcutCorrelation", "Cut Correlation", kNcuts, 0, kNcuts, kNcuts, 0, kNcuts);
196 TString cutname[kNcuts] = {"TPCclusterRatio", "TPC sigma", "TOF sigma"};
197 for(Int_t icut = 0; icut < kNcuts; icut++){
198 fHcutStatistics->GetXaxis()->SetBinLabel(fHcutStatistics->GetXaxis()->GetFirst() + icut, cutname[icut].Data());
199 fHcutCorrelation->GetXaxis()->SetBinLabel(fHcutCorrelation->GetXaxis()->GetFirst() + icut, cutname[icut].Data());
200 fHcutCorrelation->GetYaxis()->SetBinLabel(fHcutCorrelation->GetYaxis()->GetFirst() + icut, cutname[icut].Data());
202 Char_t hname[256], htitle[256];
203 for(Int_t imode = 0; imode < 2; imode++){
204 sprintf(hname, "fHclusterRatio%s", imode ? "After" : "Before");
205 sprintf(htitle, "TPC cluster Ratio %s cuts;Ratio;Entries", imode ? "after" : "before");
206 fHclusterRatio[imode] = new TH1F(hname, htitle, 20, 0., 1.);
207 for(Int_t ispec = 0; ispec < AliPID::kSPECIES; ispec++){
208 sprintf(hname, "fHnSigma%sTPC%s", AliPID::ParticleName(ispec), imode ? "after" : "before");
209 sprintf(htitle, "TPC sigma for %s %s cuts;sigma;Entries", AliPID::ParticleName(ispec), imode ? "after" : "before");
210 fHnSigmaTPC[ispec][imode] = new TH1F(hname, htitle, 200, -10., 10.);
211 sprintf(hname, "fHnSigma%sTOF%s", AliPID::ParticleName(ispec), imode ? "after" : "before");
212 sprintf(htitle, "TOF sigma for %s %s cuts;sigma;Entries", AliPID::ParticleName(ispec), imode ? "after" : "before");
213 fHnSigmaTOF[ispec][imode] = new TH1F(hname, htitle, 200, -10., 10.);
218 //_____________________________________________________________________
219 Bool_t AliESDpidCuts::AcceptTrack(const AliESDtrack *track){
221 // Check whether the tracks survived the cuts
228 Long64_t cutRequired=0, cutFullfiled = 0;
229 Double_t clusterRatio = track->GetTPCNclsF() ? static_cast<Float_t>(track->GetTPCNcls())/static_cast<Float_t>(track->GetTPCNclsF()) : 1.;
230 if(fCutTPCclusterRatio > 0.){
231 SETBIT(cutRequired, kCutClusterRatioTPC);
232 if(clusterRatio >= fCutTPCclusterRatio)
233 SETBIT(cutFullfiled, kCutClusterRatioTPC);
235 // check TPC nSigma cut
236 Float_t nsigmaTPC[AliPID::kSPECIES], nsigma; // need all sigmas for QA plotting
237 for(Int_t ispec = 0; ispec < AliPID::kSPECIES; ispec++){
238 nsigmaTPC[ispec] = nsigma = fTPCpid->GetNumberOfSigmas(track, static_cast<AliPID::EParticleType>(ispec));
239 if(!(fTPCsigmaCutRequired & 1 << ispec)) continue;
240 SETBIT(cutRequired, kCutNsigmaTPC); // We found at least one species where the n-Sigma Cut is required
241 if(nsigma >= fCutTPCnSigma[2*ispec] && nsigma <= fCutTPCnSigma[2*ispec+1]) SETBIT(cutFullfiled, kCutNsigmaTPC); // Fullfiled for at least one species
243 // check TOF nSigma cut
244 Float_t nsigmaTOF[AliPID::kSPECIES]; // see above
245 Bool_t hasTOFpid = track->GetStatus() & AliESDtrack::kTOFpid; // only apply TOF n-sigma cut when PID Status Bit is set
246 for(Int_t ispec = 0; ispec < AliPID::kSPECIES; ispec++){
247 if(hasTOFpid) nsigmaTOF[ispec] = nsigma = fTOFpid->GetNumberOfSigmas(track, static_cast<AliPID::EParticleType>(ispec));
248 if(!(fTOFsigmaCutRequired && 1 << ispec)) continue;
249 SETBIT(cutRequired, kCutNsigmaTOF);
250 if(track->GetOuterParam()->P() >= fMinMomentumTOF){
251 if(hasTOFpid && nsigma <= fCutTOFnSigma[2*ispec] && nsigma >= fCutTOFnSigma[2*ispec+1]) SETBIT(cutFullfiled, kCutNsigmaTOF);
255 // Fill Histograms before cuts
257 fHclusterRatio[0]->Fill(clusterRatio);
258 for(Int_t ispec = 0; ispec < AliPID::kSPECIES; ispec++){
259 fHnSigmaTPC[ispec][0]->Fill(nsigmaTPC[ispec]);
260 if(hasTOFpid) fHnSigmaTOF[ispec][0]->Fill(nsigmaTOF[ispec]);
263 if(cutRequired != cutFullfiled){
264 // Fill cut statistics
266 for(Int_t icut = 0; icut < kNcuts; icut++){
267 if(TESTBIT(cutRequired, icut) && !TESTBIT(cutFullfiled, icut)){
269 fHcutStatistics->Fill(icut);
270 for(Int_t jcut = 0; jcut <= icut; jcut++)
271 if(TESTBIT(cutRequired, jcut) && !TESTBIT(cutFullfiled, jcut)) fHcutCorrelation->Fill(jcut, icut);
275 return kFALSE; // At least one cut is not fullfiled
278 // Fill Histograms after cuts
280 fHclusterRatio[1]->Fill(clusterRatio);
281 for(Int_t ispec = 0; ispec < AliPID::kSPECIES; ispec++){
282 fHnSigmaTPC[ispec][1]->Fill(nsigmaTPC[ispec]);
283 if(hasTOFpid) fHnSigmaTOF[ispec][1]->Fill(nsigmaTOF[ispec]);
290 //_____________________________________________________________________
291 void AliESDpidCuts::SaveHistograms(const Char_t * location){
293 // Save the histograms to a file
295 if(!HasHistograms()){
296 AliError("Histograms not on - Exiting");
299 if(!location) location = GetName();
300 gDirectory->mkdir(location);
301 gDirectory->cd(location);
302 fHcutStatistics->Write();
303 fHcutCorrelation->Write();
305 gDirectory->mkdir("before_cuts");
306 gDirectory->mkdir("after_cuts");
308 gDirectory->cd("before_cuts");
309 fHclusterRatio[0]->Write();
310 for(Int_t ispec = 0; ispec < AliPID::kSPECIES; ispec++){
311 fHnSigmaTPC[ispec][0]->Write();
312 fHnSigmaTOF[ispec][0]->Write();
315 gDirectory->cd("../after_cuts");
316 fHclusterRatio[1]->Write();
317 for(Int_t ispec = 0; ispec < AliPID::kSPECIES; ispec++){
318 fHnSigmaTPC[ispec][1]->Write();
319 fHnSigmaTOF[ispec][1]->Write();
322 gDirectory->cd("..");
325 //_____________________________________________________________________
326 void AliESDpidCuts::DrawHistograms(){
328 // Draw the Histograms
330 TCanvas *stat = new TCanvas("cutStat", "Cut Statistics", 640, 480);
332 fHcutStatistics->SetStats(kFALSE);
333 fHcutStatistics->Draw();
334 stat->SaveAs(Form("%s_%s.gif", GetName(), stat->GetName()));
336 TCanvas *correl = new TCanvas("cutCorrelation", "Cut Correlation", 640, 480);
338 fHcutCorrelation->SetStats(kFALSE);
339 fHcutCorrelation->Draw("colz");
340 correl->SaveAs(Form("%s_%s.gif", GetName(), correl->GetName()));
342 TCanvas *cRatio = new TCanvas("ClusterRatioTPC", "TPC cluster Ratio", 640, 480);
344 fHclusterRatio[0]->SetLineColor(kRed);
345 fHclusterRatio[0]->SetStats(kFALSE);
346 fHclusterRatio[0]->Draw();
347 fHclusterRatio[1]->SetLineColor(kBlue);
348 fHclusterRatio[1]->SetStats(kFALSE);
349 fHclusterRatio[1]->Draw("same");
350 cRatio->SaveAs(Form("%s_%s.gif", GetName(), cRatio->GetName()));
352 TCanvas *cNsigTPC = new TCanvas("NsigmaTPC", "TPC n-sigma", 640, 480);
353 cNsigTPC->Divide(3,2);
354 for(Int_t ispec = 0; ispec < AliPID::kSPECIES; ispec++){
355 cNsigTPC->cd(ispec + 1);
356 fHnSigmaTPC[ispec][0]->SetLineColor(kRed);
357 fHnSigmaTPC[ispec][0]->SetStats(kFALSE);
358 fHnSigmaTPC[ispec][0]->Draw();
359 fHnSigmaTPC[ispec][1]->SetLineColor(kBlue);
360 fHnSigmaTPC[ispec][1]->SetStats(kFALSE);
361 fHnSigmaTPC[ispec][1]->Draw("same");
363 cNsigTPC->SaveAs(Form("%s_%s.gif", GetName(), cNsigTPC->GetName()));
365 TCanvas *cNsigTOF = new TCanvas("NsigmaTOF", "TOF n-sigma", 640, 480);
366 cNsigTOF->Divide(3,2);
367 for(Int_t ispec = 0; ispec < AliPID::kSPECIES; ispec++){
368 cNsigTOF->cd(ispec + 1);
369 fHnSigmaTOF[ispec][0]->SetLineColor(kRed);
370 fHnSigmaTOF[ispec][0]->SetStats(kFALSE);
371 fHnSigmaTOF[ispec][0]->Draw();
372 fHnSigmaTOF[ispec][1]->SetLineColor(kBlue);
373 fHnSigmaTOF[ispec][1]->SetStats(kFALSE);
374 fHnSigmaTOF[ispec][1]->Draw("same");
376 cNsigTOF->SaveAs(Form("%s_%s.gif", GetName(), cNsigTOF->GetName()));