]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/GammaConv/AliConversionTrackCuts.cxx
-Move cent, z axis into sparse
[u/mrichter/AliRoot.git] / PWGGA / GammaConv / AliConversionTrackCuts.cxx
1 /**************************************************************************
2  * This file is property of and copyright by the ALICE HLT Project        *
3  * ALICE Experiment at CERN, All rights reserved.                         *
4  *                                                                        *
5  * Primary Author: Svein Lindal <slindal@fys.uio.no>                      *
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 /// @file   AliConversionTrackCuts.cxx
17 /// @author Svein Lindal
18 /// @brief  Base class for analysation of conversion particle - track correlations
19
20
21 #include "AliConversionTrackCuts.h"
22 //#include "AliAODTrack.h"
23 #include "AliAODEvent.h"
24 #include <TFormula.h>
25 #include <iostream>
26 #include "TH2F.h"
27 #include "AliESDtrackCuts.h"
28
29 using namespace std;
30 ClassImp(AliConversionTrackCuts)
31
32
33 const char* AliConversionTrackCuts::fgkCutNames[AliConversionTrackCuts::kNCuts] = {
34   "nClusTPC", 
35   "FoundFindable", 
36   "Chi2PerNDF", 
37   "Kink", 
38   "DCA_Z", 
39   "DCA_XY", 
40   "TPCRefit"
41   "kAccTracks"
42 };
43
44
45
46 //________________________________________________________________________
47 AliConversionTrackCuts::AliConversionTrackCuts() : 
48   AliAnalysisCuts(),
49   fFlagsOn(0x0),
50   fFlagsOff(0x0),
51   fRejectKinkDaughters(kTRUE),
52   fDCARfixed(kTRUE),
53   fDCARptFormula(""),
54   fDCARmax(1E20),
55   fDCAZfixed(kTRUE),
56   fDCAZptFormula(""),
57   fDCAZmax(1E20),
58   fDCAXYmax(1E20),
59   fSPDminNClusters(0),
60   fITSminNClusters(0),
61   fITSmaxChi2(1E20),
62   fTPCminNClusters(0),
63   fTPCClusOverFindable(0.0),
64   fTPCmaxChi2(1E20),
65   fAODTestFilterBit(-1),
66   fRequireTPCRefit(kFALSE),
67   fESDCuts(NULL),
68   fhPhi(NULL),
69   fhPt(NULL),
70   fhPhiPt(NULL),
71   fhdcaxyPt(NULL),
72   fhdcazPt(NULL),
73   fhdca(NULL),
74   fhnclpt(NULL),
75   fhnclsfpt(NULL),
76   fHistograms(NULL) {
77   //Constructor
78   //SetUpAxes();
79 }
80 //________________________________________________________________________
81 AliConversionTrackCuts::AliConversionTrackCuts(TString name, TString title = "title") : 
82   AliAnalysisCuts(name, title),
83   fFlagsOn(0x0),
84   fFlagsOff(0x0),
85   fRejectKinkDaughters(kTRUE),
86   fDCARfixed(kTRUE),
87   fDCARptFormula(""),
88   fDCARmax(1E20),
89   fDCAZfixed(kTRUE),
90   fDCAZptFormula(""),
91   fDCAZmax(1E20),
92   fDCAXYmax(1E20),
93   fSPDminNClusters(0),
94   fITSminNClusters(0),
95   fITSmaxChi2(1E20),
96   fTPCminNClusters(0),
97   fTPCClusOverFindable(0.0),
98   fTPCmaxChi2(1E20),
99   fAODTestFilterBit(-1),
100   fRequireTPCRefit(kFALSE),
101   fESDCuts(NULL),
102   fhPhi(NULL),  
103   fhPt(NULL),
104   fhPhiPt(NULL),
105   fhdcaxyPt(NULL),
106   fhdcazPt(NULL),
107   fhdca(NULL),
108   fhnclpt(NULL),
109   fhnclsfpt(NULL),
110   fHistograms(NULL)
111 {
112   //Constructor
113 //  SetUpAxes();
114 }
115
116
117 //________________________________________________________________________________
118  AliConversionTrackCuts::~AliConversionTrackCuts() {
119    ///destructor
120    // if(fHistograms)
121    //    delete fHistograms;
122    // fHistograms = NULL;
123
124    if(fESDCuts)
125      delete fESDCuts;
126    fESDCuts = NULL;
127 }
128
129 TList * AliConversionTrackCuts::CreateHistograms() {
130   //Create the histograms
131
132   if(!fHistograms) fHistograms = new TList();
133
134   fHistograms->SetOwner(kTRUE);
135   fHistograms->SetName("trackCuts");
136
137   fhPhi = new TH2F("phi", "phi", kNCuts+2, kPreCut -0.5, kNCuts + 0.5, 
138                                    32, 0, TMath::TwoPi());
139   TAxis * xax = fhPhi->GetXaxis();
140   for(Int_t i = 0; i < kNCuts; i++){
141         xax->SetBinLabel(xax->FindFixBin(i), fgkCutNames[i]);
142   }
143   fHistograms->Add(fhPhi);
144   
145
146
147   fhPt = new TH2F("pt", "pt", kNCuts+2, kPreCut -0.5, kNCuts + 0.5, 
148                                   20, 0., 20.);
149   xax = fhPt->GetXaxis();
150   for(Int_t i = 0; i < kNCuts; i++){
151         xax->SetBinLabel(xax->FindFixBin(i), fgkCutNames[i]);
152   }
153   fHistograms->Add(fhPt);
154
155   //  fhPhiPt = new TH2F("phipt", "phipt", 100, 0, 100, 64, 0, TMath::TwoPi());
156   //fHistograms->Add(fhPhiPt);
157
158   fhdcaxyPt = new TH2F("dcaxypt", "dcaxypt", 20, 0, 20, 50, 0, 5);
159   fHistograms->Add(fhdcaxyPt);
160
161   fhdcazPt = new TH2F("dcazpt", "dcazpt", 20, 0, 20, 50, 0, 5);
162   fHistograms->Add(fhdcazPt);
163
164   fhdca = new TH2F("dca", "dca", 60, -3, 3, 60, -3, 3);
165   fHistograms->Add(fhdca);
166
167   fhnclpt = new TH2F("nclstpcvspt", "nclstpcvspt", 20, 0, 20, 50, 0, 100);
168   fHistograms->Add(fhnclpt);
169
170   fhnclsfpt = new TH2F("nclsfpt", "nclsfpt", 20, 0, 20, 60, 0, 1.2);
171   fHistograms->Add(fhnclsfpt);
172   
173   return fHistograms;
174 }
175
176
177 void AliConversionTrackCuts::FillHistograms(Int_t cutIndex, AliVTrack * track, Bool_t passed = kFALSE) {
178   //Fill histograms
179   (void)passed;
180   fhPhi->Fill(cutIndex, track->Phi());
181   fhPt->Fill(cutIndex, track->Pt());
182   //if(passed) fhPhiPt->Fill(track->Pt(), track->Phi());
183
184 }
185
186 void AliConversionTrackCuts::FillDCAHist(Float_t dcaz, Float_t dcaxy, AliVTrack * track) {
187   fhdcaxyPt->Fill(track->Pt(), dcaxy);
188   fhdcazPt->Fill(track->Pt(), dcaz);
189   fhdca->Fill(dcaz, dcaxy);
190 }
191
192
193 Bool_t AliConversionTrackCuts::AcceptTrack(AliESDtrack * track) {
194   //Check esd track
195   if(!fESDCuts)
196     fESDCuts = AliESDtrackCuts::GetStandardTPCOnlyTrackCuts();
197   
198   FillHistograms(kPreCut, track);
199
200   if( !fESDCuts->IsSelected(track)) {
201     return kFALSE;
202     
203   }
204
205   fhnclpt->Fill(track->Pt(), track->GetTPCNcls());
206   if(track->GetTPCNclsF() > 0) fhnclsfpt->Fill(track->Pt(), ((Double_t) track->GetTPCNcls())/track->GetTPCNclsF());
207   FillHistograms(kPreCut + 1, track);
208
209   ///Get impact parameters
210   Double_t extCov[15];
211   track->GetExternalCovariance(extCov);
212   Float_t b[2];
213   Float_t bCov[3];
214   track->GetImpactParameters(b,bCov);
215   if (bCov[0]<=0 || bCov[2]<=0) {
216     AliDebug(1, "Estimated b resolution lower or equal zero!");
217     bCov[0]=0; bCov[2]=0;
218   }
219   
220   Float_t dcaToVertexXY = b[0];
221   Float_t dcaToVertexZ = b[1];
222   FillDCAHist(dcaToVertexZ, dcaToVertexXY, track);
223   return kTRUE;
224 }
225
226 Bool_t AliConversionTrackCuts::AcceptTrack(AliAODTrack * track) {
227   //Check aod track
228   FillHistograms(kPreCut, track);
229   if(!track->TestFilterBit(128)) {
230     return kFALSE;
231   }
232
233   if (track->GetTPCNcls() < fTPCminNClusters) return kFALSE;
234   FillHistograms(kCutNcls, track);
235
236   if (track->Chi2perNDF() > fTPCmaxChi2) return kFALSE;
237   FillHistograms(kCutNDF, track);
238
239   AliAODVertex *vertex = track->GetProdVertex();
240   if (vertex && fRejectKinkDaughters) {
241     if (vertex->GetType() == AliAODVertex::kKink) {
242       return kFALSE;
243     }
244   }
245   FillHistograms(kCutKinc, track);
246   
247   if(TMath::Abs(track->ZAtDCA()) > fDCAZmax) {
248     return kFALSE;
249   }
250   FillHistograms(kCutDCAZ, track);
251
252   Float_t xatdca = track->XAtDCA();
253   Float_t yatdca = track->YAtDCA();
254   Float_t xy = xatdca*xatdca + yatdca*yatdca;
255   if(xy > fDCAXYmax) {
256     return kFALSE;
257   }
258
259   FillHistograms(kCutDCAXY, track);
260
261
262   fhnclpt->Fill(track->Pt(), track->GetTPCNcls());
263   if(track->GetTPCNclsF() > 0) fhnclsfpt->Fill(track->Pt(), ((Double_t) track->GetTPCNcls())/track->GetTPCNclsF());
264   FillDCAHist(track->ZAtDCA(), TMath::Sqrt(track->XAtDCA()*track->XAtDCA() + track->YAtDCA()*track->YAtDCA()), track);
265
266
267   return kTRUE;
268 }
269
270 //_________________________________________________________________________________________________
271 void AliConversionTrackCuts::Print(const Option_t *) const
272 {
273 //
274 // Print information on this cut
275 //
276
277   printf("Cut name                : %s \n", GetName());
278   printf("Kink daughters are      : %s \n", (fRejectKinkDaughters ? "rejected" : "accepted"));
279   printf("TPC requirements        : clusters/findable %f, min. cluster = %d, max chi2 = %f, %s require refit\n", fTPCClusOverFindable, fTPCminNClusters, fTPCmaxChi2, (fRequireTPCRefit) ? "" : "Don't");
280   printf("ITS requirements        : min. cluster = %d (all), %d (SPD), max chi2 = %f \n", fITSminNClusters, fSPDminNClusters, fITSmaxChi2);
281   printf("DCA z cut               : fixed to %f cm \n", fDCAZmax);
282   printf("DCA xy cut              : fixed to %f cm \n", fDCAXYmax);
283 }
284  
285 //_________________________________________________________________________________________________
286
287 Bool_t AliConversionTrackCuts::IsSelected(TObject * object ) {
288   AliAODTrack * aodtrack = dynamic_cast<AliAODTrack*>(object);
289   if (aodtrack) {
290     return AcceptTrack(aodtrack);
291   } else {
292     AliESDtrack * track = dynamic_cast<AliESDtrack*>(object);
293     if (track)
294       return AcceptTrack(track);
295   }
296
297 return kFALSE;
298
299
300
301
302