2 /**************************************************************************
3 * Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
5 * Author: The ALICE Off-line Project. *
6 * Contributors are mentioned in the code where appropriate. *
8 * Permission to use, copy, modify and distribute this software and its *
9 * documentation strictly for non-commercial purposes is hereby granted *
10 * without fee, provided that the above copyright notice appears in all *
11 * copies and that both the copyright notice and this permission notice *
12 * appear in the supporting documentation. The authors make no claims *
13 * about the suitability of this software for any purpose. It is *
14 * provided "as is" without express or implied warranty. *
15 **************************************************************************/
17 //-----------------------------------------------------------------
18 // AliSpectraAODTrackCuts class
19 //-----------------------------------------------------------------
28 #include "AliAnalysisTask.h"
29 #include "AliAnalysisManager.h"
30 #include "AliAODTrack.h"
31 #include "AliPIDResponse.h"
32 #include "AliExternalTrackParam.h"
33 #include "AliAODMCParticle.h"
34 #include "AliAODEvent.h"
35 #include "AliAODInputHandler.h"
36 #include "AliAnalysisTaskESDfilter.h"
37 #include "AliAnalysisDataContainer.h"
38 #include "AliSpectraAODTrackCuts.h"
43 const char * AliSpectraAODTrackCuts::kBinLabel[] ={"TrkBit",
57 ClassImp(AliSpectraAODTrackCuts)
60 AliSpectraAODTrackCuts::AliSpectraAODTrackCuts(const char *name) : TNamed(name, "AOD Track Cuts"), fIsSelected(0), fTrackBits(0), fMinTPCcls(0), fRequestSPDcls(0), fEtaCutMin(0), fEtaCutMax(0), fDCACut(0), fPCut(0), fPtCut(0), fYCut(0),
61 fPtCutTOFMatching(0), fHistoCuts(0), fHistoNSelectedPos(0), fHistoNSelectedNeg(0), fHistoNMatchedPos(0), fHistoNMatchedNeg(0), fHistoEtaPhiHighPt(0), fTrack(0), fPIDResponse(0)
64 fHistoCuts = new TH1I("fTrkCuts", "Track Cuts", kNTrkCuts, -0.5, kNTrkCuts - 0.5);
65 for(Int_t ibin=1;ibin<=kNTrkCuts;ibin++)fHistoCuts->GetXaxis()->SetBinLabel(ibin,kBinLabel[ibin-1]);
67 const Double_t templBins[] = {0.20,0.30,0.4,0.5,0.6,0.7,0.8,0.9,1.0,1.2,1.4,1.6,1.8,2.0,2.4,2.8,3.2,3.6,4.0,5.0,6.0,7.0,8.0,9.0,10.,12.,15.};
68 const Int_t nbinsTempl=26;
70 fHistoNSelectedPos=new TH1F("fHistoNSelectedPos","fHistoNSelectedPos",nbinsTempl,templBins);
71 fHistoNSelectedPos->GetXaxis()->SetTitle("P_{T} (GeV / c)");
72 fHistoNSelectedNeg=new TH1F("fHistoNSelectedNeg","fHistoNSelectedNeg",nbinsTempl,templBins);
73 fHistoNSelectedNeg->GetXaxis()->SetTitle("P_{T} (GeV / c)");
74 fHistoNMatchedPos=new TH1F("fHistoNMatchedPos","fHistoNMatchedPos",nbinsTempl,templBins);
75 fHistoNMatchedPos->GetXaxis()->SetTitle("P_{T} (GeV / c)");
76 fHistoNMatchedNeg=new TH1F("fHistoNMatchedNeg","fHistoNMatchedNeg",nbinsTempl,templBins);
77 fHistoNMatchedNeg->GetXaxis()->SetTitle("P_{T} (GeV / c)");
78 fHistoEtaPhiHighPt=new TH2F("fHistoEtaPhiHighPt","fHistoEtaPhiHighPt",200,-1,1,400,0,7);
79 fHistoEtaPhiHighPt->SetXTitle("eta");
80 fHistoEtaPhiHighPt->SetYTitle("phi");
82 fEtaCutMin = -100000.0; // default value of eta cut ~ no cut
83 fEtaCutMax = 100000.0; // default value of eta cut ~ no cut
84 fDCACut = 100000.0; // default value of dca cut ~ no cut
85 fPCut = 100000.0; // default value of p cut ~ no cut
86 fPtCut = 100000.0; // default value of pt cut ~ no cut
87 fPtCutTOFMatching=0.6; //default value fot matching with TOF
88 fYCut = 100000.0; // default value of y cut ~ no cut
89 fMinTPCcls=70; // ncls in TPC
90 fRequestSPDcls=kFALSE; //request a hit in the SPD
93 //_______________________________________________________
94 Bool_t AliSpectraAODTrackCuts::IsSelected(AliAODTrack * track,Bool_t FillHistStat)
96 // Returns true if Track Cuts are selected and applied
99 printf("ERROR: Could not receive track");
104 if(!CheckTrackType()){
107 if(FillHistStat)fHistoCuts->Fill(kTrkBit);
109 if(!CheckTrackCuts()){
112 if(FillHistStat)fHistoCuts->Fill(kTrkCuts);
116 if(FillHistStat)fHistoCuts->Fill(kTrkEta);
120 if(FillHistStat)fHistoCuts->Fill(kTrkDCA);
124 if(FillHistStat)fHistoCuts->Fill(kTrkP);
128 if(FillHistStat)fHistoCuts->Fill(kTrkPt);
129 if(!CheckTOFMatching(FillHistStat)){
132 if(FillHistStat)fHistoCuts->Fill(kAccepted);
135 //_________________________________________________________
137 Bool_t AliSpectraAODTrackCuts::CheckTrackType()
140 if (fTrack->TestFilterBit(fTrackBits)) return kTRUE;
143 //_________________________________________________________
145 Bool_t AliSpectraAODTrackCuts::CheckTrackCuts()
147 // Check additional track Cuts
148 Bool_t PassTrackCuts=kTRUE;
149 if (!fTrack->HasPointOnITSLayer(0) && !fTrack->HasPointOnITSLayer(1) && fRequestSPDcls)PassTrackCuts=kFALSE; //FIXME 1 SPD for the moment
150 if (fTrack->GetTPCNcls()<fMinTPCcls)PassTrackCuts=kFALSE;
151 return PassTrackCuts;
153 //________________________________________________________
154 Bool_t AliSpectraAODTrackCuts::CheckEtaCut()
157 if (fTrack->Eta() < fEtaCutMax && fTrack->Eta() > fEtaCutMin) return kTRUE;
161 Bool_t AliSpectraAODTrackCuts::CheckYCut(Double_t mass)
163 // check if the rapidity is within the set range
165 if (mass > 0.) { y = fTrack->Y(mass); }//negative mass for unidentified particles
166 if (TMath::Abs(y) > fYCut || y < -998.) return kFALSE;
169 //_______________________________________________________
170 Bool_t AliSpectraAODTrackCuts::CheckDCACut()
173 if (TMath::Abs(fTrack->DCA()) < fDCACut) return kTRUE; //FIXME for newest AOD fTrack->DCA() always gives -999
176 //________________________________________________________
177 Bool_t AliSpectraAODTrackCuts::CheckPCut()
180 if (fTrack->P() < fPCut) return kTRUE;
183 //_______________________________________________________
184 Bool_t AliSpectraAODTrackCuts::CheckPtCut()
187 if (fTrack->Pt() < fPtCut) return kTRUE;
191 //_______________________________________________________
192 Bool_t AliSpectraAODTrackCuts::CheckTOFMatching(Bool_t FillHistStat)
194 if (fTrack->Pt() < fPtCutTOFMatching) return kTRUE;
196 if(FillHistStat)fHistoCuts->Fill(kTrkPtTOF);
197 if(fTrack->Charge()>0)fHistoNSelectedPos->Fill(fTrack->Pt());
198 else fHistoNSelectedNeg->Fill(fTrack->Pt());
200 // Get PID response object, if needed
202 AliAnalysisManager *man = AliAnalysisManager::GetAnalysisManager();
203 AliInputEventHandler* inputHandler = (AliInputEventHandler*)(man->GetInputEventHandler());
204 fPIDResponse = inputHandler->GetPIDResponse();
207 AliFatal("Cannot get pid response");
211 if(fPIDResponse->CheckPIDStatus(AliPIDResponse::kTOF,fTrack)==0)return kFALSE;
213 //check the bits of the selected particles
215 status=fTrack->GetStatus();
216 if((status&AliAODTrack::kTOFout)&&FillHistStat)fHistoCuts->Fill(kTrTOFout);
217 if((status&AliAODTrack::kTIME)&&FillHistStat)fHistoCuts->Fill(kTrTIME);
218 if((status&AliAODTrack::kTOFpid)&&FillHistStat)fHistoCuts->Fill(kTrTOFpid);
221 if(FillHistStat)fHistoCuts->Fill(kTOFMatching);
222 if(fTrack->Charge()>0)fHistoNMatchedPos->Fill(fTrack->Pt());
223 else fHistoNMatchedNeg->Fill(fTrack->Pt());
224 if(fTrack->Pt()>1.5){
225 //fHistoEtaPhiHighPt->Fill(fTrack->GetOuterParam()->Eta(),fTrack->GetOuterParam()->Phi());
226 //Printf("AliExternalTrackParam * extpar=(AliExternalTrackParam*)fTrack->GetOuterParam();");
227 //AliExternalTrackParam * extpar=(AliExternalTrackParam*)fTrack->GetOuterParam();
228 fHistoEtaPhiHighPt->Fill(fTrack->Eta(),fTrack->Phi());
229 //Printf("fHistoEtaPhiHighPt->Fill(extpar->Eta(),extpar->Phi());");
230 //fHistoEtaPhiHighPt->Fill(extpar->Eta(),extpar->Phi());
236 //_______________________________________________________
237 void AliSpectraAODTrackCuts::PrintCuts() const
240 cout << "Track Cuts" << endl;
241 cout << " > TrackBit\t" << fTrackBits << endl;
242 cout << " > Eta cut\t" << fEtaCutMin <<","<< fEtaCutMax << endl;
243 cout << " > DCA cut\t" << fDCACut << endl;
244 cout << " > P cut\t" << fPCut << endl;
245 cout << " > Pt cut \t" << fPtCut << endl;
246 cout << " > TPC cls \t" << fMinTPCcls << endl;
248 //_______________________________________________________
249 void AliSpectraAODTrackCuts::SetTrackType(UInt_t bit)
251 // Set the type of track to be used. The argument should be the bit number. The mask is produced automatically.
252 fTrackBits = (0x1 << (bit - 1));
254 //_______________________________________________________
256 Long64_t AliSpectraAODTrackCuts::Merge(TCollection* list)
258 // Merge a list of AliSpectraAODTrackCuts objects with this.
259 // Returns the number of merged objects (including this).
261 // AliInfo("Merging");
269 TIterator* iter = list->MakeIterator();
272 // collections of all histograms
273 TList collections;//FIXME we should only 1 collection
274 TList collections_histoNSelectedPos;
275 TList collections_histoNSelectedNeg;
276 TList collections_histoNMatchedPos;
277 TList collections_histoNMatchedNeg;
278 TList collections_histoEtaPhiHighPt;
282 while ((obj = iter->Next())) {
283 AliSpectraAODTrackCuts* entry = dynamic_cast<AliSpectraAODTrackCuts*> (obj);
287 TH1I * histo = entry->GetHistoCuts();
288 collections.Add(histo);
289 TH1F * histoNSelectedPos = entry->GetHistoNSelectedPos();
290 collections_histoNSelectedPos.Add(histoNSelectedPos);
291 TH1F * histoNSelectedNeg = entry->GetHistoNSelectedNeg();
292 collections_histoNSelectedNeg.Add(histoNSelectedNeg);
293 TH1F * histoNMatchedPos = entry->GetHistoNMatchedPos();
294 collections_histoNMatchedPos.Add(histoNMatchedPos);
295 TH1F * histoNMatchedNeg = entry->GetHistoNMatchedNeg();
296 collections_histoNMatchedNeg.Add(histoNMatchedNeg);
297 TH2F * histoEtaPhiHighPt = entry->GetHistoEtaPhiHighPt();
298 collections_histoEtaPhiHighPt.Add(histoEtaPhiHighPt);
302 fHistoCuts->Merge(&collections);
303 fHistoNSelectedPos->Merge(&collections_histoNSelectedPos);
304 fHistoNSelectedNeg->Merge(&collections_histoNSelectedNeg);
305 fHistoNMatchedPos->Merge(&collections_histoNMatchedPos);
306 fHistoNMatchedNeg->Merge(&collections_histoNMatchedNeg);
307 fHistoEtaPhiHighPt->Merge(&collections_histoEtaPhiHighPt);