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 ///////////////////////////////////////////////////////////////////////////////
18 // TRD trigger class //
23 ///////////////////////////////////////////////////////////////////////////////
28 #include <TClonesArray.h>
29 #include <TObjArray.h>
33 #include "AliLoader.h"
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"
51 ClassImp(AliTRDtrigger)
53 //_____________________________________________________________________________
54 AliTRDtrigger::AliTRDtrigger()
80 // AliTRDtrigger default constructor
85 //_____________________________________________________________________________
86 AliTRDtrigger::AliTRDtrigger(const Text_t *name, const Text_t *title)
94 ,fDigitsManager(new AliTRDdigitsManager())
96 ,fTracklets(new TObjArray(400))
109 ,fTracks(new TClonesArray("AliTRDgtuTrack",1000))
112 // AliTRDtrigger constructor
117 //_____________________________________________________________________________
118 AliTRDtrigger::AliTRDtrigger(const AliTRDtrigger &p)
126 ,fDigitsManager(NULL)
135 ,fNtracklets(p.fNtracklets)
140 ,fNPrimary(p.fNPrimary)
144 // AliTRDtrigger copy constructor
149 ///_____________________________________________________________________________
150 AliTRDtrigger::~AliTRDtrigger()
153 // AliTRDtrigger destructor
157 fTracklets->Delete();
168 //_____________________________________________________________________________
169 AliTRDtrigger &AliTRDtrigger::operator=(const AliTRDtrigger &p)
172 // Assignment operator
175 if (this != &p) ((AliTRDtrigger &) p).Copy(*this);
180 //_____________________________________________________________________________
181 void AliTRDtrigger::Copy(TObject &) const
187 AliFatal("Not implemented");
191 //_____________________________________________________________________________
192 void AliTRDtrigger::Init()
195 fModule = new AliTRDmodule(fTrigParam);
198 fField = fTrigParam->GetField();
199 fGeo = (AliTRDgeometry*)AliTRDgeometry::GetGeometry(fRunLoader);
201 fCalib = AliTRDcalibDB::Instance();
203 AliError("No instance of AliTRDcalibDB.");
207 fCParam = AliTRDCommonParam::Instance();
209 AliError("No common parameters.");
215 //_____________________________________________________________________________
216 Bool_t AliTRDtrigger::Open(const Char_t *name, Int_t nEvent)
219 // Opens the AliROOT file.
222 TString evfoldname = AliConfig::GetDefaultEventFolderName();
223 fRunLoader = AliRunLoader::GetRunLoader(evfoldname);
226 fRunLoader = AliRunLoader::Open(name);
229 AliError(Form("Can not open session for file %s.",name));
233 // Import the Trees for the event nEvent in the file
234 fRunLoader->GetEvent(nEvent);
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);
248 //_____________________________________________________________________________
249 Bool_t AliTRDtrigger::ReadDigits()
252 // Reads the digits arrays from the input aliroot file
256 AliError("Can not find the Run Loader");
260 AliLoader* loader = fRunLoader->GetLoader("TRDLoader");
261 if (!loader->TreeD()) {
262 loader->LoadDigits();
264 if (!loader->TreeD()) {
268 return (fDigitsManager->ReadDigits(loader->TreeD()));
272 //_____________________________________________________________________________
273 Bool_t AliTRDtrigger::ReadDigits(AliRawReader* rawReader)
276 // Reads the digits arrays from the ddl file
279 AliTRDrawData *raw = new AliTRDrawData();
280 fDigitsManager = raw->Raw2Digits(rawReader);
286 //_____________________________________________________________________________
287 Bool_t AliTRDtrigger::ReadTracklets(AliRunLoader *rl)
290 // Reads the tracklets find the tracks
295 AliLoader *loader = rl->GetLoader("TRDLoader");
296 loader->LoadTracks();
297 fTrackletTree = loader->TreeT();
299 TBranch *branch = fTrackletTree->GetBranch("TRDmcmTracklet");
301 AliError("Can't get the branch !");
304 TObjArray *tracklets = new TObjArray(400);
305 branch->SetAddress(&tracklets);
307 Int_t nEntries = (Int_t) fTrackletTree->GetEntries();
311 Int_t iStackPrev = -1;
313 for (iEntry = 0; iEntry < nEntries; iEntry++) {
315 fTrackletTree->GetEvent(iEntry);
317 for (itrk = 0; itrk < tracklets->GetEntriesFast(); itrk++) {
319 fTrk = (AliTRDmcmTracklet*)tracklets->UncheckedAt(itrk);
320 idet = fTrk->GetDetector();
321 iStack = idet / (AliTRDgeometry::Nplan());
323 if (iStackPrev != iStack) {
324 if (iStackPrev == -1) {
328 MakeTracks(idet - AliTRDgeometry::Nplan());
334 Tracklets()->Add(fTrk);
336 if ((iEntry == (nEntries-1)) &&
337 (itrk == (tracklets->GetEntriesFast() - 1))) {
339 MakeTracks(idet-AliTRDgeometry::Nplan());
347 loader->UnloadTracks();
353 //_____________________________________________________________________________
354 Bool_t AliTRDtrigger::MakeTracklets(Bool_t makeTracks)
357 // Create tracklets from digits
361 Int_t chamEnd = AliTRDgeometry::Ncham();
363 Int_t planEnd = AliTRDgeometry::Nplan();
365 Int_t sectEnd = AliTRDgeometry::Nsect();
367 fTrkTest = new AliTRDmcmTracklet(0,0,0);
368 fMCM = new AliTRDmcm(fTrigParam,0);
377 Int_t iStackPrev = -1;
380 for (Int_t isect = sectBeg; isect < sectEnd; isect++) {
382 for (Int_t icham = chamBeg; icham < chamEnd; icham++) {
384 // Number of ROBs in the chamber
392 for (Int_t iplan = planBeg; iplan < planEnd; iplan++) {
394 idet = fGeo->GetDetector(iplan,icham,isect);
398 iStack = idet / (AliTRDgeometry::Nplan());
399 if (iStackPrev != iStack) {
400 if (iStackPrev == -1) {
404 MakeTracks(idet-AliTRDgeometry::Nplan());
411 Int_t nRowMax = fCParam->GetRowMax(iplan,icham,isect);
412 Int_t nColMax = fCParam->GetColMax(iplan);
413 Int_t nTimeTotal = fCalib->GetNumberOfTimeBins();
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) {
423 fTrack0 = fDigitsManager->GetDictionary(idet,0);
424 if (!fTrack0) return kFALSE;
426 fTrack1 = fDigitsManager->GetDictionary(idet,1);
427 if (!fTrack1) return kFALSE;
429 fTrack2 = fDigitsManager->GetDictionary(idet,2);
430 if (!fTrack2) return kFALSE;
433 for (Int_t iRob = 0; iRob < fNROB; iRob++) {
435 for (Int_t iMcm = 0; iMcm < kNMCM; iMcm++) {
438 fMCM->SetRobId(iRob);
439 fMCM->SetChaId(idet);
441 SetMCMcoordinates(iMcm);
443 row = fMCM->GetRow();
445 if ((row < 0) || (row >= nRowMax)) {
446 AliError("MCM row number out of range.");
450 fMCM->GetColRange(col1,col2);
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));
460 fMCM->SetADC(col-col1,time,amp);
464 if (fTrigParam->GetTailCancelation()) {
465 fMCM->Filter(fTrigParam->GetNexponential(),fTrigParam->GetFilterType());
470 for (Int_t iSeed = 0; iSeed < kMaxTrackletsPerMCM; iSeed++) {
472 if (fMCM->GetSeedCol()[iSeed] < 0) {
476 if (fTrigParam->GetDebugLevel() > 1) {
477 AliInfo(Form("Add tracklet %d in col %02d \n",fNtracklets,fMCM->GetSeedCol()[iSeed]));
480 if (TestTracklet(idet,row,iSeed,0)) {
481 AddTracklet(idet,row,iSeed,fNtracklets++);
493 // Compress the arrays
494 fDigits->Compress(1,0);
495 fTrack0->Compress(1,0);
496 fTrack1->Compress(1,0);
497 fTrack2->Compress(1,0);
499 WriteTracklets(idet);
507 MakeTracks(idet - AliTRDgeometry::Nplan());
515 //_____________________________________________________________________________
516 void AliTRDtrigger::SetMCMcoordinates(Int_t imcm)
519 // Configure MCM position in the pad plane
522 Int_t robid = fMCM->GetRobId();
524 // setting the Row and Col range
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
529 Int_t mcmid = imcm%(kNmcmRob*kNmcmRob);
531 if (robid%kNcolRob == 0) {
533 if (mcmid%kNmcmRob == 0) {
534 fMCM->SetColRange(18*0-1,18*1-1+2+1);
536 if (mcmid%kNmcmRob == 1) {
537 fMCM->SetColRange(18*1-1,18*2-1+2+1);
539 if (mcmid%kNmcmRob == 2) {
540 fMCM->SetColRange(18*2-1,18*3-1+2+1);
542 if (mcmid%kNmcmRob == 3) {
543 fMCM->SetColRange(18*3-1,18*4-1+2+1);
549 if (mcmid%kNmcmRob == 0) {
550 fMCM->SetColRange(18*4-1,18*5-1+2+1);
552 if (mcmid%kNmcmRob == 1) {
553 fMCM->SetColRange(18*5-1,18*6-1+2+1);
555 if (mcmid%kNmcmRob == 2) {
556 fMCM->SetColRange(18*6-1,18*7-1+2+1);
558 if (mcmid%kNmcmRob == 3) {
559 fMCM->SetColRange(18*7-1,18*8-1+2+1);
564 fMCM->SetRow(kNmcmRob*(robid/kNcolRob)+mcmid/kNmcmRob);
568 //_____________________________________________________________________________
569 Bool_t AliTRDtrigger::TestTracklet(Int_t det, Int_t row, Int_t seed, Int_t n)
572 // Check first the tracklet pt
575 Int_t nTimeTotal = fCalib->GetNumberOfTimeBins();
577 // Calibration fill 2D
578 AliTRDCalibra *calibra = AliTRDCalibra::Instance();
580 AliInfo("Could not get Calibra instance\n");
585 fTrkTest->SetDetector(det);
586 fTrkTest->SetRow(row);
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
594 for (Int_t iTime = 0; iTime < nTimeTotal; iTime++) {
596 amp[0] = fMCM->GetADC(iCol-1,iTime);
597 amp[1] = fMCM->GetADC(iCol ,iTime);
598 amp[2] = fMCM->GetADC(iCol+1,iTime);
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);
605 if (fMCM->IsCluster(iCol,iTime)) {
607 fTrkTest->AddCluster(iCol+iCol1,iTime,amp,track);
610 else if ((iCol+1+1) < kMcmCol) {
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);
616 if (fMCM->IsCluster(iCol+1,iTime)) {
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);
623 fTrkTest->AddCluster(iCol+1+iCol1,iTime,amp,track);
631 fTrkTest->CookLabel(0.8);
633 if (fTrkTest->GetLabel() >= fNPrimary) {
634 Info("AddTracklet","Only primaries are stored!");
639 fTrkTest->MakeTrackletGraph(fGeo,fField);
641 // TRD Online calibration
642 if (calibra->GetMcmTracking()) {
643 calibra->UpdateHistogramcm(fTrkTest);
646 fTrkTest->MakeClusAmpGraph();
648 if (TMath::Abs(fTrkTest->GetPt()) < fTrigParam->GetLtuPtCut()) {
656 //_____________________________________________________________________________
657 void AliTRDtrigger::AddTracklet(Int_t det, Int_t row, Int_t seed, Int_t n)
660 // Add a found tracklet
663 Int_t nTimeTotal = fCalib->GetNumberOfTimeBins();
665 fTrk = new AliTRDmcmTracklet(det,row,n);
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
672 for (Int_t iTime = 0; iTime < nTimeTotal; iTime++) {
674 amp[0] = fMCM->GetADC(iCol-1,iTime);
675 amp[1] = fMCM->GetADC(iCol ,iTime);
676 amp[2] = fMCM->GetADC(iCol+1,iTime);
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);
683 if (fMCM->IsCluster(iCol,iTime)) {
685 fTrk->AddCluster(iCol+iCol1,iTime,amp,track);
688 else if ((iCol+1+1) < kMcmCol) {
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);
694 if (fMCM->IsCluster(iCol+1,iTime)) {
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);
701 fTrk->AddCluster(iCol+1+iCol1,iTime,amp,track);
709 fTrk->CookLabel(0.8);
711 if (fTrk->GetLabel() >= fNPrimary) {
712 Info("AddTracklet","Only primaries are stored!");
717 fTrk->MakeTrackletGraph(fGeo,fField);
718 fTrk->MakeClusAmpGraph();
719 if (TMath::Abs(fTrk->GetPt()) < fTrigParam->GetLtuPtCut()) {
723 Tracklets()->Add(fTrk);
727 //_____________________________________________________________________________
728 Bool_t AliTRDtrigger::WriteTracklets(Int_t det)
731 // Fills TRDmcmTracklet branch in the tree with the Tracklets
732 // found in detector = det. For det=-1 writes the tree.
735 if ((det < -1) || (det >= AliTRDgeometry::Ndet())) {
736 AliError(Form("Unexpected detector index %d.",det));
740 TBranch *branch = fTrackletTree->GetBranch("TRDmcmTracklet");
742 TObjArray *ioArray = 0;
743 branch = fTrackletTree->Branch("TRDmcmTracklet","TObjArray",&ioArray,32000,0);
746 if ((det >= 0) && (det < AliTRDgeometry::Ndet())) {
748 Int_t nTracklets = Tracklets()->GetEntriesFast();
749 TObjArray *detTracklets = new TObjArray(400);
751 for (Int_t i = 0; i < nTracklets; i++) {
753 AliTRDmcmTracklet *trk = (AliTRDmcmTracklet *) Tracklets()->UncheckedAt(i);
755 if (det == trk->GetDetector()) {
756 detTracklets->AddLast(trk);
761 branch->SetAddress(&detTracklets);
762 fTrackletTree->Fill();
772 AliInfo(Form("Writing the Tracklet tree %s for event %d."
773 ,fTrackletTree->GetName(),fRunLoader->GetEventNumber()));
775 AliLoader* loader = fRunLoader->GetLoader("TRDLoader");
776 loader->WriteTracks("OVERWRITE");
786 //_____________________________________________________________________________
787 void AliTRDtrigger::MakeTracks(Int_t det)
790 // Create GTU tracks per module (stack of 6 chambers)
795 Int_t nRowMax, iplan, icham, isect, row;
797 if ((det < 0) || (det >= AliTRDgeometry::Ndet())) {
798 AliError(Form("Unexpected detector index %d.",det));
802 Int_t nTracklets = Tracklets()->GetEntriesFast();
804 AliTRDmcmTracklet *trk;
805 for (Int_t i = 0; i < nTracklets; i++) {
807 trk = (AliTRDmcmTracklet *) Tracklets()->UncheckedAt(i);
809 iplan = fGeo->GetPlane(trk->GetDetector());
810 icham = fGeo->GetChamber(trk->GetDetector());
811 isect = fGeo->GetSector(trk->GetDetector());
813 nRowMax = fCParam->GetRowMax(iplan,icham,isect);
816 fModule->AddTracklet(trk->GetDetector(),
828 fModule->SortTracklets();
829 fModule->RemoveMultipleTracklets();
830 fModule->SortZ((Int_t)fGeo->GetChamber(det));
831 fModule->FindTracks();
832 fModule->SortTracks();
833 fModule->RemoveMultipleTracks();
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;
842 AddTrack(gtutrk,det);
847 //_____________________________________________________________________________
848 void AliTRDtrigger::AddTrack(const AliTRDgtuTrack *t, Int_t det)
851 // Add a track to the list
854 AliTRDgtuTrack *track = new(fTracks->operator[](fTracks->GetEntriesFast()))
856 track->SetDetector(det);
860 //_____________________________________________________________________________
861 TObjArray* AliTRDtrigger::Tracklets()
864 // Returns list of tracklets
868 fTracklets = new TObjArray(400);
874 //_____________________________________________________________________________
875 void AliTRDtrigger::ResetTracklets()
878 // Resets the list of tracklets
882 fTracklets->Delete();
887 //_____________________________________________________________________________
888 Int_t AliTRDtrigger::GetNumberOfTracks() const
891 // Returns number of tracks
894 return fTracks->GetEntriesFast();
898 //_____________________________________________________________________________
899 AliTRDgtuTrack* AliTRDtrigger::GetTrack(Int_t i) const
902 // Returns a given track from the list
905 return (AliTRDgtuTrack *) fTracks->UncheckedAt(i);