1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 /* $Id: AliTRDgtuTMU.cxx 28397 2008-09-02 09:33:00Z cblume $ */
18 ////////////////////////////////////////////////////////////////////////////
20 // Track Matching Unit (TMU) simulation //
22 // Author: J. Klein (Jochen.Klein@cern.ch) //
24 ////////////////////////////////////////////////////////////////////////////
31 #include "AliESDEvent.h"
32 #include "AliESDTrdTrack.h"
35 #include "AliTRDgeometry.h"
36 #include "AliTRDpadPlane.h"
38 #include "AliTRDgtuParam.h"
39 #include "AliTRDgtuTMU.h"
40 #include "AliTRDtrackGTU.h"
42 ClassImp(AliTRDgtuTMU)
44 AliTRDgtuTMU::AliTRDgtuTMU(Int_t stack, Int_t sector) :
47 fZChannelTracklets(0x0),
53 fGtuParam = AliTRDgtuParam::Instance();
54 fTracklets = new TObjArray*[fGtuParam->GetNLayers()];
55 fZChannelTracklets = new TList*[fGtuParam->GetNLayers()];
56 for (Int_t layer = 0; layer < fGtuParam->GetNLayers(); layer++) {
57 fTracklets[layer] = new TObjArray();
58 fZChannelTracklets[layer] = new TList[fGtuParam->GetNZChannels()];
60 fTracks = new TList*[fGtuParam->GetNZChannels()];
61 for (Int_t zch = 0; zch < fGtuParam->GetNZChannels(); zch++) {
62 fTracks[zch] = new TList[fGtuParam->GetNRefLayers()];
70 AliTRDgtuTMU::~AliTRDgtuTMU()
72 for (Int_t zch = 0; zch < fGtuParam->GetNZChannels(); zch++) {
73 delete [] fTracks[zch];
76 for (Int_t layer = 0; layer < fGtuParam->GetNLayers(); layer++) {
77 fTracklets[layer]->Delete();
78 delete [] fZChannelTracklets[layer];
79 delete fTracklets[layer];
81 delete [] fZChannelTracklets;
85 Bool_t AliTRDgtuTMU::SetSector(Int_t sector)
89 if (sector > -1 && sector < fGtuParam->GetGeo()->Nsector() ) {
94 AliError(Form("Invalid sector given: %i", sector));
98 Bool_t AliTRDgtuTMU::SetStack(Int_t stack)
100 // Set the stack (necessary for tracking)
102 if (stack > -1 && stack < fGtuParam->GetGeo()->Nstack() ) {
107 AliError(Form("Invalid stack given: %i", stack));
111 Bool_t AliTRDgtuTMU::Reset()
115 for (Int_t zch = 0; zch < fGtuParam->GetNZChannels(); zch++) {
116 for (Int_t reflayeridx = 0; reflayeridx < fGtuParam->GetNRefLayers(); reflayeridx++) {
117 fTracks[zch][reflayeridx].Clear();
121 // delete all added tracklets
122 for (Int_t layer = 0; layer < fGtuParam->GetNLinks()/2; layer++) {
123 for (Int_t zch = 0; zch < fGtuParam->GetNZChannels(); zch++)
124 fZChannelTracklets[layer][zch].Clear();
125 fTracklets[layer]->Delete();
134 Bool_t AliTRDgtuTMU::AddTracklet(AliTRDtrackletBase *tracklet, Int_t link)
136 // add a tracklet on the given link
138 AliTRDtrackletGTU *trkl = new AliTRDtrackletGTU(tracklet);
139 fTracklets[(Int_t) link/2]->Add(trkl);
144 Bool_t AliTRDgtuTMU::RunTMU(TList *ListOfTracks, AliESDEvent *esd)
146 // performs the analysis as in a TMU module of the GTU, i. e.
148 // calculation of track parameteres (pt, deflection, ???)
150 if (fStack < 0 || fSector < 0) {
151 AliError("No valid stack/sector set for this TMU! Tracking aborted!");
155 // ----- Input units -----
156 AliInfo("--------- Running Input units ----------");
157 for (Int_t layer = 0; layer < fGtuParam->GetNLayers(); layer++) {
158 if (!RunInputUnit(layer)) {
159 AliError(Form("Input unit in layer %i failed", layer));
164 // ----- Z-channel units -----
165 AliInfo("--------- Running Z-channel units ----------");
166 for (Int_t layer = 0; layer < fGtuParam->GetNLayers(); layer++) {
167 fZChannelTracklets[layer] = new TList[fGtuParam->GetNZChannels()];
168 if (!RunZChannelUnit(layer)) {
169 AliError(Form("Z-Channel unit in layer %i failed", layer));
174 // ----- track finding -----
175 AliInfo("--------- Running tracking units ----------");
176 for (Int_t zch = 0; zch < fGtuParam->GetNZChannels(); zch++) {
177 AliInfo(Form("Track finder for Zchannel: %i", zch));
178 if (!RunTrackFinder(zch, ListOfTracks)) {
179 AliError(Form("Track Finder in z-channel %i failed", zch));
184 // ----- Track Merging -----
185 if (!RunTrackMerging(ListOfTracks)) {
186 AliError("Track merging failed");
190 // ----- track reconstruction -----
191 if (!RunTrackReconstruction(ListOfTracks)) {
192 AliError("Track reconstruction failed");
197 TIter next(ListOfTracks);
198 while (AliTRDtrackGTU *trk = (AliTRDtrackGTU*) next()) {
199 AliESDTrdTrack *trdtrack = trk->CreateTrdTrack();
200 esd->AddTrdTrack(trdtrack);
208 Bool_t AliTRDgtuTMU::RunInputUnit(Int_t layer)
210 // precalculations for the tracking and reconstruction
212 fTracklets[layer]->Sort();
213 TIter next(fTracklets[layer]);
215 while ( AliTRDtrackletGTU *trk = (AliTRDtrackletGTU*) next() ) {
216 trk->SetIndex(fTracklets[layer]->IndexOf(trk));
218 Int_t alpha = (trk->GetYbin() >> fGtuParam->GetBitExcessY()) * fGtuParam->GetCiAlpha(layer);
219 alpha = ( 2 * trk->GetdY() - (alpha >> fGtuParam->GetBitExcessAlpha()) + 1 ) >> 1;
220 trk->SetAlpha(alpha);
222 Int_t yproj = trk->GetdY() * (fGtuParam->GetCiYProj(layer));
223 yproj = ( ( ( (yproj >> fGtuParam->GetBitExcessYProj()) + trk->GetYbin() ) >> 2) + 1) >> 1;
224 trk->SetYProj(yproj);
226 trk->SetYPrime(trk->GetYbin() + fGtuParam->GetYt(fStack, layer, trk->GetZbin()));
228 // printf("InputUnit : GetIndex(): %3i, GetZbin(): %2i, GetY() : %5i, GetdY() : %3i, GetYPrime() : %5i, GetYProj() : %5i, GetAlpha() : %3i \n",
229 // trk->GetIndex(), trk->GetZbin(), trk->GetYbin(), trk->GetdY(), trk->GetYPrime(), trk->GetYProj(), trk->GetAlpha() );
234 Bool_t AliTRDgtuTMU::RunZChannelUnit(Int_t layer)
236 // run the z-channel unit
238 TIter next(fTracklets[layer]);
240 while (AliTRDtrackletGTU *trk = (AliTRDtrackletGTU*) next()) {
241 printf("*TMU* Tracklet in stack %d, layer %2d: 0x%08x ", fStack, layer, trk->GetTrackletWord());
242 for (Int_t zch = 0; zch < fGtuParam->GetNZChannels(); zch++) {
243 if (fGtuParam->IsInZChannel(fStack, layer, zch, trk->GetZbin()) ) {
244 trk->SetSubChannel(zch, fGtuParam->GetZSubchannel(fStack, layer, zch, trk->GetZbin()) );
245 printf("Z%i(%i) ", zch, trk->GetSubChannel(zch));
247 TIter nexttrkl(&fZChannelTracklets[layer][zch], kIterBackward);
248 AliTRDtrackletGTU *t = 0x0;
249 while ((t = (AliTRDtrackletGTU*) nexttrkl.Next())) {
250 if (t->GetSubChannel(zch) < trk->GetSubChannel(zch) ||
251 (t->GetSubChannel(zch) == trk->GetSubChannel(zch) && t->GetYProj() < trk->GetYProj()) )
254 fZChannelTracklets[layer][zch].AddAfter(t, trk);
264 Bool_t AliTRDgtuTMU::RunTrackFinder(Int_t zch, TList* /* ListOfTracks */)
266 // run the track finding
268 Int_t *notr = new Int_t[fGtuParam->GetNLayers()];
269 Int_t *ptrA = new Int_t[fGtuParam->GetNLayers()];
270 Int_t *ptrB = new Int_t[fGtuParam->GetNLayers()];
272 Bool_t *bHitA = new Bool_t[fGtuParam->GetNLayers()];
273 Bool_t *bHitB = new Bool_t[fGtuParam->GetNLayers()];
274 Bool_t *bAligned = new Bool_t[fGtuParam->GetNLayers()];
275 Bool_t *bAlignedA = new Bool_t[fGtuParam->GetNLayers()];
276 Bool_t *bAlignedB = new Bool_t[fGtuParam->GetNLayers()];
277 Bool_t *bDone = new Bool_t[fGtuParam->GetNLayers()];
280 Int_t *inc = new Int_t[fGtuParam->GetNLayers()];
281 Int_t *incprime = new Int_t[fGtuParam->GetNLayers()];
283 // ----- signals within current layer -----
294 AliTRDtrackletGTU *trkRA = 0x0; // reference tracklet A
295 AliTRDtrackletGTU *trkRB = 0x0; // reference tracklet B
296 AliTRDtrackletGTU *trkA = 0x0; // tracklet A in current layer
297 AliTRDtrackletGTU *trkB = 0x0; // tracklet B in current layer
299 AliTRDtrackletGTU *trkEnd = new AliTRDtrackletGTU();
300 for (Int_t i = 0; i < fGtuParam->GetNZChannels(); i++)
301 trkEnd->SetSubChannel(i, 7);
306 for (Int_t refLayerIdx = 0; refLayerIdx < fGtuParam->GetNRefLayers(); refLayerIdx++) {
307 Int_t reflayer = fGtuParam->GetRefLayer(refLayerIdx);
308 AliInfo(Form("~~~~~ Reflayer: %i", reflayer));
310 ready = kFALSE; // ready if all channels done
312 for (Int_t layer = 0; layer < fGtuParam->GetNLayers(); layer++) {
313 notr[layer] = fZChannelTracklets[layer][zch].GetEntries();
314 ptrA[layer] = 0; // notr[layer] > 0 ? 0 : -1;
315 ptrB[layer] = 1; // notr[layer] > 1 ? 1 : -1;
317 // not necessary here
318 // bDone[layer] = (ptrB[layer] >= notr[layer] - 1); // trkB is last one
319 // bDone[layer] = (notr[layer] <= 0);
320 // ready = ready && bDone[layer];
323 AliInfo(Form("in layer: %i (zchannel = %i) there are: %i tracklets", layer, zch, notr[layer]));
326 if (ptrA[reflayer] < 0 && ptrB[reflayer] < 0)
330 // ----- reference tracklets -----
333 if (0 <= ptrA[reflayer] && ptrA[reflayer] < notr[reflayer])
334 trkRA = (AliTRDtrackletGTU*) fZChannelTracklets[reflayer][zch].At(ptrA[reflayer]);
336 AliInfo(Form("No valid tracklet in the reference at ptr: %i! Aborting!", ptrA[reflayer]));
340 if (0 <= ptrB[reflayer] && ptrB[reflayer] < notr[reflayer])
341 trkRB = (AliTRDtrackletGTU*) fZChannelTracklets[reflayer][zch].At(ptrB[reflayer]);
343 AliInfo(Form("ptrRA: %i, ptrRB: %i", ptrA[reflayer], ptrB[reflayer]));
344 Yplus = trkRA->GetYProj() + fGtuParam->GetDeltaY();
345 Yminus = trkRA->GetYProj() - fGtuParam->GetDeltaY();
346 Alphaplus = trkRA->GetAlpha() + fGtuParam->GetDeltaAlpha();
347 Alphaminus = trkRA->GetAlpha() - fGtuParam->GetDeltaAlpha();
349 YBplus = trkRB->GetYProj() + fGtuParam->GetDeltaY();
350 YBminus = trkRB->GetYProj() - fGtuParam->GetDeltaY();
352 else { // irrelevant (should be, is it?)
353 YBplus = trkRA->GetYProj() + fGtuParam->GetDeltaY();
354 YBminus = trkRA->GetYProj() - fGtuParam->GetDeltaY();
361 for (Int_t layer = 0; layer < fGtuParam->GetNLayers(); layer++) {
362 bHitA[layer] = bHitB[layer] = bAligned[layer] = kFALSE;
363 inc[layer] = incprime[layer] = 0;
365 if (layer == reflayer) {
366 bHitA[layer] = kTRUE;
367 bAligned[layer] = kTRUE;
374 if (0 <= ptrA[layer] && ptrA[layer] < notr[layer])
375 trkA = (AliTRDtrackletGTU*) fZChannelTracklets[layer][zch].At(ptrA[layer]);
376 if (0 <= ptrB[layer] && ptrB[layer] < notr[layer])
377 trkB = (AliTRDtrackletGTU*) fZChannelTracklets[layer][zch].At(ptrB[layer]);
379 bAlignedA[layer] = kFALSE;
380 bAlignedB[layer] = kFALSE;
383 bHitA[layer] = ( !(trkA->GetSubChannel(zch) < trkRA->GetSubChannel(zch) || (trkA->GetSubChannel(zch) == trkRA->GetSubChannel(zch) && trkA->GetYProj() < Yminus) ) &&
384 !(trkA->GetSubChannel(zch) > trkRA->GetSubChannel(zch) || (trkA->GetSubChannel(zch) == trkRA->GetSubChannel(zch) && trkA->GetYProj() > Yplus) ) &&
385 !(trkA->GetAlpha() < Alphaminus) &&
386 !(trkA->GetAlpha() > Alphaplus) );
387 bAlignedA[layer] = !(trkA->GetSubChannel(zch) < trkRA->GetSubChannel(zch) || (trkA->GetSubChannel(zch) == trkRA->GetSubChannel(zch) && trkA->GetYProj() < Yminus) );
391 bAlignedA[layer] = kTRUE;
395 bHitB[layer] = ( !(trkB->GetSubChannel(zch) < trkRA->GetSubChannel(zch) || (trkB->GetSubChannel(zch) == trkRA->GetSubChannel(zch) && trkB->GetYProj() < Yminus) ) &&
396 !(trkB->GetSubChannel(zch) > trkRA->GetSubChannel(zch) || (trkB->GetSubChannel(zch) == trkRA->GetSubChannel(zch) && trkB->GetYProj() > Yplus) ) &&
397 !(Alphaminus > trkB->GetAlpha()) &&
398 !(Alphaplus > trkB->GetAlpha()) );
399 bAlignedB[layer] = (trkB->GetSubChannel(zch) > trkRA->GetSubChannel(zch) || (trkB->GetSubChannel(zch) == trkRA->GetSubChannel(zch) && trkB->GetYProj() > Yplus) );
403 bAlignedB[layer] = kTRUE;
406 bAligned[layer] = bAlignedA[layer] || bAlignedB[layer]; //???
407 // bAligned[layer] = bAlignedA[layer]; //???
409 if (bAligned[layer] && (bHitA[layer] || bHitB[layer]) )
411 else if (!bAligned[layer] )
415 if ((trkA->GetSubChannel(zch) > trkRB->GetSubChannel(zch)) || (trkA->GetSubChannel(zch) == trkRB->GetSubChannel(zch) && trkA->GetYProj() > YBplus) )
422 // pre-calculation for the layer shifting (alignment w. r. t. trkRB)
425 if ((trkA->GetSubChannel(zch) < trkRB->GetSubChannel(zch)) || (trkA->GetSubChannel(zch) == trkRB->GetSubChannel(zch) && trkA->GetYProj() < YBminus )) // could trkA be aligned for trkRB
435 if ((trkB->GetSubChannel(zch) < trkRB->GetSubChannel(zch)) || (trkB->GetSubChannel(zch) == trkRB->GetSubChannel(zch) && trkB->GetYProj() < YBminus )) // could trkB be aligned for trkRB
442 } // end of loop over layers
444 AliInfo(Form("logic calculation finished, Nhits: %i", NHits));
447 // ----- track registration -----
448 AliInfo("***** TMU: Track found *****");
449 AliTRDtrackGTU *track = new AliTRDtrackGTU();
450 track->SetSector(fSector);
451 track->SetStack(fStack);
452 for (Int_t layer = 0; layer < fGtuParam->GetNLayers(); layer++ ) {
453 if (bHitA[layer] || layer == reflayer)
454 track->AddTracklet((AliTRDtrackletGTU* ) fZChannelTracklets[layer][zch].At(ptrA[layer]), layer );
455 else if (bHitB[layer])
456 track->AddTracklet((AliTRDtrackletGTU* ) fZChannelTracklets[layer][zch].At(ptrB[layer]), layer );
459 Bool_t registerTrack = kTRUE;
460 for (Int_t layerIdx = refLayerIdx; layerIdx > 0; layerIdx--) {
461 if (track->IsTrackletInLayer(fGtuParam->GetRefLayer(layerIdx)))
462 registerTrack = kFALSE;
465 track->SetZChannel(zch);
466 track->SetRefLayerIdx(refLayerIdx);
467 fTracks[zch][refLayerIdx].Add(track);
471 if ( (NUnc != 0) && (NUnc + NHits >= 4) ) // could this position of the reference layer give some track //??? special check in case of hit?
473 else if (NWayBeyond > 2) // no track possible for both reference tracklets
478 if (inc[reflayer] != 0) // reflayer is shifted
479 for (Int_t layer = 0; layer < fGtuParam->GetNLayers(); layer++) {
480 if (layer == reflayer)
482 inc[layer] = incprime[layer];
484 else { // reflayer is not shifted
485 for (Int_t layer = 0; layer < fGtuParam->GetNLayers(); layer++) {
486 if (layer == reflayer || notr[layer] == 0)
491 if (0 <= ptrA[layer] && ptrA[layer] < notr[layer])
492 trkA = (AliTRDtrackletGTU*) fZChannelTracklets[layer][zch].At(ptrA[layer]);
494 if (0 <= ptrB[layer] && ptrB[layer] < notr[layer])
495 trkB = (AliTRDtrackletGTU*) fZChannelTracklets[layer][zch].At(ptrB[layer]);
498 if ( !(trkA->GetSubChannel(zch) < trkRA->GetSubChannel(zch) || (trkA->GetSubChannel(zch) == trkRA->GetSubChannel(zch) && trkA->GetYProj() < Yminus) ) &&
499 !(trkA->GetSubChannel(zch) > trkRA->GetSubChannel(zch) || (trkA->GetSubChannel(zch) == trkRA->GetSubChannel(zch) && trkA->GetYProj() > Yplus ) ) ) // trkA could hit trkRA
502 if ( trkB->GetSubChannel(zch) < trkRA->GetSubChannel(zch) || (trkB->GetSubChannel(zch) == trkRA->GetSubChannel(zch) && trkB->GetYProj() < Yminus) )
504 else if ( !(trkB->GetSubChannel(zch) > trkRA->GetSubChannel(zch) || (trkB->GetSubChannel(zch) == trkRA->GetSubChannel(zch) && trkB->GetYProj() > Yplus) ) )
507 inc[layer] = incprime[layer];
510 inc[layer] = incprime[layer];
516 for (Int_t layer = 0; layer < fGtuParam->GetNLayers(); layer++) {
518 bDone[layer] = ptrB[layer] < 0 || ptrB[layer] >= notr[layer];
519 ready = ready && bDone[layer];
521 if (inc[layer] != 0 && ptrA[layer] >= notr[layer])
522 AliError(Form("Invalid increment: %i at ptrA: %i, notr: %i", inc[layer], ptrA[layer], notr[layer]));
524 // AliInfo(Form("Shifting layer: %i, notr: %i, ptrA: %i, ptrB: %i, inc: %i", layer, notr[layer], ptrA[layer], ptrB[layer], inc[layer]));
525 AliInfo(Form(" -- Layer: %i %i %i +%i %s (no: %i)", layer, ptrA[layer], ptrB[layer], inc[layer], bDone[layer] ? "done" : " ", notr[layer]));
526 ptrA[layer] += inc[layer];
527 ptrB[layer] += inc[layer];
531 } // end of loop over reflayer
548 Bool_t AliTRDgtuTMU::RunTrackMerging(TList* ListOfTracks)
550 TList **tracksRefMerged = new TList*[fGtuParam->GetNZChannels()];
551 TList **tracksRefUnique = new TList*[fGtuParam->GetNZChannels()];
553 for (Int_t zch = 0; zch < fGtuParam->GetNZChannels(); zch++) {
554 tracksRefMerged[zch] = new TList;
555 tracksRefUnique[zch] = new TList;
558 TList *tracksZMergedStage0 = new TList;
559 TList *tracksZUniqueStage0 = new TList;
561 TList **tracksZSplitted = new TList*[2];
562 for (Int_t i = 0; i < 2; i++)
563 tracksZSplitted[i] = new TList;
565 TList *tracksZMergedStage1 = new TList;
567 AliTRDtrackGTU **trkInRefLayer = new AliTRDtrackGTU*[fGtuParam->GetNRefLayers()];
569 Bool_t done = kFALSE;
571 AliTRDtrackGTU *trkStage0 = 0x0;
573 for (Int_t zch = 0; zch < fGtuParam->GetNZChannels(); zch++) {
574 // ----- Merging and Unification in Reflayers -----
578 for (Int_t refLayerIdx = 0; refLayerIdx < fGtuParam->GetNRefLayers(); refLayerIdx++) {
579 trkInRefLayer[refLayerIdx] = (AliTRDtrackGTU*) fTracks[zch][refLayerIdx].First();
580 if (trkInRefLayer[refLayerIdx] == 0) {
583 else if (trkStage0 == 0x0 ) {
584 trkStage0 = trkInRefLayer[refLayerIdx];
585 minIdx = refLayerIdx;
588 else if (trkInRefLayer[refLayerIdx]->GetZSubChannel() < trkStage0->GetZSubChannel() ||
589 (trkInRefLayer[refLayerIdx]->GetZSubChannel() == trkStage0->GetZSubChannel() && trkInRefLayer[refLayerIdx]->GetYapprox() < trkStage0->GetYapprox()) ) {
590 minIdx = refLayerIdx;
591 trkStage0 = trkInRefLayer[refLayerIdx];
596 tracksRefMerged[zch]->Add(trkStage0);
597 fTracks[zch][minIdx].RemoveFirst();
598 } while (trkStage0 != 0);
600 Uniquifier(tracksRefMerged[zch], tracksRefUnique[zch]);
603 // ----- Merging in zchannels - 1st stage -----
608 for (Int_t zch = 0; zch < fGtuParam->GetNZChannels(); zch++) {
609 AliTRDtrackGTU *trk = (AliTRDtrackGTU*) tracksRefUnique[zch]->First();
613 else if (trkStage0 == 0x0 ) {
618 else if ( ((trk->GetZChannel() + 3 * trk->GetZSubChannel()) / 2 - 1) < ((trkStage0->GetZChannel() + 3 * trkStage0->GetZSubChannel()) / 2 -1 ) ) {
626 tracksZMergedStage0->Add(trkStage0);
627 tracksRefUnique[minIdx]->RemoveFirst();
628 } while (trkStage0 != 0);
630 Uniquifier(tracksZMergedStage0, tracksZUniqueStage0);
632 // ----- Splitting in z -----
634 TIter next(tracksZUniqueStage0);
635 while (AliTRDtrackGTU *trk = (AliTRDtrackGTU*) next()) {
636 tracksZSplitted[(trk->GetZChannel() + 3 * (trk->GetZSubChannel() - 1)) % 2]->Add(trk);
639 // ----- Merging in zchanels - 2nd stage -----
644 for (Int_t i = 0; i < 2; i++) {
645 AliTRDtrackGTU *trk = (AliTRDtrackGTU*) tracksZSplitted[i]->First();
649 else if (trkStage0 == 0x0 ) {
654 else if ( ((trk->GetZChannel() + 3 * (trk->GetZSubChannel() - 1)) / 2) < ((trkStage0->GetZChannel() + 3 * (trkStage0->GetZSubChannel() - 1)) / 2) ) {
661 tracksZMergedStage1->Add(trkStage0);
662 tracksZSplitted[minIdx]->RemoveFirst();
663 } while (trkStage0 != 0);
665 Uniquifier(tracksZMergedStage1, ListOfTracks);
670 Bool_t AliTRDgtuTMU::RunTrackReconstruction(TList* ListOfTracks)
672 TIter next(ListOfTracks);
673 while (AliTRDtrackGTU *track = (AliTRDtrackGTU*) next()) {
674 CalculateTrackParams(track);
679 Bool_t AliTRDgtuTMU::CalculateTrackParams(AliTRDtrackGTU *track)
681 // calculate the track parameters
684 AliError("No track to calculate!");
694 AliInfo(Form("There are %i tracklets in this track.", track->GetNTracklets()));
696 for (Int_t layer = 0; layer < fGtuParam->GetNLayers(); layer++) {
697 if (!track->IsTrackletInLayer(layer)) {
700 AliTRDtrackletGTU *trk = track->GetTracklet(layer);
702 AliError(Form("Could not get tracklet in layer %i\n", layer));
705 AliInfo(Form("trk yprime: %i", trk->GetYPrime()));
706 a += (((Int_t) (2048 * fGtuParam->GetAki(track->GetTrackletMask(), layer))) * trk->GetYPrime() + 1) >> 8;
707 b += fGtuParam->GetBki(track->GetTrackletMask(), layer) * trk->GetYPrime() * fGtuParam->GetBinWidthY();
708 c += fGtuParam->GetCki(track->GetTrackletMask(), layer) * trk->GetYPrime() * fGtuParam->GetBinWidthY();
714 fGtuParam->GetIntersectionPoints(track->GetTrackletMask(), x1, x2);
715 AliInfo(Form("Intersection points: %f, %f", x1, x2));
716 AliInfo(Form("Sum: a = %5i, b = %9.2f, c = %9.2f\n", a, b, c));
717 track->SetFitParams(a, b, c);
719 Float_t r = fGtuParam->GetRadius(a, b, x1, x2);
720 Int_t pt = (Int_t) (2 * r);
727 AliInfo(Form("Track parameters: a = %i, b = %f, c = %f, x1 = %f, x2 = %f, r = %f, pt = %f (trkl mask: %i)", a, b, c, x1, x2, r, track->GetPt(), track->GetTrackletMask()));
731 Bool_t AliTRDgtuTMU::WriteTrackletsToTree(TTree *trklTree)
734 AliError("No tree given");
737 AliTRDtrackletGTU *trkl = 0x0;
738 TBranch *branch = trklTree->GetBranch("gtutracklets");
740 branch = trklTree->Branch("gtutracklets", "AliTRDtrackletGTU", &trkl, 32000, 99);
743 AliInfo(Form("---------- Writing tracklets to tree (not yet) ----------"));
744 for (Int_t layer = 0; layer < fGtuParam->GetNLayers(); layer++) {
745 TIter next(fTracklets[layer]);
746 while ((trkl = (AliTRDtrackletGTU*) next())) {
747 AliInfo(Form("InputUnit : GetIndex(): %3i, GetZbin(): %2i, GetY() : %5i, GetdY() : %3i, GetYPrime() : %5i, GetYProj() : %5i, GetAlpha() : %3i, Zidx(2..0): %i %i %i", trkl->GetIndex(), trkl->GetZbin(), trkl->GetYbin(), trkl->GetdY(), trkl->GetYPrime(), trkl->GetYProj(), trkl->GetAlpha(), trkl->GetSubChannel(2), trkl->GetSubChannel(1), trkl->GetSubChannel(0) ));
748 branch->SetAddress(&trkl);
755 Bool_t AliTRDgtuTMU::Uniquifier(TList *inlist, TList *outlist)
758 AliTRDtrackGTU *trkStage0 = 0x0;
759 AliTRDtrackGTU *trkStage1 = 0x0;
762 trkStage0 = (AliTRDtrackGTU*) next();
764 Bool_t tracks_equal = kFALSE;
765 if (trkStage0 != 0 && trkStage1 != 0) {
766 for (Int_t layer = 0; layer < fGtuParam->GetNLayers(); layer++) {
767 if (trkStage0->GetTrackletIndex(layer) != -1 && trkStage0->GetTrackletIndex(layer) == trkStage1->GetTrackletIndex(layer)) {
768 tracks_equal = kTRUE;
775 if (trkStage0->GetNTracklets() > trkStage1->GetNTracklets())
776 trkStage1 = trkStage0;
779 if (trkStage1 != 0x0)
780 outlist->Add(trkStage1);
781 trkStage1 = trkStage0;
784 } while (trkStage1 != 0x0);