Line data Source code
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 : // A GTU track //
19 : // //
20 : // Author: J. Klein (Jochen.Klein@cern.ch) //
21 : // //
22 : ////////////////////////////////////////////////////////////////////////////
23 :
24 : /* $Id: AliTRDtrackGTU.cxx 27566 2008-07-24 15:31:08Z cblume $ */
25 :
26 : #include "TObject.h"
27 : #include "TObjArray.h"
28 : #include "TClass.h"
29 : #include "TMath.h"
30 : #include "TH1F.h"
31 :
32 : #include "AliESDTrdTrack.h"
33 : #include "AliLog.h"
34 : #include "AliTRDgtuParam.h"
35 : #include "AliTRDtrackGTU.h"
36 : #include "AliTRDtrackletGTU.h"
37 : #include "AliTRDtrackletMCM.h"
38 : #include "AliESDTrdTrack.h"
39 :
40 48 : ClassImp(AliTRDtrackGTU)
41 :
42 : AliTRDtrackGTU::AliTRDtrackGTU() :
43 88 : TObject(),
44 88 : fStack(-1),
45 88 : fSector(-1),
46 88 : fPID(0),
47 88 : fTracklets(0x0),
48 88 : fTrackletMask(0),
49 88 : fNTracklets(0),
50 88 : fRefLayerIdx(-1),
51 88 : fZChannel(-1),
52 88 : fZSubChannel(-1),
53 88 : fA(0),
54 88 : fB(0),
55 88 : fC(0),
56 88 : fInvPtDev(0),
57 88 : fLabel(-1)
58 440 : {
59 : // default ctor
60 :
61 264 : fTracklets = new TClonesArray("AliTRDtrackletGTU", 6);
62 1232 : for (Int_t iTracklet = 0; iTracklet < 6; iTracklet++)
63 1584 : new ((*fTracklets)[iTracklet]) AliTRDtrackletGTU();
64 : // fTracklets->BypassStreamer(kFALSE);
65 176 : }
66 :
67 : AliTRDtrackGTU::AliTRDtrackGTU(const AliTRDtrackGTU &rhs) :
68 0 : TObject(),
69 0 : fStack(rhs.fStack),
70 0 : fSector(rhs.fSector),
71 0 : fPID(rhs.fPID),
72 0 : fTracklets(0x0),
73 0 : fTrackletMask(rhs.fTrackletMask),
74 0 : fNTracklets(rhs.fNTracklets),
75 0 : fRefLayerIdx(rhs.fRefLayerIdx),
76 0 : fZChannel(rhs.fZChannel),
77 0 : fZSubChannel(rhs.fZSubChannel),
78 0 : fA(rhs.fA),
79 0 : fB(rhs.fB),
80 0 : fC(rhs.fC),
81 0 : fInvPtDev(rhs.fInvPtDev),
82 0 : fLabel(rhs.fLabel)
83 0 : {
84 0 : fTracklets = new TClonesArray("AliTRDtrackletGTU", 6);
85 0 : for (Int_t iTracklet = 0; iTracklet < 6; iTracklet++)
86 0 : new ((*fTracklets)[iTracklet]) AliTRDtrackletGTU(*((AliTRDtrackletGTU*)(*(rhs.fTracklets))[iTracklet]));
87 0 : }
88 :
89 : AliTRDtrackGTU& AliTRDtrackGTU::operator=(const AliTRDtrackGTU &rhs)
90 : {
91 0 : if (&rhs != this) {
92 0 : TObject::operator=(rhs);
93 0 : fStack = rhs.fStack;
94 0 : fSector = rhs.fSector;
95 0 : fPID = rhs.fPID;
96 0 : fTrackletMask = rhs.fTrackletMask;
97 0 : fNTracklets = rhs.fNTracklets;
98 0 : fRefLayerIdx = rhs.fRefLayerIdx;
99 0 : fZChannel = rhs.fZChannel;
100 0 : fZSubChannel = rhs.fZSubChannel;
101 0 : fA = rhs.fA;
102 0 : fB = rhs.fB;
103 0 : fC = rhs.fC;
104 0 : fInvPtDev = rhs.fInvPtDev;
105 0 : fLabel = rhs.fLabel;
106 0 : for (Int_t iTracklet = 0; iTracklet < 6; iTracklet++)
107 0 : new ((*fTracklets)[iTracklet]) AliTRDtrackletGTU(*((AliTRDtrackletGTU*)(*(rhs.fTracklets))[iTracklet]));
108 0 : }
109 :
110 0 : return *this;
111 0 : }
112 :
113 : AliTRDtrackGTU& AliTRDtrackGTU::operator=(const AliESDTrdTrack &rhs)
114 : {
115 0 : if ((void*) &rhs != (void*) this) {
116 0 : TObject::operator=(rhs);
117 0 : fStack = rhs.GetStack();
118 0 : fSector = rhs.GetSector();
119 0 : fPID = rhs.GetPID();
120 0 : fTrackletMask = rhs.GetLayerMask();
121 0 : fNTracklets = 0;
122 0 : fRefLayerIdx = -1;
123 0 : fZChannel = -1;
124 0 : fZSubChannel = -1;
125 0 : fA = rhs.GetA();
126 0 : fB = rhs.GetB();
127 0 : fC = rhs.GetC();
128 0 : fLabel = rhs.GetLabel();
129 0 : for (Int_t iTracklet = 0; iTracklet < 6; iTracklet++) {
130 0 : AliTRDtrackletGTU *trkl = new ((*fTracklets)[iTracklet]) AliTRDtrackletGTU();
131 0 : if (fTrackletMask & (1 << iTracklet)) {
132 0 : ++fNTracklets;
133 0 : trkl->SetIndex(rhs.GetTrackletIndex(iTracklet));
134 0 : }
135 : else
136 0 : trkl->SetIndex(-1);
137 : }
138 0 : }
139 :
140 0 : return *this;
141 0 : }
142 :
143 : AliTRDtrackGTU::~AliTRDtrackGTU()
144 456 : {
145 : // dtor
146 :
147 76 : fTracklets->Delete();
148 152 : delete fTracklets;
149 228 : }
150 :
151 : void AliTRDtrackGTU::AddTracklet(const AliTRDtrackletGTU * const tracklet, Int_t layer)
152 : {
153 : // add a tracklet to this track
154 :
155 668 : if ( (fTrackletMask & (1 << layer)) != 0 ) {
156 0 : AliError(Form("Only one tracklet per layer (%i) possible! Mask: 0x%02x", layer, fTrackletMask));
157 0 : return;
158 : }
159 :
160 334 : new ((*fTracklets)[layer]) AliTRDtrackletGTU(*tracklet);
161 334 : fNTracklets++;
162 334 : fTrackletMask |= (1 << layer);
163 668 : }
164 :
165 : AliTRDtrackletGTU* AliTRDtrackGTU::GetTracklet(Int_t layer) const
166 : {
167 : // get a pointer to the tracklet in the layer specified
168 :
169 816 : if (IsTrackletInLayer(layer))
170 387 : return ((AliTRDtrackletGTU*) (*fTracklets)[layer]);
171 : else
172 21 : return 0x0;
173 408 : }
174 :
175 : Int_t AliTRDtrackGTU::GetNTracklets() const
176 : {
177 : // returns the number of tracklets in this track
178 :
179 68 : return fNTracklets;
180 : }
181 :
182 : Bool_t AliTRDtrackGTU::IsTrackletInLayer(Int_t layer) const
183 : {
184 : // checks for a tracklet in the given layer
185 :
186 1732 : if ( (GetTrackletMask() & (1 << layer)) != 0)
187 751 : return kTRUE;
188 : else
189 115 : return kFALSE;
190 866 : }
191 :
192 : void AliTRDtrackGTU::SetFitParams(Float_t a, Float_t b, Float_t c)
193 : {
194 : // set the fit parameters
195 :
196 30 : fA = a;
197 15 : fB = b;
198 15 : fC = c;
199 15 : }
200 :
201 : Int_t AliTRDtrackGTU::GetZSubChannel()
202 : {
203 : // returns the z-subchannel
204 :
205 142 : if (fZSubChannel < 0) {
206 350 : for (Int_t layer = 0; layer < AliTRDgtuParam::GetNLayers(); layer++) {
207 150 : if (IsTrackletInLayer(layer)) {
208 117 : AliTRDtrackletGTU *trkl = (AliTRDtrackletGTU*) (*fTracklets)[layer];
209 117 : if (trkl) {
210 209 : if ((fZSubChannel > -1) &&
211 92 : (fZSubChannel != trkl->GetSubChannel(GetZChannel())))
212 0 : AliError(Form("found inconsistent z-subchannels: track = %i/%i, trkl = %i",
213 : GetZChannel(), fZSubChannel, trkl->GetSubChannel(GetZChannel())));
214 117 : fZSubChannel = trkl->GetSubChannel(GetZChannel());
215 117 : }
216 : else {
217 0 : AliError("no tracklet where one should be according to layer mask");
218 : }
219 117 : }
220 : }
221 25 : }
222 71 : return fZSubChannel;
223 : }
224 :
225 : Int_t AliTRDtrackGTU::GetYapprox()
226 : {
227 : // returns an approximated y-position for the track
228 : // taken from the projected y-position of the tracklet in the reference layer
229 : // in which the track was found
230 :
231 30 : if ((fRefLayerIdx > -1) && (fRefLayerIdx < AliTRDgtuParam::GetNRefLayers()))
232 10 : return ((AliTRDtrackletGTU*) (*fTracklets)[AliTRDgtuParam::GetRefLayer(fRefLayerIdx)])->GetYProj();
233 : else
234 0 : return 0;
235 10 : }
236 :
237 : AliESDTrdTrack* AliTRDtrackGTU::CreateTrdTrack() const
238 : {
239 : // creates an AliESDTrdTrack to be added to the ESD
240 :
241 30 : AliESDTrdTrack *trk = new AliESDTrdTrack();
242 15 : trk->SetA((Int_t) fA);
243 15 : trk->SetB((Int_t) fB);
244 15 : trk->SetC((Int_t) fC);
245 15 : trk->SetLayerMask(fTrackletMask);
246 15 : trk->SetPID(fPID);
247 15 : trk->SetStack(fStack);
248 15 : trk->SetSector(fSector);
249 15 : trk->SetLabel(fLabel);
250 :
251 15 : if ( TMath::Abs(fInvPtDev) > AliTRDgtuParam::GetInvPtDevCut() )
252 4 : trk->SetFlags(1 << 5);
253 :
254 210 : for (Int_t iLayer = 0; iLayer < AliTRDgtuParam::GetNLayers(); iLayer++) {
255 90 : AliTRDtrackletGTU *trklGTU = GetTracklet(iLayer);
256 90 : if (trklGTU) {
257 73 : trk->SetTrackletIndex(trklGTU->GetIndex(), iLayer);
258 73 : AliESDTrdTracklet *trkl = trklGTU->GetTrackletESD();
259 73 : if (trkl)
260 0 : trk->AddTrackletReference(trkl, iLayer);
261 73 : }
262 : }
263 :
264 15 : return trk;
265 0 : }
266 :
267 : Bool_t AliTRDtrackGTU::CookLabel()
268 : {
269 : // assign label from tracklets according to frequency
270 :
271 : Int_t nLabels = 0;
272 30 : Int_t label[6] = { 0 };
273 15 : Int_t count[6] = { 0 };
274 :
275 210 : for (Int_t iTracklet = 0; iTracklet < 6; iTracklet++) {
276 90 : if (!IsTrackletInLayer(iTracklet))
277 : continue;
278 :
279 73 : Int_t currLabel = GetTracklet(iTracklet)->GetLabel();
280 :
281 : Bool_t assigned = kFALSE;
282 176 : for (Int_t iLabel = 0; iLabel < nLabels; iLabel++) {
283 61 : if (currLabel == label[iLabel]) {
284 55 : count[iLabel]++;
285 : assigned = kTRUE;
286 55 : break;
287 : }
288 : }
289 :
290 73 : if (!assigned) {
291 18 : label[nLabels] = currLabel;
292 18 : count[nLabels] = 1;
293 18 : nLabels++;
294 18 : }
295 73 : }
296 :
297 15 : Int_t index[32];
298 15 : TMath::Sort(6, count, index);
299 15 : fLabel = label[index[0]];
300 :
301 15 : return kTRUE;
302 15 : }
|