]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/AliTPCcalibTracksCuts.cxx
Fixed sprintf, so that the pointer address is printed correctly
[u/mrichter/AliRoot.git] / TPC / AliTPCcalibTracksCuts.cxx
CommitLineData
10757ee9 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
17//////////////////////////////////////////////////////
18// //
19// Class to specify cuts for track analysis //
20// with AliTPCcalibTracks //
21// //
22//////////////////////////////////////////////////////
23
24#include <iostream>
25#include <TString.h>
26#include <TChain.h>
27#include <TList.h>
7eac86b0 28#include "AliTPCseed.h"
29#include "AliESDtrack.h"
10757ee9 30#include "AliTPCcalibTracksCuts.h"
31
32ClassImp(AliTPCcalibTracksCuts)
33
9389f9a4 34
35AliTPCcalibTracksCuts::AliTPCcalibTracksCuts():
36 TNamed("calibTracksCuts", "calibTracksCuts"),
37 fMinClusters(0), // number of clusters
38 fMinRatio(0), // kMinRratio = 0.4
39 fMax1pt(0), // kMax1pt = 0.5
40 fEdgeYXCutNoise(0), // kEdgeYXCutNoise = 0.13
7eac86b0 41 fEdgeThetaCutNoise(0) // kEdgeThetaCutNoise = 0.018
9389f9a4 42{
43 //
44 // default constructor
45 //
46}
47
48
10757ee9 49AliTPCcalibTracksCuts::AliTPCcalibTracksCuts(Int_t minClusters, Float_t minRatio, Float_t max1pt,
7eac86b0 50 Float_t edgeXZCutNoise, Float_t edgeThetaCutNoise):
9389f9a4 51 TNamed("calibTracksCuts", "calibTracksCuts"),
52 fMinClusters(minClusters), // number of clusters
53 fMinRatio(minRatio), // kMinRratio = 0.4
54 fMax1pt(max1pt), // kMax1pt = 0.5
55 fEdgeYXCutNoise(edgeXZCutNoise), // kEdgeYXCutNoise = 0.13
7eac86b0 56 fEdgeThetaCutNoise(edgeThetaCutNoise) // kEdgeThetaCutNoise = 0.018
9389f9a4 57{
10757ee9 58 //
59 // Constructor for AliTPCcalibTracksCuts
60 // specify the cuts to be set on the processed tracks
61 // default cuts are for comics
62 //
9389f9a4 63}
64
65AliTPCcalibTracksCuts::AliTPCcalibTracksCuts(AliTPCcalibTracksCuts *cuts):
66 TNamed(cuts->GetName(), cuts->GetTitle()),
67 fMinClusters(cuts->GetMinClusters()), // number of clusters
68 fMinRatio(cuts->GetMinRatio()), // kMinRratio = 0.4
69 fMax1pt( cuts->GetMax1pt()), // kMax1pt = 0.5
70 fEdgeYXCutNoise(cuts->GetEdgeYXCutNoise()), // kEdgeYXCutNoise = 0.13
7eac86b0 71 fEdgeThetaCutNoise( cuts->GetEdgeThetaCutNoise()) // kEdgeThetaCutNoise = 0.018
9389f9a4 72{
73 //
74 // copy constructor
75 //
10757ee9 76}
77
9389f9a4 78
79
10757ee9 80AliTPCcalibTracksCuts::~AliTPCcalibTracksCuts(){
81 //
82 // Destructor
83 //
84 cout << "AliTPCcalibTracksCuts destructor called, nothing happend." << endl;
85}
86
87
10757ee9 88
10757ee9 89
7eac86b0 90 AliTPCcalibTracksCuts * AliTPCcalibTracksCuts::CreateCuts(char* ctype){
10757ee9 91 //
7eac86b0 92 // Create predefined cuts
10757ee9 93 // (creates AliTPCcalibTracksCuts object)
10757ee9 94 //
95 // The following predefined sets of cuts can be selected:
96 // laser: 20, 0.4, 0.5, 0.13, 0.018
97 // cosmic: 20, 0.4, 0.5, 0.13, 0.018
98 // lowflux: 20, 0.4, 5, 0.2, 0.0001
99 // highflux: 20, 0.4, 5, 0.2, 0.0001
100 //
101
102 TString cutType(ctype);
103 cutType.ToUpper();
104 AliTPCcalibTracksCuts *cuts = 0;
105 if (cutType == "LASER")
7eac86b0 106 cuts = new AliTPCcalibTracksCuts(20, 0.4, 5, 0.13, 0.018);
10757ee9 107 else if (cutType == "COSMIC")
7eac86b0 108 cuts = new AliTPCcalibTracksCuts(20, 0.4, 0.5, 0.13, 0.018);
10757ee9 109 else if (cutType == "LOWFLUX")
7eac86b0 110 cuts = new AliTPCcalibTracksCuts(20, 0.4, 5, 0.2, 0.0001);
10757ee9 111 else if (cutType == "HIGHFLUX")
7eac86b0 112 cuts = new AliTPCcalibTracksCuts(20, 0.4, 5, 0.2, 0.0001);
10757ee9 113 else {
7eac86b0 114 cuts = new AliTPCcalibTracksCuts(20, 0.4, 5, 0.2, 0.0001);
115 cerr << "WARNING! unknown type '" << ctype << "', cuts set to default values for cosmics." << endl;
116 cutType = "COSMIC";
10757ee9 117 }
10757ee9 118 cout << "Cuts were set to predefined set: " << cutType << endl;
7eac86b0 119 return cuts;
10757ee9 120}
121
7eac86b0 122
123
124Int_t AliTPCcalibTracksCuts::AcceptTrack(const AliTPCseed * track) const {
125 //
126 // Function, that decides wheather a given track is accepted for
127 // the analysis or not.
128 // Returns 0 if a track is accepted or an integer different from 0
129 // to indicate the failed cut
130 //
131
132 //
133 // edge induced noise tracks - NEXT RELEASE will be removed during tracking
134 if ( TMath::Abs(track->GetY() / track->GetX()) > fEdgeYXCutNoise )
135 if ( TMath::Abs(track->GetTgl()) < fEdgeThetaCutNoise ) return 1;
136 if (track->GetNumberOfClusters() < fMinClusters) return 2;
137 Float_t ratio = track->GetNumberOfClusters() / (track->GetNFoundable() + 1.);
138 if (ratio < fMinRatio) return 3;
139 // Float_t mpt = track->Get1Pt(); // Get1Pt() doesn't exist any more
140 Float_t mpt = track->GetSigned1Pt();
141 if (TMath::Abs(mpt) > fMax1pt) return 4;
142
143 return 0;
144}
145
146Int_t AliTPCcalibTracksCuts::AcceptTrack(const AliESDtrack * track) const {
147 //
148 // Function, that decides wheather a given track is accepted for
149 // the analysis or not.
150 // Returns 0 if a track is accepted or an integer different from 0
151 // to indicate the failed cut
152 //
153
154 //
155 // edge induced noise tracks - NEXT RELEASE will be removed during tracking
156 if ( TMath::Abs(track->GetY() / track->GetX()) > fEdgeYXCutNoise )
157 if ( TMath::Abs(track->GetTgl()) < fEdgeThetaCutNoise ) return 1;
158 if (track->GetTPCNcls() < fMinClusters) return 2;
159 Float_t ratio = track->GetTPCNcls() / (track->GetTPCNclsF() + 1.);
160 if (ratio < fMinRatio) return 3;
161 // Float_t mpt = track->Get1Pt(); // Get1Pt() doesn't exist any more
162 Float_t mpt = track->GetSigned1Pt();
163 if (TMath::Abs(mpt) > fMax1pt) return 4;
164
165 return 0;
10757ee9 166}
167
7eac86b0 168void AliTPCcalibTracksCuts::Print(Option_t*) const {
10757ee9 169 //
170 // Print the cut contents
171 //
7eac86b0 172 cout << "<AliTPCcalibTracksCuts>: The following cuts are specified: " << endl;
173 cout << "fMinClusters: " << fMinClusters << endl;
174 cout << "fMinRatio: " << fMinRatio << endl;
175 cout << "fMax1pt: " << fMax1pt << endl;
176 cout << "fEdgeYXCutNoise: " << fEdgeYXCutNoise << endl;
177 cout << "fEdgeThetaCutNoise: " << fEdgeThetaCutNoise << endl;
10757ee9 178} // Prints out the specified cuts