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