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