]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDtrigger.cxx
First upload of a macro to create raw data manually
[u/mrichter/AliRoot.git] / TRD / AliTRDtrigger.cxx
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 //  TRD trigger class                                                        //
19 //                                                                           //
20 //  Author:                                                                  //
21 //    Bogdan Vulpescu                                                        //
22 //                                                                           //
23 ///////////////////////////////////////////////////////////////////////////////
24
25 #include <TTree.h>
26 #include <TBranch.h>
27 #include <TMatrixD.h>
28 #include <TClonesArray.h>
29 #include <TObjArray.h>
30
31 #include "AliLog.h"
32 #include "AliRun.h"
33 #include "AliLoader.h"
34
35 #include "AliTRDdigitsManager.h"
36 #include "AliTRDgeometry.h"
37 #include "AliTRDdataArrayI.h"
38 #include "AliTRDcalibDB.h"
39 #include "AliTRDCommonParam.h"
40 #include "AliTRDrawData.h"
41 #include "AliTRDtrigger.h"
42 #include "AliTRDmodule.h"
43 #include "AliTRDmcmTracklet.h"
44 #include "AliTRDgtuTrack.h"
45 #include "AliTRDtrigParam.h"
46 #include "AliTRDmcm.h"
47 #include "AliTRDzmaps.h"
48 #include "AliTRDCalibra.h"
49 #include "Cal/AliTRDCalPIDLQ.h"
50
51 ClassImp(AliTRDtrigger)
52
53 //_____________________________________________________________________________
54 AliTRDtrigger::AliTRDtrigger()
55   :TNamed()
56   ,fField(0)
57   ,fGeo(NULL)
58   ,fCalib(NULL)
59   ,fCParam(NULL)
60   ,fTrigParam(NULL)
61   ,fRunLoader(NULL)
62   ,fDigitsManager(NULL)
63   ,fTrackletTree(NULL)
64   ,fTracklets(NULL)
65   ,fNROB(0)
66   ,fMCM(NULL)
67   ,fTrk(NULL)
68   ,fTrkTest(NULL)
69   ,fModule(NULL)
70   ,fGTUtrk(NULL)
71   ,fNtracklets(0)
72   ,fDigits(NULL)
73   ,fTrack0(NULL)
74   ,fTrack1(NULL)
75   ,fTrack2(NULL)
76   ,fNPrimary(0)
77   ,fTracks(NULL)
78 {
79   //
80   // AliTRDtrigger default constructor
81   //
82
83 }
84
85 //_____________________________________________________________________________
86 AliTRDtrigger::AliTRDtrigger(const Text_t *name, const Text_t *title)
87   :TNamed(name,title)
88   ,fField(0)
89   ,fGeo(NULL)
90   ,fCalib(NULL)
91   ,fCParam(NULL)
92   ,fTrigParam(NULL)
93   ,fRunLoader(NULL)
94   ,fDigitsManager(new AliTRDdigitsManager())
95   ,fTrackletTree(NULL)
96   ,fTracklets(new TObjArray(400))
97   ,fNROB(0)
98   ,fMCM(NULL)
99   ,fTrk(NULL)
100   ,fTrkTest(NULL)
101   ,fModule(NULL)
102   ,fGTUtrk(NULL)
103   ,fNtracklets(0)
104   ,fDigits(NULL)
105   ,fTrack0(NULL)
106   ,fTrack1(NULL)
107   ,fTrack2(NULL)
108   ,fNPrimary(0)
109   ,fTracks(new TClonesArray("AliTRDgtuTrack",1000))
110 {
111   //
112   // AliTRDtrigger constructor
113   //
114
115 }
116
117 //_____________________________________________________________________________
118 AliTRDtrigger::AliTRDtrigger(const AliTRDtrigger &p)
119   :TNamed(p)
120   ,fField(p.fField)
121   ,fGeo(NULL)
122   ,fCalib(NULL)
123   ,fCParam(NULL)
124   ,fTrigParam(NULL)
125   ,fRunLoader(NULL)
126   ,fDigitsManager(NULL)
127   ,fTrackletTree(NULL)
128   ,fTracklets(NULL)
129   ,fNROB(p.fNROB)
130   ,fMCM(NULL)
131   ,fTrk(NULL)
132   ,fTrkTest(NULL)
133   ,fModule(NULL)
134   ,fGTUtrk(NULL)
135   ,fNtracklets(p.fNtracklets)
136   ,fDigits(NULL)
137   ,fTrack0(NULL)
138   ,fTrack1(NULL)
139   ,fTrack2(NULL)
140   ,fNPrimary(p.fNPrimary)
141   ,fTracks(NULL)
142 {
143   //
144   // AliTRDtrigger copy constructor
145   //
146
147 }
148
149 ///_____________________________________________________________________________
150 AliTRDtrigger::~AliTRDtrigger()
151 {
152   //
153   // AliTRDtrigger destructor
154   //
155
156   if (fTracklets) {
157     fTracklets->Delete();
158     delete fTracklets;
159   }
160
161   if (fTracks) {
162     fTracks->Delete();
163     delete fTracks;
164   }
165
166 }
167
168 //_____________________________________________________________________________
169 AliTRDtrigger &AliTRDtrigger::operator=(const AliTRDtrigger &p)
170 {
171   //
172   // Assignment operator
173   //
174
175   if (this != &p) ((AliTRDtrigger &) p).Copy(*this);
176   return *this;
177
178 }
179
180 //_____________________________________________________________________________
181 void AliTRDtrigger::Copy(TObject &) const
182 {
183   //
184   // Copy function
185   //
186
187   AliFatal("Not implemented");
188
189 }
190
191 //_____________________________________________________________________________
192 void AliTRDtrigger::Init()
193 {
194
195   fModule = new AliTRDmodule(fTrigParam); 
196   fTracks->Clear();
197
198   fField  = fTrigParam->GetField();
199   fGeo    = (AliTRDgeometry*)AliTRDgeometry::GetGeometry(fRunLoader);
200
201   fCalib  = AliTRDcalibDB::Instance();
202   if (!fCalib) {
203     AliError("No instance of AliTRDcalibDB.");
204     return;  
205   }
206
207   fCParam = AliTRDCommonParam::Instance();
208   if (!fCParam) {
209     AliError("No common parameters.");
210     return;
211   }
212
213 }
214
215 //_____________________________________________________________________________
216 Bool_t AliTRDtrigger::Open(const Char_t *name, Int_t nEvent)
217 {
218   //
219   // Opens the AliROOT file.
220   //
221
222   TString evfoldname = AliConfig::GetDefaultEventFolderName();
223   fRunLoader         = AliRunLoader::GetRunLoader(evfoldname);
224
225   if (!fRunLoader) {
226     fRunLoader = AliRunLoader::Open(name);
227   }
228   if (!fRunLoader) {
229     AliError(Form("Can not open session for file %s.",name));
230     return kFALSE;
231   }
232
233   // Import the Trees for the event nEvent in the file
234   fRunLoader->GetEvent(nEvent);
235
236   // Open output
237   TObjArray *ioArray = 0;
238   AliLoader* loader  = fRunLoader->GetLoader("TRDLoader");
239   loader->MakeTree("T");
240   fTrackletTree = loader->TreeT();
241   fTrackletTree->Branch("TRDmcmTracklet","TObjArray",&ioArray,32000,0);
242   Init();
243
244   return kTRUE;
245
246 }
247
248 //_____________________________________________________________________________
249 Bool_t AliTRDtrigger::ReadDigits() 
250 {
251   //
252   // Reads the digits arrays from the input aliroot file
253   //
254
255   if (!fRunLoader) {
256     AliError("Can not find the Run Loader");
257     return kFALSE;
258   }
259
260   AliLoader* loader = fRunLoader->GetLoader("TRDLoader");
261   if (!loader->TreeD()) {
262     loader->LoadDigits();
263   }
264   if (!loader->TreeD()) {
265      return kFALSE;
266   }
267
268   return (fDigitsManager->ReadDigits(loader->TreeD()));
269
270 }
271
272 //_____________________________________________________________________________
273 Bool_t AliTRDtrigger::ReadDigits(AliRawReader* rawReader)
274 {
275   //
276   // Reads the digits arrays from the ddl file
277   //
278
279   AliTRDrawData *raw = new AliTRDrawData();
280   fDigitsManager     = raw->Raw2Digits(rawReader);
281
282   return kTRUE;
283
284 }
285
286 //_____________________________________________________________________________
287 Bool_t AliTRDtrigger::ReadTracklets(AliRunLoader *rl) 
288 {
289   //
290   // Reads the tracklets find the tracks
291   //
292
293   Int_t idet;
294
295   AliLoader *loader = rl->GetLoader("TRDLoader");
296   loader->LoadTracks();
297   fTrackletTree     = loader->TreeT();
298
299   TBranch *branch   = fTrackletTree->GetBranch("TRDmcmTracklet");
300   if (!branch) {
301     AliError("Can't get the branch !");
302     return kFALSE;
303   }
304   TObjArray *tracklets = new TObjArray(400);
305   branch->SetAddress(&tracklets);
306
307   Int_t nEntries   = (Int_t) fTrackletTree->GetEntries();
308   Int_t iEntry;
309   Int_t itrk;
310   Int_t iStack;
311   Int_t iStackPrev = -1;
312   
313   for (iEntry = 0; iEntry < nEntries; iEntry++) {    
314
315     fTrackletTree->GetEvent(iEntry);
316     
317     for (itrk = 0; itrk < tracklets->GetEntriesFast(); itrk++) {
318
319       fTrk   = (AliTRDmcmTracklet*)tracklets->UncheckedAt(itrk);
320       idet   = fTrk->GetDetector();
321       iStack = idet / (AliTRDgeometry::Nplan());
322
323       if (iStackPrev != iStack) {
324         if (iStackPrev == -1) {
325           iStackPrev = iStack;
326         } 
327         else {
328           MakeTracks(idet - AliTRDgeometry::Nplan());
329           ResetTracklets();
330           iStackPrev = iStack;
331         }
332       }
333       
334       Tracklets()->Add(fTrk);
335
336       if ((iEntry == (nEntries-1)) && 
337           (itrk   == (tracklets->GetEntriesFast() - 1))) {
338         idet++;
339         MakeTracks(idet-AliTRDgeometry::Nplan());
340         ResetTracklets();
341       }
342
343     }
344
345   }
346
347   loader->UnloadTracks();
348
349   return kTRUE;
350
351 }
352
353 //_____________________________________________________________________________
354 Bool_t AliTRDtrigger::MakeTracklets(Bool_t makeTracks)
355 {
356   //
357   // Create tracklets from digits
358   //
359
360   Int_t chamBeg = 0;
361   Int_t chamEnd = AliTRDgeometry::Ncham();
362   Int_t planBeg = 0;
363   Int_t planEnd = AliTRDgeometry::Nplan();
364   Int_t sectBeg = 0;
365   Int_t sectEnd = AliTRDgeometry::Nsect();
366
367   fTrkTest = new AliTRDmcmTracklet(0,0,0);
368   fMCM     = new AliTRDmcm(fTrigParam,0);
369
370   Int_t   time;
371   Int_t   col;
372   Int_t   row;
373   Int_t   col1;
374   Int_t   col2;
375   Int_t   idet       = -1;
376   Int_t   iStack     = -1;
377   Int_t   iStackPrev = -1;
378   Float_t amp;
379
380   for (Int_t isect = sectBeg; isect < sectEnd; isect++) {
381
382     for (Int_t icham = chamBeg; icham < chamEnd; icham++) {
383
384       // Number of ROBs in the chamber
385       if(icham == 2) {
386         fNROB = 6;
387       } 
388       else {
389         fNROB = 8;
390       }
391
392       for (Int_t iplan = planBeg; iplan < planEnd; iplan++) {
393
394         idet = fGeo->GetDetector(iplan,icham,isect);
395         ResetTracklets();
396         
397         if (makeTracks) {
398           iStack = idet / (AliTRDgeometry::Nplan());
399           if (iStackPrev != iStack) {
400             if (iStackPrev == -1) {
401               iStackPrev = iStack;
402             } 
403             else {
404               MakeTracks(idet-AliTRDgeometry::Nplan());
405               ResetTracklets();
406               iStackPrev = iStack;
407             }
408           }
409         }
410
411         Int_t nRowMax    = fCParam->GetRowMax(iplan,icham,isect);
412         Int_t nColMax    = fCParam->GetColMax(iplan);
413         Int_t nTimeTotal = fCalib->GetNumberOfTimeBins();
414
415         // Get the digits
416         fDigits = fDigitsManager->GetDigits(idet);
417         if (!fDigits) return kFALSE;
418         // This is to take care of switched off super modules
419         if (fDigits->GetNtime() == 0) {
420           continue;
421         }
422         fDigits->Expand();
423         fTrack0 = fDigitsManager->GetDictionary(idet,0);
424         if (!fTrack0) return kFALSE;
425         fTrack0->Expand();
426         fTrack1 = fDigitsManager->GetDictionary(idet,1);
427         if (!fTrack1) return kFALSE;
428         fTrack1->Expand();
429         fTrack2 = fDigitsManager->GetDictionary(idet,2); 
430         if (!fTrack2) return kFALSE;
431         fTrack2->Expand();
432
433         for (Int_t iRob = 0; iRob < fNROB; iRob++) {
434
435           for (Int_t iMcm = 0; iMcm < kNMCM; iMcm++) {
436
437             fMCM->Reset();
438             fMCM->SetRobId(iRob);
439             fMCM->SetChaId(idet);
440
441             SetMCMcoordinates(iMcm);
442
443             row = fMCM->GetRow();
444
445             if ((row < 0) || (row >= nRowMax)) {
446               AliError("MCM row number out of range.");
447               continue;
448             }
449
450             fMCM->GetColRange(col1,col2);
451             
452             for (time = 0; time < nTimeTotal; time++) {
453               for (col = col1; col < col2; col++) {
454                 if ((col >= 0) && (col < nColMax)) {
455                   amp = TMath::Abs(fDigits->GetDataUnchecked(row,col,time));
456                 } 
457                 else {
458                   amp = 0.0;
459                 }
460                 fMCM->SetADC(col-col1,time,amp);
461               }
462             }
463
464             if (fTrigParam->GetTailCancelation()) {
465               fMCM->Filter(fTrigParam->GetNexponential(),fTrigParam->GetFilterType());
466             }
467             
468             if (fMCM->Run()) {
469
470               for (Int_t iSeed = 0; iSeed < kMaxTrackletsPerMCM; iSeed++) {
471                 
472                 if (fMCM->GetSeedCol()[iSeed] < 0) {
473                   continue;
474                 }
475
476                 if (fTrigParam->GetDebugLevel()   > 1) { 
477                   AliInfo(Form("Add tracklet %d in col %02d \n",fNtracklets,fMCM->GetSeedCol()[iSeed]));
478                 }
479
480                 if (TestTracklet(idet,row,iSeed,0)) {
481                   AddTracklet(idet,row,iSeed,fNtracklets++);
482                 }
483
484               }
485
486             }
487
488           }
489
490       
491         }
492
493         // Compress the arrays
494         fDigits->Compress(1,0);
495         fTrack0->Compress(1,0);
496         fTrack1->Compress(1,0);
497         fTrack2->Compress(1,0);
498
499         WriteTracklets(idet);
500
501      }
502     }
503   }
504
505   if (makeTracks) {
506     idet++;
507     MakeTracks(idet - AliTRDgeometry::Nplan());
508     ResetTracklets();
509   }
510
511   return kTRUE;
512
513 }
514
515 //_____________________________________________________________________________
516 void AliTRDtrigger::SetMCMcoordinates(Int_t imcm)
517 {
518   //
519   // Configure MCM position in the pad plane
520   //
521
522   Int_t robid = fMCM->GetRobId();
523
524   // setting the Row and Col range
525
526   const Int_t kNcolRob = 2;  // number of ROBs per chamber in column direction
527   const Int_t kNmcmRob = 4;  // number of MCMs per ROB in column/row direction
528
529   Int_t mcmid = imcm%(kNmcmRob*kNmcmRob);
530
531   if (robid%kNcolRob == 0) {
532
533     if (mcmid%kNmcmRob == 0) {
534       fMCM->SetColRange(18*0-1,18*1-1+2+1);
535     }
536     if (mcmid%kNmcmRob == 1) {
537       fMCM->SetColRange(18*1-1,18*2-1+2+1);
538     }
539     if (mcmid%kNmcmRob == 2) {
540       fMCM->SetColRange(18*2-1,18*3-1+2+1);
541     }
542     if (mcmid%kNmcmRob == 3) {
543       fMCM->SetColRange(18*3-1,18*4-1+2+1);
544     }
545
546   } 
547   else {
548
549     if (mcmid%kNmcmRob == 0) {
550       fMCM->SetColRange(18*4-1,18*5-1+2+1);
551     }
552     if (mcmid%kNmcmRob == 1) {
553       fMCM->SetColRange(18*5-1,18*6-1+2+1);
554     }
555     if (mcmid%kNmcmRob == 2) {
556       fMCM->SetColRange(18*6-1,18*7-1+2+1);
557     }
558     if (mcmid%kNmcmRob == 3) {
559       fMCM->SetColRange(18*7-1,18*8-1+2+1);
560     }
561
562   } 
563
564   fMCM->SetRow(kNmcmRob*(robid/kNcolRob)+mcmid/kNmcmRob);
565
566 }
567
568 //_____________________________________________________________________________
569 Bool_t AliTRDtrigger::TestTracklet(Int_t det, Int_t row, Int_t seed, Int_t n)
570 {
571   //
572   // Check first the tracklet pt
573   //
574
575   Int_t nTimeTotal  = fCalib->GetNumberOfTimeBins();
576
577   // Calibration fill 2D
578   AliTRDCalibra *calibra = AliTRDCalibra::Instance();
579   if (!calibra) {
580     AliInfo("Could not get Calibra instance\n");
581   }
582
583   fTrkTest->Reset();
584
585   fTrkTest->SetDetector(det);
586   fTrkTest->SetRow(row);
587   fTrkTest->SetN(n);
588
589   Int_t iCol, iCol1, iCol2, track[3];
590   iCol = fMCM->GetSeedCol()[seed];  // 0....20 (MCM)
591   fMCM->GetColRange(iCol1,iCol2);   // range in the pad plane
592             
593   Float_t amp[3];
594   for (Int_t iTime = 0; iTime < nTimeTotal; iTime++) {
595
596     amp[0] = fMCM->GetADC(iCol-1,iTime);
597     amp[1] = fMCM->GetADC(iCol  ,iTime);
598     amp[2] = fMCM->GetADC(iCol+1,iTime);
599
600     // extract track contribution only from the central pad
601     track[0] = fTrack0->GetDataUnchecked(row,iCol+iCol1,iTime);
602     track[1] = fTrack1->GetDataUnchecked(row,iCol+iCol1,iTime);
603     track[2] = fTrack2->GetDataUnchecked(row,iCol+iCol1,iTime);
604
605     if      (fMCM->IsCluster(iCol,iTime)) {
606
607       fTrkTest->AddCluster(iCol+iCol1,iTime,amp,track);
608
609     } 
610     else if ((iCol+1+1) < kMcmCol) {
611
612       amp[0] = fMCM->GetADC(iCol-1+1,iTime);
613       amp[1] = fMCM->GetADC(iCol  +1,iTime);
614       amp[2] = fMCM->GetADC(iCol+1+1,iTime);
615
616       if (fMCM->IsCluster(iCol+1,iTime)) {
617
618         // extract track contribution only from the central pad
619         track[0] = fTrack0->GetDataUnchecked(row,iCol+1+iCol1,iTime);
620         track[1] = fTrack1->GetDataUnchecked(row,iCol+1+iCol1,iTime);
621         track[2] = fTrack2->GetDataUnchecked(row,iCol+1+iCol1,iTime);
622
623         fTrkTest->AddCluster(iCol+1+iCol1,iTime,amp,track);
624
625       }
626
627     } 
628
629   }
630
631   fTrkTest->CookLabel(0.8);  
632   /*
633   if (fTrkTest->GetLabel() >= fNPrimary) {
634     Info("AddTracklet","Only primaries are stored!");
635     return;
636   }
637   */
638   // LTU Pt cut
639   fTrkTest->MakeTrackletGraph(fGeo,fField);
640
641   // TRD Online calibration
642   if (calibra->GetMcmTracking()) {
643     calibra->UpdateHistogramcm(fTrkTest);
644   }
645
646   fTrkTest->MakeClusAmpGraph();
647
648   if (TMath::Abs(fTrkTest->GetPt()) < fTrigParam->GetLtuPtCut()) {
649     return kFALSE;
650   }
651   
652   return kTRUE;  
653
654 }
655
656 //_____________________________________________________________________________
657 void AliTRDtrigger::AddTracklet(Int_t det, Int_t row, Int_t seed, Int_t n)
658 {
659   //
660   // Add a found tracklet
661   //
662
663   Int_t nTimeTotal  = fCalib->GetNumberOfTimeBins();
664
665   fTrk = new AliTRDmcmTracklet(det,row,n);
666
667   Int_t iCol, iCol1, iCol2, track[3];
668   iCol = fMCM->GetSeedCol()[seed];  // 0....20 (MCM)
669   fMCM->GetColRange(iCol1,iCol2);   // range in the pad plane
670             
671   Float_t amp[3];
672   for (Int_t iTime = 0; iTime < nTimeTotal; iTime++) {
673
674     amp[0] = fMCM->GetADC(iCol-1,iTime);
675     amp[1] = fMCM->GetADC(iCol  ,iTime);
676     amp[2] = fMCM->GetADC(iCol+1,iTime);
677
678     // extract track contribution only from the central pad
679     track[0] = fTrack0->GetDataUnchecked(row,iCol+iCol1,iTime);
680     track[1] = fTrack1->GetDataUnchecked(row,iCol+iCol1,iTime);
681     track[2] = fTrack2->GetDataUnchecked(row,iCol+iCol1,iTime);
682
683     if      (fMCM->IsCluster(iCol,iTime)) {
684
685       fTrk->AddCluster(iCol+iCol1,iTime,amp,track);
686
687     } 
688     else if ((iCol+1+1) < kMcmCol) {
689
690       amp[0] = fMCM->GetADC(iCol-1+1,iTime);
691       amp[1] = fMCM->GetADC(iCol  +1,iTime);
692       amp[2] = fMCM->GetADC(iCol+1+1,iTime);
693
694       if (fMCM->IsCluster(iCol+1,iTime)) {
695
696         // extract track contribution only from the central pad
697         track[0] = fTrack0->GetDataUnchecked(row,iCol+1+iCol1,iTime);
698         track[1] = fTrack1->GetDataUnchecked(row,iCol+1+iCol1,iTime);
699         track[2] = fTrack2->GetDataUnchecked(row,iCol+1+iCol1,iTime);
700
701         fTrk->AddCluster(iCol+1+iCol1,iTime,amp,track);
702
703       }
704
705     }
706
707   }
708
709   fTrk->CookLabel(0.8);  
710   /*
711   if (fTrk->GetLabel() >= fNPrimary) {
712     Info("AddTracklet","Only primaries are stored!");
713     return;
714   }
715   */
716   // LTU Pt cut
717   fTrk->MakeTrackletGraph(fGeo,fField);
718   fTrk->MakeClusAmpGraph();
719   if (TMath::Abs(fTrk->GetPt()) < fTrigParam->GetLtuPtCut()) {
720     return;
721   }
722       
723   Tracklets()->Add(fTrk);
724
725 }
726
727 //_____________________________________________________________________________
728 Bool_t AliTRDtrigger::WriteTracklets(Int_t det) 
729 {
730   //
731   // Fills TRDmcmTracklet branch in the tree with the Tracklets 
732   // found in detector = det. For det=-1 writes the tree. 
733   //
734
735   if ((det < -1) || (det >= AliTRDgeometry::Ndet())) {
736     AliError(Form("Unexpected detector index %d.",det));
737     return kFALSE;
738   }
739
740   TBranch *branch = fTrackletTree->GetBranch("TRDmcmTracklet");
741   if (!branch) {
742     TObjArray *ioArray = 0;
743     branch = fTrackletTree->Branch("TRDmcmTracklet","TObjArray",&ioArray,32000,0);
744   }
745
746   if ((det >= 0) && (det < AliTRDgeometry::Ndet())) {
747
748     Int_t nTracklets = Tracklets()->GetEntriesFast();
749     TObjArray *detTracklets = new TObjArray(400);
750
751     for (Int_t i = 0; i < nTracklets; i++) {
752
753       AliTRDmcmTracklet *trk = (AliTRDmcmTracklet *) Tracklets()->UncheckedAt(i);
754       
755       if (det == trk->GetDetector()) {
756         detTracklets->AddLast(trk);
757       }
758
759     }
760
761     branch->SetAddress(&detTracklets);
762     fTrackletTree->Fill();
763
764     delete detTracklets;
765
766     return kTRUE;
767
768   }
769
770   if (det == -1) {
771
772     AliInfo(Form("Writing the Tracklet tree %s for event %d."
773                 ,fTrackletTree->GetName(),fRunLoader->GetEventNumber()));
774
775     AliLoader* loader = fRunLoader->GetLoader("TRDLoader");
776     loader->WriteTracks("OVERWRITE");
777     
778     return kTRUE;  
779
780   }
781
782   return kFALSE;
783
784 }
785
786 //_____________________________________________________________________________
787 void AliTRDtrigger::MakeTracks(Int_t det)
788 {
789   //
790   // Create GTU tracks per module (stack of 6 chambers)
791   //
792   
793   fModule->Reset();
794
795   Int_t nRowMax, iplan, icham, isect, row;
796
797   if ((det < 0) || (det >= AliTRDgeometry::Ndet())) {
798     AliError(Form("Unexpected detector index %d.",det));
799     return;
800   }
801   
802   Int_t nTracklets = Tracklets()->GetEntriesFast();
803   
804   AliTRDmcmTracklet *trk;
805   for (Int_t i = 0; i < nTracklets; i++) {
806     
807     trk = (AliTRDmcmTracklet *) Tracklets()->UncheckedAt(i);
808     
809     iplan = fGeo->GetPlane(trk->GetDetector());
810     icham = fGeo->GetChamber(trk->GetDetector());
811     isect = fGeo->GetSector(trk->GetDetector());
812
813     nRowMax = fCParam->GetRowMax(iplan,icham,isect);
814     row = trk->GetRow();
815
816     fModule->AddTracklet(trk->GetDetector(),
817                          row,
818                          trk->GetRowz(),
819                          trk->GetSlope(),
820                          trk->GetOffset(),
821                          trk->GetTime0(),
822                          trk->GetNclusters(),
823                          trk->GetLabel(),
824                          trk->GetdQdl());
825     
826   }
827
828   fModule->SortTracklets();
829   fModule->RemoveMultipleTracklets();
830   fModule->SortZ((Int_t)fGeo->GetChamber(det));
831   fModule->FindTracks();
832   fModule->SortTracks();
833   fModule->RemoveMultipleTracks();
834
835   Int_t nModTracks = fModule->GetNtracks();
836   AliTRDgtuTrack *gtutrk;
837   for (Int_t i = 0; i < nModTracks; i++) {
838     gtutrk = (AliTRDgtuTrack*)fModule->GetTrack(i);
839     if (TMath::Abs(gtutrk->GetPt()) < fTrigParam->GetGtuPtCut()) continue;
840     gtutrk->CookLabel();
841     gtutrk->MakePID();
842     AddTrack(gtutrk,det);
843   }
844   
845 }
846
847 //_____________________________________________________________________________
848 void AliTRDtrigger::AddTrack(const AliTRDgtuTrack *t, Int_t det)  
849
850   //
851   // Add a track to the list
852   //
853
854   AliTRDgtuTrack *track = new(fTracks->operator[](fTracks->GetEntriesFast())) 
855                           AliTRDgtuTrack(*t);
856   track->SetDetector(det); 
857
858 }
859
860 //_____________________________________________________________________________
861 TObjArray* AliTRDtrigger::Tracklets()
862
863   //
864   // Returns list of tracklets
865   //
866
867   if (!fTracklets) {
868     fTracklets = new TObjArray(400); 
869   }
870   return fTracklets;       
871
872 }
873
874 //_____________________________________________________________________________
875 void AliTRDtrigger::ResetTracklets()                              
876 {
877   //
878   // Resets the list of tracklets
879   //
880  
881   if (fTracklets) {
882     fTracklets->Delete();             
883   }
884
885 }
886
887 //_____________________________________________________________________________
888 Int_t AliTRDtrigger::GetNumberOfTracks() const                     
889
890   //
891   // Returns number of tracks
892   //
893
894   return fTracks->GetEntriesFast(); 
895                 
896 }
897
898 //_____________________________________________________________________________
899 AliTRDgtuTrack* AliTRDtrigger::GetTrack(Int_t i) const               
900
901   //
902   // Returns a given track from the list
903   //
904
905   return (AliTRDgtuTrack *) fTracks->UncheckedAt(i); 
906
907 }