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 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
20 // TRD cluster finder //
22 ///////////////////////////////////////////////////////////////////////////////
24 #include <TClonesArray.h>
25 #include <TObjArray.h>
27 #include "AliRunLoader.h"
28 #include "AliLoader.h"
29 #include "AliTreeLoader.h"
30 #include "AliAlignObj.h"
32 #include "AliTRDclusterizer.h"
33 #include "AliTRDcluster.h"
34 #include "AliTRDReconstructor.h"
35 #include "AliTRDgeometry.h"
36 #include "AliTRDarrayDictionary.h"
37 #include "AliTRDarrayADC.h"
38 #include "AliTRDdigitsManager.h"
39 #include "AliTRDdigitsParam.h"
40 #include "AliTRDrawData.h"
41 #include "AliTRDcalibDB.h"
42 #include "AliTRDtransform.h"
43 #include "AliTRDSignalIndex.h"
44 #include "AliTRDrawStreamBase.h"
45 #include "AliTRDrawStream.h"
46 #include "AliTRDfeeParam.h"
47 #include "AliTRDtrackletWord.h"
49 #include "TTreeStream.h"
51 #include "Cal/AliTRDCalROC.h"
52 #include "Cal/AliTRDCalDet.h"
53 #include "Cal/AliTRDCalSingleChamberStatus.h"
55 ClassImp(AliTRDclusterizer)
57 //_____________________________________________________________________________
58 AliTRDclusterizer::AliTRDclusterizer(const AliTRDReconstructor *const rec)
66 ,fDigitsManager(new AliTRDdigitsManager())
67 ,fTrackletContainer(NULL)
69 ,fTransform(new AliTRDtransform(0))
76 ,fMinLeftRightCutSigma(0)
82 ,fCalGainFactorROC(NULL)
83 ,fCalGainFactorDetValue(0)
86 ,fCalPadStatusROC(NULL)
94 // AliTRDclusterizer default constructor
97 SetBit(kLabels, kTRUE);
98 SetBit(knewDM, kFALSE);
100 fRawVersion = AliTRDfeeParam::Instance()->GetRAWversion();
102 // Initialize debug stream
104 if(fReconstructor->GetRecoParam()->GetStreamLevel(AliTRDrecoParam::kClusterizer) > 1){
105 TDirectory *savedir = gDirectory;
106 //fgGetDebugStream = new TTreeSRedirector("TRD.ClusterizerDebug.root");
113 //_____________________________________________________________________________
114 AliTRDclusterizer::AliTRDclusterizer(const Text_t *name, const Text_t *title, const AliTRDReconstructor *const rec)
122 ,fDigitsManager(new AliTRDdigitsManager())
123 ,fTrackletContainer(NULL)
125 ,fTransform(new AliTRDtransform(0))
132 ,fMinLeftRightCutSigma(0)
138 ,fCalGainFactorROC(NULL)
139 ,fCalGainFactorDetValue(0)
141 ,fCalNoiseDetValue(0)
142 ,fCalPadStatusROC(NULL)
150 // AliTRDclusterizer constructor
153 SetBit(kLabels, kTRUE);
154 SetBit(knewDM, kFALSE);
156 AliTRDcalibDB *trd = 0x0;
157 if (!(trd = AliTRDcalibDB::Instance())) {
158 AliFatal("Could not get calibration object");
161 fDigitsManager->CreateArrays();
163 fRawVersion = AliTRDfeeParam::Instance()->GetRAWversion();
169 //_____________________________________________________________________________
170 AliTRDclusterizer::AliTRDclusterizer(const AliTRDclusterizer &c)
172 ,fReconstructor(c.fReconstructor)
178 ,fDigitsManager(NULL)
179 ,fTrackletContainer(NULL)
188 ,fMinLeftRightCutSigma(0)
194 ,fCalGainFactorROC(NULL)
195 ,fCalGainFactorDetValue(0)
197 ,fCalNoiseDetValue(0)
198 ,fCalPadStatusROC(NULL)
206 // AliTRDclusterizer copy constructor
209 SetBit(kLabels, kTRUE);
210 SetBit(knewDM, kFALSE);
216 //_____________________________________________________________________________
217 AliTRDclusterizer::~AliTRDclusterizer()
220 // AliTRDclusterizer destructor
223 if (fRecPoints/* && IsClustersOwner()*/){
224 fRecPoints->Delete();
229 fTracklets->Delete();
233 if (fDigitsManager) {
234 delete fDigitsManager;
235 fDigitsManager = NULL;
238 if (fTrackletContainer){
239 delete [] fTrackletContainer[0];
240 delete [] fTrackletContainer[1];
241 delete [] fTrackletContainer;
242 fTrackletContainer = NULL;
257 //_____________________________________________________________________________
258 AliTRDclusterizer &AliTRDclusterizer::operator=(const AliTRDclusterizer &c)
261 // Assignment operator
266 ((AliTRDclusterizer &) c).Copy(*this);
273 //_____________________________________________________________________________
274 void AliTRDclusterizer::Copy(TObject &c) const
280 ((AliTRDclusterizer &) c).fClusterTree = NULL;
281 ((AliTRDclusterizer &) c).fRecPoints = NULL;
282 ((AliTRDclusterizer &) c).fTrackletTree = NULL;
283 ((AliTRDclusterizer &) c).fDigitsManager = NULL;
284 ((AliTRDclusterizer &) c).fTrackletContainer = NULL;
285 ((AliTRDclusterizer &) c).fRawVersion = fRawVersion;
286 ((AliTRDclusterizer &) c).fTransform = NULL;
287 ((AliTRDclusterizer &) c).fDigits = NULL;
288 ((AliTRDclusterizer &) c).fIndexes = NULL;
289 ((AliTRDclusterizer &) c).fMaxThresh = 0;
290 ((AliTRDclusterizer &) c).fMaxThreshTest = 0;
291 ((AliTRDclusterizer &) c).fSigThresh = 0;
292 ((AliTRDclusterizer &) c).fMinMaxCutSigma= 0;
293 ((AliTRDclusterizer &) c).fMinLeftRightCutSigma = 0;
294 ((AliTRDclusterizer &) c).fLayer = 0;
295 ((AliTRDclusterizer &) c).fDet = 0;
296 ((AliTRDclusterizer &) c).fVolid = 0;
297 ((AliTRDclusterizer &) c).fColMax = 0;
298 ((AliTRDclusterizer &) c).fTimeTotal = 0;
299 ((AliTRDclusterizer &) c).fCalGainFactorROC = NULL;
300 ((AliTRDclusterizer &) c).fCalGainFactorDetValue = 0;
301 ((AliTRDclusterizer &) c).fCalNoiseROC = NULL;
302 ((AliTRDclusterizer &) c).fCalNoiseDetValue = 0;
303 ((AliTRDclusterizer &) c).fCalPadStatusROC = NULL;
304 ((AliTRDclusterizer &) c).fClusterROC = 0;
305 ((AliTRDclusterizer &) c).firstClusterROC= 0;
306 ((AliTRDclusterizer &) c).fNoOfClusters = 0;
307 ((AliTRDclusterizer &) c).fBaseline = 0;
308 ((AliTRDclusterizer &) c).fRawStream = NULL;
312 //_____________________________________________________________________________
313 Bool_t AliTRDclusterizer::Open(const Char_t *name, Int_t nEvent)
316 // Opens the AliROOT file. Output and input are in the same file
319 TString evfoldname = AliConfig::GetDefaultEventFolderName();
320 fRunLoader = AliRunLoader::GetRunLoader(evfoldname);
323 fRunLoader = AliRunLoader::Open(name);
327 AliError(Form("Can not open session for file %s.",name));
338 //_____________________________________________________________________________
339 Bool_t AliTRDclusterizer::OpenOutput()
342 // Open the output file
345 if (!fReconstructor->IsWritingClusters()) return kTRUE;
347 TObjArray *ioArray = 0x0;
349 AliLoader* loader = fRunLoader->GetLoader("TRDLoader");
350 loader->MakeTree("R");
352 fClusterTree = loader->TreeR();
353 fClusterTree->Branch("TRDcluster", "TObjArray", &ioArray, 32000, 0);
359 //_____________________________________________________________________________
360 Bool_t AliTRDclusterizer::OpenOutput(TTree *const clusterTree)
363 // Connect the output tree
367 if (fReconstructor->IsWritingClusters()){
368 TObjArray *ioArray = 0x0;
369 fClusterTree = clusterTree;
370 fClusterTree->Branch("TRDcluster", "TObjArray", &ioArray, 32000, 0);
375 //_____________________________________________________________________________
376 Bool_t AliTRDclusterizer::OpenInput(Int_t nEvent)
379 // Opens a ROOT-file with TRD-hits and reads in the digits-tree
382 // Import the Trees for the event nEvent in the file
383 fRunLoader->GetEvent(nEvent);
389 //_____________________________________________________________________________
390 Bool_t AliTRDclusterizer::WriteClusters(Int_t det)
393 // Fills TRDcluster branch in the tree with the clusters
394 // found in detector = det. For det=-1 writes the tree.
398 (det >= AliTRDgeometry::Ndet())) {
399 AliError(Form("Unexpected detector index %d.\n",det));
403 TObjArray *ioArray = new TObjArray(400);
404 TBranch *branch = fClusterTree->GetBranch("TRDcluster");
406 branch = fClusterTree->Branch("TRDcluster","TObjArray",&ioArray,32000,0);
407 } else branch->SetAddress(&ioArray);
409 Int_t nRecPoints = RecPoints()->GetEntriesFast();
411 for (Int_t i = 0; i < nRecPoints; i++) {
412 AliTRDcluster *c = (AliTRDcluster *) RecPoints()->UncheckedAt(i);
413 if(det != c->GetDetector()) continue;
416 fClusterTree->Fill();
419 Int_t detOld = -1, nw(0);
420 for (Int_t i = 0; i < nRecPoints; i++) {
421 AliTRDcluster *c = (AliTRDcluster *) RecPoints()->UncheckedAt(i);
422 if(c->GetDetector() != detOld){
423 nw += ioArray->GetEntriesFast();
424 fClusterTree->Fill();
426 detOld = c->GetDetector();
430 if(ioArray->GetEntriesFast()){
431 nw += ioArray->GetEntriesFast();
432 fClusterTree->Fill();
435 AliDebug(2, Form("Clusters FOUND[%d] WRITTEN[%d] STATUS[%s]", nRecPoints, nw, nw==nRecPoints?"OK":"FAILED"));
442 //_____________________________________________________________________________
443 Bool_t AliTRDclusterizer::WriteTracklets(Int_t det)
446 // Write the raw data tracklets into seperate file
449 UInt_t **leaves = new UInt_t *[2];
450 for (Int_t i=0; i<2 ;i++){
451 leaves[i] = new UInt_t[258];
452 leaves[i][0] = det; // det
453 leaves[i][1] = i; // side
454 memcpy(leaves[i]+2, fTrackletContainer[i], sizeof(UInt_t) * 256);
458 AliDataLoader *dl = fRunLoader->GetLoader("TRDLoader")->GetDataLoader("tracklets");
460 fTrackletTree = dl->Tree();
463 TBranch *trkbranch = fTrackletTree->GetBranch("trkbranch");
465 trkbranch = fTrackletTree->Branch("trkbranch",leaves[0],"det/i:side/i:tracklets[256]/i");
468 for (Int_t i=0; i<2; i++){
469 if (leaves[i][2]>0) {
470 trkbranch->SetAddress(leaves[i]);
471 fTrackletTree->Fill();
475 // AliDataLoader *dl = fRunLoader->GetLoader("TRDLoader")->GetDataLoader("tracklets"); //jkl: wrong here
476 // dl->WriteData("OVERWRITE"); //jkl: wrong here
484 //_____________________________________________________________________________
485 Bool_t AliTRDclusterizer::ReadDigits()
488 // Reads the digits arrays from the input aliroot file
492 AliError("No run loader available");
496 AliLoader* loader = fRunLoader->GetLoader("TRDLoader");
497 if (!loader->TreeD()) {
498 loader->LoadDigits();
501 // Read in the digit arrays
502 return (fDigitsManager->ReadDigits(loader->TreeD()));
506 //_____________________________________________________________________________
507 Bool_t AliTRDclusterizer::ReadDigits(TTree *digitsTree)
510 // Reads the digits arrays from the input tree
513 // Read in the digit arrays
514 return (fDigitsManager->ReadDigits(digitsTree));
518 //_____________________________________________________________________________
519 Bool_t AliTRDclusterizer::ReadDigits(AliRawReader *rawReader)
522 // Reads the digits arrays from the ddl file
526 fDigitsManager = raw.Raw2Digits(rawReader);
532 //_____________________________________________________________________________
533 Bool_t AliTRDclusterizer::MakeClusters()
536 // Creates clusters from digits
539 // Propagate info from the digits manager
540 if (TestBit(kLabels)){
541 SetBit(kLabels, fDigitsManager->UsesDictionaries());
544 Bool_t fReturn = kTRUE;
545 for (Int_t i = 0; i < AliTRDgeometry::kNdet; i++){
547 AliTRDarrayADC *digitsIn = (AliTRDarrayADC*) fDigitsManager->GetDigits(i); //mod
548 // This is to take care of switched off super modules
549 if (!digitsIn->HasData()) continue;
551 digitsIn->DeleteNegatives(); // Restore digits array to >=0 values
552 AliTRDSignalIndex* indexes = fDigitsManager->GetIndexes(i);
553 if (indexes->IsAllocated() == kFALSE){
554 fDigitsManager->BuildIndexes(i);
558 if (indexes->HasEntry()){
559 if (TestBit(kLabels)){
560 for (Int_t iDict = 0; iDict < AliTRDdigitsManager::kNDict; iDict++){
561 AliTRDarrayDictionary *tracksIn = 0; //mod
562 tracksIn = (AliTRDarrayDictionary *) fDigitsManager->GetDictionary(i,iDict); //mod
566 fR = MakeClusters(i);
567 fReturn = fR && fReturn;
571 // if(IsWritingClusters()) WriteClusters(i);
575 // Clear arrays of this chamber, to prepare for next event
576 fDigitsManager->ClearArrays(i);
579 if(fReconstructor->IsWritingClusters()) WriteClusters(-1);
581 AliInfo(Form("Number of found clusters : %d", RecPoints()->GetEntriesFast()));
587 //_____________________________________________________________________________
588 Bool_t AliTRDclusterizer::Raw2Clusters(AliRawReader *rawReader)
591 // Creates clusters from raw data
594 return Raw2ClustersChamber(rawReader);
598 //_____________________________________________________________________________
599 Bool_t AliTRDclusterizer::Raw2ClustersChamber(AliRawReader *rawReader)
602 // Creates clusters from raw data
605 // Create the digits manager
606 if (!fDigitsManager){
607 SetBit(knewDM, kTRUE);
608 fDigitsManager = new AliTRDdigitsManager(kTRUE);
609 fDigitsManager->CreateArrays();
612 fDigitsManager->SetUseDictionaries(TestBit(kLabels));
614 // ----- preparing tracklet output -----
615 if (fReconstructor->IsWritingTracklets()) {
616 AliDataLoader *trklLoader = AliRunLoader::Instance()->GetLoader("TRDLoader")->GetDataLoader("tracklets");
618 //AliInfo("Could not get the tracklets data loader, adding it now!");
619 trklLoader = new AliDataLoader("TRD.Tracklets.root","tracklets", "tracklets");
620 AliRunLoader::Instance()->GetLoader("TRDLoader")->AddDataLoader(trklLoader);
622 AliTreeLoader *trklTreeLoader = dynamic_cast<AliTreeLoader*> (trklLoader->GetBaseLoader("tracklets-raw"));
623 if (!trklTreeLoader) {
624 trklTreeLoader = new AliTreeLoader("tracklets-raw", trklLoader);
625 trklLoader->AddBaseLoader(trklTreeLoader);
627 if (!trklTreeLoader->Tree())
628 trklTreeLoader->MakeTree();
631 // tracklet container for raw tracklet writing
632 if (!fTrackletContainer && ( fReconstructor->IsWritingTracklets() || fReconstructor->IsProcessingTracklets() )) {
633 // maximum tracklets for one HC
634 const Int_t kTrackletChmb=256;
635 fTrackletContainer = new UInt_t *[2];
636 fTrackletContainer[0] = new UInt_t[kTrackletChmb];
637 fTrackletContainer[1] = new UInt_t[kTrackletChmb];
638 memset(fTrackletContainer[0], 0, kTrackletChmb*sizeof(UInt_t)); //jkl
639 memset(fTrackletContainer[1], 0, kTrackletChmb*sizeof(UInt_t)); //jkl
643 fRawStream = AliTRDrawStreamBase::GetRawStream(rawReader);
645 fRawStream->SetReader(rawReader);
647 if(fReconstructor->IsHLT()){
648 if(fRawStream->InheritsFrom(AliTRDrawStream::Class()))
649 ((AliTRDrawStream*)fRawStream)->DisableErrorStorage();
651 fRawStream->SetSharedPadReadout(kFALSE);
652 fRawStream->SetNoErrorWarning();
656 AliDebug(1,Form("Stream version: %s", fRawStream->IsA()->GetName()));
659 while ((det = fRawStream->NextChamber(fDigitsManager,fTrackletContainer)) < AliTRDgeometry::kNdet){
660 if (fDigitsManager->GetIndexes(det)->HasEntry()){
662 fDigitsManager->ClearArrays(det);
665 if (!fReconstructor->IsWritingTracklets()) continue;
666 if (*(fTrackletContainer[0]) > 0 || *(fTrackletContainer[1]) > 0) WriteTracklets(det);
669 if (fReconstructor->IsWritingTracklets()) {
670 if (AliDataLoader *trklLoader = AliRunLoader::Instance()->GetLoader("TRDLoader")->GetDataLoader("tracklets")) {
672 if (AliTreeLoader *trklTreeLoader = (AliTreeLoader*) trklLoader->GetBaseLoader("tracklets-raw"))
673 trklTreeLoader->WriteData("OVERWRITE");
674 trklLoader->UnloadAll();
679 if (fTrackletContainer){
680 delete [] fTrackletContainer[0];
681 delete [] fTrackletContainer[1];
682 delete [] fTrackletContainer;
683 fTrackletContainer = NULL;
686 if(fReconstructor->IsWritingClusters()) WriteClusters(-1);
688 if(!TestBit(knewDM)){
689 delete fDigitsManager;
690 fDigitsManager = NULL;
695 AliInfo(Form("Number of found clusters : %d", fNoOfClusters));
700 //_____________________________________________________________________________
701 UChar_t AliTRDclusterizer::GetStatus(Short_t &signal)
704 // Check if a pad is masked
709 if(signal>0 && TESTBIT(signal, 10)){
711 for(int ibit=0; ibit<4; ibit++){
712 if(TESTBIT(signal, 11+ibit)){
713 SETBIT(status, ibit);
714 CLRBIT(signal, 11+ibit);
721 //_____________________________________________________________________________
722 void AliTRDclusterizer::SetPadStatus(const UChar_t status, UChar_t &out) const {
724 // Set the pad status into out
725 // First three bits are needed for the position encoding
730 //_____________________________________________________________________________
731 UChar_t AliTRDclusterizer::GetPadStatus(UChar_t encoding) const {
733 // return the staus encoding of the corrupted pad
735 return static_cast<UChar_t>(encoding >> 3);
738 //_____________________________________________________________________________
739 Int_t AliTRDclusterizer::GetCorruption(UChar_t encoding) const {
741 // Return the position of the corruption
746 //_____________________________________________________________________________
747 Bool_t AliTRDclusterizer::MakeClusters(Int_t det)
750 // Generates the cluster.
754 fDigits = (AliTRDarrayADC *) fDigitsManager->GetDigits(det); //mod
755 fBaseline = fDigitsManager->GetDigitsParam()->GetADCbaseline(det);
757 // This is to take care of switched off super modules
758 if (!fDigits->HasData()) return kFALSE;
760 fIndexes = fDigitsManager->GetIndexes(det);
761 if (fIndexes->IsAllocated() == kFALSE) {
762 AliError("Indexes do not exist!");
766 AliTRDcalibDB* const calibration = AliTRDcalibDB::Instance();
768 AliFatal("No AliTRDcalibDB instance available\n");
772 if (!fReconstructor){
773 AliError("Reconstructor not set\n");
777 const AliTRDrecoParam *const recoParam = fReconstructor->GetRecoParam();
779 fMaxThresh = (Short_t)recoParam->GetClusMaxThresh();
780 fMaxThreshTest = (Short_t)(recoParam->GetClusMaxThresh()/2+fBaseline);
781 fSigThresh = (Short_t)recoParam->GetClusSigThresh();
782 fMinMaxCutSigma = recoParam->GetMinMaxCutSigma();
783 fMinLeftRightCutSigma = recoParam->GetMinLeftRightCutSigma();
784 const Int_t iEveryNTB = recoParam->GetRecEveryNTB();
786 Int_t istack = fIndexes->GetStack();
787 fLayer = fIndexes->GetLayer();
788 Int_t isector = fIndexes->GetSM();
790 // Start clustering in the chamber
792 fDet = AliTRDgeometry::GetDetector(fLayer,istack,isector);
794 AliError("Strange Detector number Missmatch!");
798 AliDebug(2, Form("Det[%d] @ Sec[%d] Stk[%d] Ly[%d]", fDet, isector, istack, fLayer));
800 // TRD space point transformation
801 fTransform->SetDetector(det);
803 Int_t iGeoLayer = AliGeomManager::kTRD1 + fLayer;
804 Int_t iGeoModule = istack + AliTRDgeometry::Nstack() * isector;
805 fVolid = AliGeomManager::LayerToVolUID(iGeoLayer,iGeoModule);
807 if(fReconstructor->IsProcessingTracklets() && fTrackletContainer)
808 AddTrackletsToArray();
810 fColMax = fDigits->GetNcol();
811 fTimeTotal = fDigitsManager->GetDigitsParam()->GetNTimeBins(det);
813 // Check consistency between OCDB and raw data
814 Int_t nTimeOCDB = calibration->GetNumberOfTimeBinsDCS();
815 if(fReconstructor->IsHLT()){
816 if((nTimeOCDB > -1) && (fTimeTotal != nTimeOCDB)){
817 AliWarning(Form("Number of timebins does not match OCDB value (RAW[%d] OCDB[%d]), using raw value"
818 ,fTimeTotal,nTimeOCDB));
822 AliWarning("Undefined number of timebins in OCDB, using value from raw data.");
824 AliError("Number of timebins in raw data is negative, skipping chamber!");
827 }else if(nTimeOCDB == -2){
828 AliError("Mixed number of timebins in OCDB, no reconstruction of TRD data!");
830 }else if(fTimeTotal != nTimeOCDB){
831 AliError(Form("Number of timebins in raw data does not match OCDB value (RAW[%d] OCDB[%d]), skipping chamber!"
832 ,fTimeTotal,nTimeOCDB));
837 // Detector wise calibration object for the gain factors
838 const AliTRDCalDet *calGainFactorDet = calibration->GetGainFactorDet();
839 // Calibration object with pad wise values for the gain factors
840 fCalGainFactorROC = calibration->GetGainFactorROC(fDet);
841 // Calibration value for chamber wise gain factor
842 fCalGainFactorDetValue = calGainFactorDet->GetValue(fDet);
844 // Detector wise calibration object for the noise
845 const AliTRDCalDet *calNoiseDet = calibration->GetNoiseDet();
846 // Calibration object with pad wise values for the noise
847 fCalNoiseROC = calibration->GetNoiseROC(fDet);
848 // Calibration value for chamber wise noise
849 fCalNoiseDetValue = calNoiseDet->GetValue(fDet);
851 // Calibration object with the pad status
852 fCalPadStatusROC = calibration->GetPadStatusROC(fDet);
854 firstClusterROC = -1;
857 SetBit(kLUT, recoParam->UseLUT());
858 SetBit(kGAUS, recoParam->UseGAUS());
860 // Apply the gain and the tail cancelation via digital filter
861 if(recoParam->UseTailCancelation()) TailCancelation(recoParam);
863 MaxStruct curr, last;
864 Int_t nMaximas = 0, nCorrupted = 0;
866 // Here the clusterfining is happening
868 for(curr.time = 0; curr.time < fTimeTotal; curr.time+=iEveryNTB){
869 while(fIndexes->NextRCIndex(curr.row, curr.col)){
870 if(fDigits->GetData(curr.row, curr.col, curr.time) > fMaxThreshTest && IsMaximum(curr, curr.padStatus, &curr.signals[0])){
872 if(curr.col==last.col+2 && curr.row==last.row && curr.time==last.time) FivePadCluster(last, curr);
875 last=curr; curr.fivePad=kFALSE;
879 if(last.row>-1) CreateCluster(last);
881 if(recoParam->GetStreamLevel(AliTRDrecoParam::kClusterizer) > 2 && fReconstructor->IsDebugStreaming()){
882 TTreeSRedirector* fDebugStream = fReconstructor->GetDebugStream(AliTRDrecoParam::kClusterizer);
883 (*fDebugStream) << "MakeClusters"
884 << "Detector=" << det
885 << "NMaxima=" << nMaximas
886 << "NClusters=" << fClusterROC
887 << "NCorrupted=" << nCorrupted
890 if (TestBit(kLabels)) AddLabels();
896 //_____________________________________________________________________________
897 Bool_t AliTRDclusterizer::IsMaximum(const MaxStruct &Max, UChar_t &padStatus, Short_t *const Signals)
900 // Returns true if this row,col,time combination is a maximum.
901 // Gives back the padStatus and the signals of the center pad and the two neighbouring pads.
904 Float_t gain = fCalGainFactorDetValue * fCalGainFactorROC->GetValue(Max.col,Max.row);
905 Signals[1] = (Short_t)((fDigits->GetData(Max.row, Max.col, Max.time) - fBaseline) / gain + 0.5f);
906 if(Signals[1] <= fMaxThresh) return kFALSE;
908 if(Max.col < 1 || Max.col + 1 >= fColMax) return kFALSE;
910 Short_t noiseMiddleThresh = (Short_t)(fMinMaxCutSigma*fCalNoiseDetValue*fCalNoiseROC->GetValue(Max.col, Max.row));
911 if (Signals[1] <= noiseMiddleThresh) return kFALSE;
914 fCalPadStatusROC->GetStatus(Max.col-1, Max.row)
915 ,fCalPadStatusROC->GetStatus(Max.col, Max.row)
916 ,fCalPadStatusROC->GetStatus(Max.col+1, Max.row)
919 gain = fCalGainFactorDetValue * fCalGainFactorROC->GetValue(Max.col-1,Max.row);
920 Signals[0] = (Short_t)((fDigits->GetData(Max.row, Max.col-1, Max.time) - fBaseline) / gain + 0.5f);
921 gain = fCalGainFactorDetValue * fCalGainFactorROC->GetValue(Max.col+1,Max.row);
922 Signals[2] = (Short_t)((fDigits->GetData(Max.row, Max.col+1, Max.time) - fBaseline) / gain + 0.5f);
924 if(!(status[0] | status[1] | status[2])) {//all pads are good
925 if ((Signals[2] <= Signals[1]) && (Signals[0] < Signals[1])) {
926 if ((Signals[2] > fSigThresh) || (Signals[0] > fSigThresh)) {
927 if(Signals[0]<0)Signals[0]=0;
928 if(Signals[2]<0)Signals[2]=0;
929 Short_t noiseSumThresh = (Short_t)(fMinLeftRightCutSigma * fCalNoiseDetValue
930 * fCalNoiseROC->GetValue(Max.col, Max.row));
931 if ((Signals[2]+Signals[0]+Signals[1]) <= noiseSumThresh) return kFALSE;
936 } else { // at least one of the pads is bad, and reject candidates with more than 1 problematic pad
937 if(Signals[0]<0)Signals[0]=0;
938 if(Signals[2]<0)Signals[2]=0;
939 if (status[2] && (!(status[0] || status[1])) && Signals[1] > Signals[0] && Signals[0] > fSigThresh) {
941 SetPadStatus(status[2], padStatus);
944 else if (status[0] && (!(status[1] || status[2])) && Signals[1] >= Signals[2] && Signals[2] > fSigThresh) {
946 SetPadStatus(status[0], padStatus);
949 else if (status[1] && (!(status[0] || status[2])) && ((Signals[2] > fSigThresh) || (Signals[0] > fSigThresh))) {
950 Signals[1] = fMaxThresh;
951 SetPadStatus(status[1], padStatus);
958 //_____________________________________________________________________________
959 Bool_t AliTRDclusterizer::FivePadCluster(MaxStruct &ThisMax, MaxStruct &NeighbourMax)
962 // Look for 5 pad cluster with minimum in the middle
963 // Gives back the ratio
966 if (ThisMax.col >= fColMax - 3) return kFALSE;
968 if (ThisMax.col < fColMax - 5){
969 gain = fCalGainFactorDetValue * fCalGainFactorROC->GetValue(ThisMax.col+4,ThisMax.row);
970 if (fDigits->GetData(ThisMax.row, ThisMax.col+4, ThisMax.time) - fBaseline >= fSigThresh * gain)
973 if (ThisMax.col > 1) {
974 gain = fCalGainFactorDetValue * fCalGainFactorROC->GetValue(ThisMax.col-2,ThisMax.row);
975 if (fDigits->GetData(ThisMax.row, ThisMax.col-2, ThisMax.time) - fBaseline >= fSigThresh * gain)
979 const Float_t kEpsilon = 0.01;
980 Double_t padSignal[5] = {ThisMax.signals[0], ThisMax.signals[1], ThisMax.signals[2],
981 NeighbourMax.signals[1], NeighbourMax.signals[2]};
983 // Unfold the two maxima and set the signal on
984 // the overlapping pad to the ratio
985 Float_t ratio = Unfold(kEpsilon,fLayer,padSignal);
986 ThisMax.signals[2] = (Short_t)(ThisMax.signals[2]*ratio + 0.5f);
987 NeighbourMax.signals[0] = (Short_t)(NeighbourMax.signals[0]*(1-ratio) + 0.5f);
988 ThisMax.fivePad=kTRUE;
989 NeighbourMax.fivePad=kTRUE;
994 //_____________________________________________________________________________
995 void AliTRDclusterizer::CreateCluster(const MaxStruct &Max)
998 // Creates a cluster at the given position and saves it in fRecPoints
1001 Int_t nPadCount = 1;
1002 Short_t signals[7] = { 0, 0, Max.signals[0], Max.signals[1], Max.signals[2], 0, 0 };
1003 if(!fReconstructor->IsHLT()) CalcAdditionalInfo(Max, signals, nPadCount);
1005 AliTRDcluster cluster(fDet, ((UChar_t) Max.col), ((UChar_t) Max.row), ((UChar_t) Max.time), signals, fVolid);
1006 cluster.SetNPads(nPadCount);
1007 if(TestBit(kLUT)) cluster.SetRPhiMethod(AliTRDcluster::kLUT);
1008 else if(TestBit(kGAUS)) cluster.SetRPhiMethod(AliTRDcluster::kGAUS);
1009 else cluster.SetRPhiMethod(AliTRDcluster::kCOG);
1011 cluster.SetFivePad(Max.fivePad);
1012 // set pads status for the cluster
1013 UChar_t maskPosition = GetCorruption(Max.padStatus);
1015 cluster.SetPadMaskedPosition(maskPosition);
1016 cluster.SetPadMaskedStatus(GetPadStatus(Max.padStatus));
1018 cluster.SetXcorr(fReconstructor->UseClusterRadialCorrection());
1020 // Transform the local cluster coordinates into calibrated
1021 // space point positions defined in the local tracking system.
1022 // Here the calibration for T0, Vdrift and ExB is applied as well.
1023 if(!TestBit(kSkipTrafo)) if(!fTransform->Transform(&cluster)) return;
1025 // Temporarily store the Max.Row, column and time bin of the center pad
1026 // Used to later on assign the track indices
1027 cluster.SetLabel(Max.row, 0);
1028 cluster.SetLabel(Max.col, 1);
1029 cluster.SetLabel(Max.time,2);
1031 //needed for HLT reconstruction
1032 AddClusterToArray(&cluster);
1034 // Store the index of the first cluster in the current ROC
1035 if (firstClusterROC < 0) firstClusterROC = fNoOfClusters;
1041 //_____________________________________________________________________________
1042 void AliTRDclusterizer::CalcAdditionalInfo(const MaxStruct &Max, Short_t *const signals, Int_t &nPadCount)
1044 // Look to the right
1046 while (fDigits->GetData(Max.row, Max.col-ii, Max.time) >= fSigThresh) {
1049 if (Max.col < ii) break;
1053 while (fDigits->GetData(Max.row, Max.col+ii, Max.time) >= fSigThresh) {
1056 if (Max.col+ii >= fColMax) break;
1059 // Store the amplitudes of the pads in the cluster for later analysis
1060 // and check whether one of these pads is masked in the database
1061 signals[2]=Max.signals[0];
1062 signals[3]=Max.signals[1];
1063 signals[4]=Max.signals[2];
1065 for(Int_t i = 0; i<2; i++)
1068 gain = fCalGainFactorDetValue * fCalGainFactorROC->GetValue(Max.col-3+i,Max.row);
1069 signals[i] = (Short_t)((fDigits->GetData(Max.row, Max.col-3+i, Max.time) - fBaseline) / gain + 0.5f);
1071 if(Max.col+3-i < fColMax){
1072 gain = fCalGainFactorDetValue * fCalGainFactorROC->GetValue(Max.col+3-i,Max.row);
1073 signals[6-i] = (Short_t)((fDigits->GetData(Max.row, Max.col+3-i, Max.time) - fBaseline) / gain + 0.5f);
1076 /*for (Int_t jPad = Max.Col-3; jPad <= Max.Col+3; jPad++) {
1077 if ((jPad >= 0) && (jPad < fColMax))
1078 signals[jPad-Max.Col+3] = TMath::Nint(fDigits->GetData(Max.Row,jPad,Max.Time));
1082 //_____________________________________________________________________________
1083 void AliTRDclusterizer::AddClusterToArray(AliTRDcluster* cluster)
1086 // Add a cluster to the array
1089 Int_t n = RecPoints()->GetEntriesFast();
1090 if(n!=fNoOfClusters)AliError(Form("fNoOfClusters != RecPoints()->GetEntriesFast %i != %i \n", fNoOfClusters, n));
1091 new((*RecPoints())[n]) AliTRDcluster(*cluster);
1094 //_____________________________________________________________________________
1095 void AliTRDclusterizer::AddTrackletsToArray()
1098 // Add the online tracklets of this chamber to the array
1101 UInt_t* trackletword;
1102 for(Int_t side=0; side<2; side++)
1105 trackletword=fTrackletContainer[side];
1106 while(trackletword[trkl]>0){
1107 Int_t n = TrackletsArray()->GetEntriesFast();
1108 AliTRDtrackletWord tmp(trackletword[trkl]);
1109 new((*TrackletsArray())[n]) AliTRDcluster(&tmp,fDet,fVolid);
1115 //_____________________________________________________________________________
1116 Bool_t AliTRDclusterizer::AddLabels()
1119 // Add the track indices to the found clusters
1122 const Int_t kNclus = 3;
1123 const Int_t kNdict = AliTRDdigitsManager::kNDict;
1124 const Int_t kNtrack = kNdict * kNclus;
1126 Int_t iClusterROC = 0;
1133 // Temporary array to collect the track indices
1134 Int_t *idxTracks = new Int_t[kNtrack*fClusterROC];
1136 // Loop through the dictionary arrays one-by-one
1137 // to keep memory consumption low
1138 AliTRDarrayDictionary *tracksIn = 0; //mod
1139 for (Int_t iDict = 0; iDict < kNdict; iDict++) {
1141 // tracksIn should be expanded beforehand!
1142 tracksIn = (AliTRDarrayDictionary *) fDigitsManager->GetDictionary(fDet,iDict);
1144 // Loop though the clusters found in this ROC
1145 for (iClusterROC = 0; iClusterROC < fClusterROC; iClusterROC++) {
1147 AliTRDcluster *cluster = (AliTRDcluster *)
1148 RecPoints()->UncheckedAt(firstClusterROC+iClusterROC);
1149 row = cluster->GetLabel(0);
1150 col = cluster->GetLabel(1);
1151 time = cluster->GetLabel(2);
1153 for (iPad = 0; iPad < kNclus; iPad++) {
1154 Int_t iPadCol = col - 1 + iPad;
1155 Int_t index = tracksIn->GetData(row,iPadCol,time); //Modification of -1 in Track
1156 idxTracks[3*iPad+iDict + iClusterROC*kNtrack] = index;
1163 // Copy the track indices into the cluster
1164 // Loop though the clusters found in this ROC
1165 for (iClusterROC = 0; iClusterROC < fClusterROC; iClusterROC++) {
1167 AliTRDcluster *cluster = (AliTRDcluster *)
1168 RecPoints()->UncheckedAt(firstClusterROC+iClusterROC);
1169 cluster->SetLabel(-9999,0);
1170 cluster->SetLabel(-9999,1);
1171 cluster->SetLabel(-9999,2);
1173 cluster->AddTrackIndex(&idxTracks[iClusterROC*kNtrack]);
1177 delete [] idxTracks;
1183 //_____________________________________________________________________________
1184 Float_t AliTRDclusterizer::Unfold(Double_t eps, Int_t layer, const Double_t *const padSignal) const
1187 // Method to unfold neighbouring maxima.
1188 // The charge ratio on the overlapping pad is calculated
1189 // until there is no more change within the range given by eps.
1190 // The resulting ratio is then returned to the calling method.
1193 AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
1195 AliError("No AliTRDcalibDB instance available\n");
1200 Int_t itStep = 0; // Count iteration steps
1202 Double_t ratio = 0.5; // Start value for ratio
1203 Double_t prevRatio = 0.0; // Store previous ratio
1205 Double_t newLeftSignal[3] = { 0.0, 0.0, 0.0 }; // Array to store left cluster signal
1206 Double_t newRightSignal[3] = { 0.0, 0.0, 0.0 }; // Array to store right cluster signal
1207 Double_t newSignal[3] = { 0.0, 0.0, 0.0 };
1209 // Start the iteration
1210 while ((TMath::Abs(prevRatio - ratio) > eps) && (itStep < 10)) {
1215 // Cluster position according to charge ratio
1216 Double_t maxLeft = (ratio*padSignal[2] - padSignal[0])
1217 / (padSignal[0] + padSignal[1] + ratio * padSignal[2]);
1218 Double_t maxRight = (padSignal[4] - (1-ratio)*padSignal[2])
1219 / ((1.0 - ratio)*padSignal[2] + padSignal[3] + padSignal[4]);
1221 // Set cluster charge ratio
1222 irc = calibration->PadResponse(1.0, maxLeft, layer, newSignal);
1223 Double_t ampLeft = padSignal[1] / newSignal[1];
1224 irc = calibration->PadResponse(1.0, maxRight, layer, newSignal);
1225 Double_t ampRight = padSignal[3] / newSignal[1];
1227 // Apply pad response to parameters
1228 irc = calibration->PadResponse(ampLeft ,maxLeft ,layer,newLeftSignal );
1229 irc = calibration->PadResponse(ampRight,maxRight,layer,newRightSignal);
1231 // Calculate new overlapping ratio
1232 ratio = TMath::Min((Double_t) 1.0
1233 ,newLeftSignal[2] / (newLeftSignal[2] + newRightSignal[0]));
1241 //_____________________________________________________________________________
1242 void AliTRDclusterizer::TailCancelation(const AliTRDrecoParam* const recoParam)
1245 // Applies the tail cancelation
1252 TTreeSRedirector *fDebugStream = fReconstructor->GetDebugStream(AliTRDrecoParam::kClusterizer);
1253 Bool_t debugStreaming = recoParam->GetStreamLevel(AliTRDrecoParam::kClusterizer) > 7 && fReconstructor->IsDebugStreaming();
1254 Int_t nexp = recoParam->GetTCnexp();
1255 while(fIndexes->NextRCIndex(iRow, iCol))
1257 // if corrupted then don't make the tail cancallation
1258 if (fCalPadStatusROC->GetStatus(iCol, iRow)) continue;
1261 for (iTime = 0; iTime < fTimeTotal; iTime++)
1262 (*fDebugStream) << "TailCancellation"
1268 // Apply the tail cancelation via the digital filter
1269 DeConvExp(fDigits->GetDataAddress(iRow,iCol),fTimeTotal,nexp);
1271 } // while irow icol
1277 //_____________________________________________________________________________
1278 void AliTRDclusterizer::DeConvExp(Short_t *const arr, const Int_t nTime, const Int_t nexp)
1281 // Tail cancellation by deconvolution for PASA v4 TRF
1285 Float_t coefficients[2];
1287 // Initialization (coefficient = alpha, rates = lambda)
1293 if (nexp == 1) { // 1 Exponentials
1299 if (nexp == 2) { // 2 Exponentials
1301 fReconstructor->GetRecoParam()->GetTCParams(par);
1302 r1 = par[0];//1.156;
1303 r2 = par[1];//0.130;
1304 c1 = par[2];//0.114;
1305 c2 = par[3];//0.624;
1308 coefficients[0] = c1;
1309 coefficients[1] = c2;
1313 rates[0] = TMath::Exp(-dt/(r1));
1314 rates[1] = (nexp == 1) ? .0 : TMath::Exp(-dt/(r2));
1316 Float_t reminder[2] = { .0, .0 };
1317 Float_t correction = 0.0;
1318 Float_t result = 0.0;
1320 for (int i = 0; i < nTime; i++) {
1322 result = arr[i] - correction - fBaseline; // No rescaling
1323 arr[i] = (Short_t)(result + fBaseline + 0.5f);
1326 for (int k = 0; k < 2; k++) {
1327 correction += reminder[k] = rates[k] * (reminder[k] + coefficients[k] * result);
1334 //_____________________________________________________________________________
1335 void AliTRDclusterizer::ResetRecPoints()
1338 // Resets the list of rec points
1342 fRecPoints->Clear();
1344 // delete fRecPoints;
1348 //_____________________________________________________________________________
1349 TClonesArray *AliTRDclusterizer::RecPoints()
1352 // Returns the list of rec points
1356 if(!(fRecPoints = AliTRDReconstructor::GetClusters())){
1357 // determine number of clusters which has to be allocated
1358 Float_t nclusters = fReconstructor->GetRecoParam()->GetNClusters();
1360 fRecPoints = new TClonesArray("AliTRDcluster", Int_t(nclusters));
1362 //SetClustersOwner(kTRUE);
1363 AliTRDReconstructor::SetClusters(0x0);
1369 //_____________________________________________________________________________
1370 TClonesArray *AliTRDclusterizer::TrackletsArray()
1373 // Returns the list of rec points
1376 if (!fTracklets && fReconstructor->IsProcessingTracklets()) {
1377 fTracklets = new TClonesArray("AliTRDcluster", 2*MAXTRACKLETSPERHC);
1378 //SetClustersOwner(kTRUE);
1379 //AliTRDReconstructor::SetTracklets(0x0);