2 /**************************************************************************
3 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5 * Author: The ALICE Off-line Project. *
6 * Contributors are mentioned in the code where appropriate. *
8 * Permission to use, copy, modify and distribute this software and its *
9 * documentation strictly for non-commercial purposes is hereby granted *
10 * without fee, provided that the above copyright notice appears in all *
11 * copies and that both the copyright notice and this permission notice *
12 * appear in the supporting documentation. The authors make no claims *
13 * about the suitability of this software for any purpose. It is *
14 * provided "as is" without express or implied warranty. *
15 **************************************************************************/
19 ///////////////////////////////////////////////////////////////////////////////
21 // TRD cluster finder //
23 ///////////////////////////////////////////////////////////////////////////////
31 #include "AliRunLoader.h"
32 #include "AliLoader.h"
33 #include "AliRawReader.h"
36 #include "AliTRDclusterizerV1.h"
37 #include "AliTRDgeometry.h"
38 #include "AliTRDdataArrayF.h"
39 #include "AliTRDdataArrayI.h"
40 #include "AliTRDdigitsManager.h"
41 #include "AliTRDpadPlane.h"
42 #include "AliTRDrawData.h"
43 #include "AliTRDcalibDB.h"
44 #include "AliTRDSimParam.h"
45 #include "AliTRDRecParam.h"
46 #include "AliTRDCommonParam.h"
47 #include "AliTRDcluster.h"
49 ClassImp(AliTRDclusterizerV1)
51 //_____________________________________________________________________________
52 AliTRDclusterizerV1::AliTRDclusterizerV1()
57 // AliTRDclusterizerV1 default constructor
62 //_____________________________________________________________________________
63 AliTRDclusterizerV1::AliTRDclusterizerV1(const Text_t* name, const Text_t* title)
64 :AliTRDclusterizer(name,title)
65 ,fDigitsManager(new AliTRDdigitsManager())
68 // AliTRDclusterizerV1 constructor
71 fDigitsManager->CreateArrays();
75 //_____________________________________________________________________________
76 AliTRDclusterizerV1::AliTRDclusterizerV1(const AliTRDclusterizerV1 &c)
81 // AliTRDclusterizerV1 copy constructor
86 //_____________________________________________________________________________
87 AliTRDclusterizerV1::~AliTRDclusterizerV1()
90 // AliTRDclusterizerV1 destructor
94 delete fDigitsManager;
95 fDigitsManager = NULL;
100 //_____________________________________________________________________________
101 AliTRDclusterizerV1 &AliTRDclusterizerV1::operator=(const AliTRDclusterizerV1 &c)
104 // Assignment operator
107 if (this != &c) ((AliTRDclusterizerV1 &) c).Copy(*this);
112 //_____________________________________________________________________________
113 void AliTRDclusterizerV1::Copy(TObject &c) const
119 ((AliTRDclusterizerV1 &) c).fDigitsManager = 0;
121 AliTRDclusterizer::Copy(c);
125 //_____________________________________________________________________________
126 Bool_t AliTRDclusterizerV1::ReadDigits()
129 // Reads the digits arrays from the input aliroot file
133 AliError("No run loader available");
137 AliLoader* loader = fRunLoader->GetLoader("TRDLoader");
138 if (!loader->TreeD()) {
139 loader->LoadDigits();
142 // Read in the digit arrays
143 return (fDigitsManager->ReadDigits(loader->TreeD()));
147 //_____________________________________________________________________________
148 Bool_t AliTRDclusterizerV1::ReadDigits(AliRawReader* rawReader)
151 // Reads the digits arrays from the ddl file
155 fDigitsManager = raw.Raw2Digits(rawReader);
161 //_____________________________________________________________________________
162 Bool_t AliTRDclusterizerV1::MakeClusters()
165 // Generates the cluster.
176 AliTRDdataArrayI *digitsIn;
177 AliTRDdataArrayI *track0;
178 AliTRDdataArrayI *track1;
179 AliTRDdataArrayI *track2;
182 AliTRDgeometry *geo = AliTRDgeometry::GetGeometry(fRunLoader);
183 AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
185 AliError("No AliTRDcalibDB instance available\n");
189 AliTRDSimParam *simParam = AliTRDSimParam::Instance();
191 AliError("No AliTRDSimParam instance available\n");
195 AliTRDRecParam *recParam = AliTRDRecParam::Instance();
197 AliError("No AliTRDRecParam instance available\n");
201 AliTRDCommonParam *commonParam = AliTRDCommonParam::Instance();
203 AliError("Could not get common parameters\n");
208 Float_t ADCthreshold = simParam->GetADCthreshold();
209 // Threshold value for the maximum
210 Float_t maxThresh = recParam->GetClusMaxThresh();
211 // Threshold value for the digit signal
212 Float_t sigThresh = recParam->GetClusSigThresh();
214 // Iteration limit for unfolding procedure
215 const Float_t kEpsilon = 0.01;
216 const Int_t kNclus = 3;
217 const Int_t kNsig = 5;
218 const Int_t kNtrack = 3 * kNclus;
222 Double_t ratioLeft = 1.0;
223 Double_t ratioRight = 1.0;
225 Double_t padSignal[kNsig];
226 Double_t clusterSignal[kNclus];
227 Double_t clusterPads[kNclus];
228 Int_t clusterTracks[kNtrack];
231 Int_t chamEnd = AliTRDgeometry::Ncham();
233 Int_t planEnd = AliTRDgeometry::Nplan();
235 Int_t sectEnd = AliTRDgeometry::Nsect();
236 Int_t nTimeTotal = calibration->GetNumberOfTimeBins();
238 AliDebug(1,Form("Number of Time Bins = %d.\n",nTimeTotal));
240 // Start clustering in every chamber
241 for (icham = chamBeg; icham < chamEnd; icham++) {
242 for (iplan = planBeg; iplan < planEnd; iplan++) {
243 for (isect = sectBeg; isect < sectEnd; isect++) {
245 Int_t idet = geo->GetDetector(iplan,icham,isect);
247 Int_t nRowMax = commonParam->GetRowMax(iplan,icham,isect);
248 Int_t nColMax = commonParam->GetColMax(iplan);
250 AliTRDpadPlane *padPlane = commonParam->GetPadPlane(iplan,icham);
253 Int_t nClusters2pad = 0;
254 Int_t nClusters3pad = 0;
255 Int_t nClusters4pad = 0;
256 Int_t nClusters5pad = 0;
257 Int_t nClustersLarge = 0;
259 AliDebug(1,Form("Analyzing chamber %d, plane %d, sector %d.\n"
260 ,icham,iplan,isect));
263 digitsIn = fDigitsManager->GetDigits(idet);
264 // This is to take care of switched off super modules
265 if (digitsIn->GetNtime() == 0) {
269 track0 = fDigitsManager->GetDictionary(idet,0);
271 track1 = fDigitsManager->GetDictionary(idet,1);
273 track2 = fDigitsManager->GetDictionary(idet,2);
276 AliTRDdataArrayF *digitsOut = new AliTRDdataArrayF(digitsIn->GetNrow()
278 ,digitsIn->GetNtime());
279 Transform(digitsIn, digitsOut,idet,nRowMax,nColMax,nTimeTotal,ADCthreshold);
281 // Loop through the chamber and find the maxima
282 for ( row = 0; row < nRowMax; row++) {
283 for ( col = 2; col < nColMax; col++) {
284 for (time = 0; time < nTimeTotal; time++) {
286 Float_t signalL = TMath::Abs(digitsOut->GetDataUnchecked(row,col ,time));
287 Float_t signalM = TMath::Abs(digitsOut->GetDataUnchecked(row,col-1,time));
288 Float_t signalR = TMath::Abs(digitsOut->GetDataUnchecked(row,col-2,time));
290 // Look for the maximum
291 if (signalM >= maxThresh) {
292 if ((TMath::Abs(signalL) <= signalM) &&
293 (TMath::Abs(signalR) <= signalM) &&
294 ((TMath::Abs(signalL) + TMath::Abs(signalR)) > sigThresh)) {
295 // Maximum found, mark the position by a negative signal
296 digitsOut->SetDataUnchecked(row,col-1,time,-signalM);
304 // Now check the maxima and calculate the cluster position
305 for ( row = 0; row < nRowMax ; row++) {
306 for (time = 0; time < nTimeTotal; time++) {
307 for ( col = 1; col < nColMax-1; col++) {
310 if (digitsOut->GetDataUnchecked(row,col,time) < 0) {
312 for (iPad = 0; iPad < kNclus; iPad++) {
313 Int_t iPadCol = col - 1 + iPad;
314 clusterSignal[iPad] = TMath::Abs(digitsOut->GetDataUnchecked(row
317 clusterTracks[3*iPad ] = track0->GetDataUnchecked(row,iPadCol,time) - 1;
318 clusterTracks[3*iPad+1] = track1->GetDataUnchecked(row,iPadCol,time) - 1;
319 clusterTracks[3*iPad+2] = track2->GetDataUnchecked(row,iPadCol,time) - 1;
322 // Count the number of pads in the cluster
325 while (TMath::Abs(digitsOut->GetDataUnchecked(row,col-ii ,time)) >= sigThresh) {
328 if (col-ii < 0) break;
331 while (TMath::Abs(digitsOut->GetDataUnchecked(row,col+ii+1,time)) >= sigThresh) {
334 if (col+ii+1 >= nColMax) break;
361 // Look for 5 pad cluster with minimum in the middle
362 Bool_t fivePadCluster = kFALSE;
363 if (col < (nColMax - 3)) {
364 if (digitsOut->GetDataUnchecked(row,col+2,time) < 0) {
365 fivePadCluster = kTRUE;
367 if ((fivePadCluster) && (col < (nColMax - 5))) {
368 if (digitsOut->GetDataUnchecked(row,col+4,time) >= sigThresh) {
369 fivePadCluster = kFALSE;
372 if ((fivePadCluster) && (col > 1)) {
373 if (digitsOut->GetDataUnchecked(row,col-2,time) >= sigThresh) {
374 fivePadCluster = kFALSE;
380 // Modify the signal of the overlapping pad for the left part
381 // of the cluster which remains from a previous unfolding
383 clusterSignal[0] *= ratioLeft;
388 // Unfold the 5 pad cluster
389 if (fivePadCluster) {
390 for (iPad = 0; iPad < kNsig; iPad++) {
391 padSignal[iPad] = TMath::Abs(digitsOut->GetDataUnchecked(row
395 // Unfold the two maxima and set the signal on
396 // the overlapping pad to the ratio
397 ratioRight = Unfold(kEpsilon,iplan,padSignal);
398 ratioLeft = 1.0 - ratioRight;
399 clusterSignal[2] *= ratioRight;
404 Double_t clusterCharge = clusterSignal[0]
408 // The position of the cluster
409 clusterPads[0] = row + 0.5;
410 // Take the shift of the additional time bins into account
411 clusterPads[2] = time + 0.5;
413 if (recParam->LUTOn()) {
414 // Calculate the position of the cluster by using the
415 // lookup table method
416 clusterPads[1] = recParam->LUTposition(iplan,clusterSignal[0]
421 // Calculate the position of the cluster by using the
422 // center of gravity method
423 for (Int_t i = 0; i < 5; i++) {
426 padSignal[2] = TMath::Abs(digitsOut->GetDataUnchecked(row,col ,time)); // central pad
427 padSignal[1] = TMath::Abs(digitsOut->GetDataUnchecked(row,col-1,time)); // left pad
428 padSignal[3] = TMath::Abs(digitsOut->GetDataUnchecked(row,col+1,time)); // right pad
430 (TMath::Abs(digitsOut->GetDataUnchecked(row,col-2,time)) < padSignal[1])) {
431 padSignal[0] = TMath::Abs(digitsOut->GetDataUnchecked(row,col-2,time));
433 if ((col < nColMax - 3) &&
434 (TMath::Abs(digitsOut->GetDataUnchecked(row,col+2,time)) < padSignal[3])) {
435 padSignal[4] = TMath::Abs(digitsOut->GetDataUnchecked(row,col+2,time));
437 clusterPads[1] = GetCOG(padSignal);
440 Double_t q0 = clusterSignal[0];
441 Double_t q1 = clusterSignal[1];
442 Double_t q2 = clusterSignal[2];
443 Double_t clusterSigmaY2 = (q1*(q0+q2)+4*q0*q2) /
444 (clusterCharge*clusterCharge);
447 // Calculate the position and the error
451 Int_t clusterTimeBin = TMath::Nint(time - calibration->GetT0(idet, col, row));
453 Double_t colSize = padPlane->GetColSize(col);
454 Double_t rowSize = padPlane->GetRowSize(row);
456 Double_t clusterPos[3];
457 clusterPos[0] = padPlane->GetColPos(col) - (clusterPads[1]+0.5)*colSize;
458 clusterPos[1] = padPlane->GetRowPos(row) - 0.5*rowSize;
459 clusterPos[2] = CalcXposFromTimebin(clusterPads[2],idet,col,row);
460 Double_t clusterSig[2];
461 clusterSig[0] = (clusterSigmaY2 + 1./12.) * colSize*colSize;
462 clusterSig[1] = rowSize * rowSize / 12.;
465 // Add the cluster to the output array
466 AliTRDcluster * cluster = AddCluster(clusterPos
475 printf("Add a cluster: q=%f, det=%d, x=%f, y=%f, z=%f\n",clusterCharge
476 ,idet,clusterPos[0],clusterPos[1],clusterPos[2]);
478 Short_t signals[7]={ 0, 0, 0, 0, 0, 0, 0 };
479 for (Int_t jPad = col-3; jPad <= col+3; jPad++) {
480 if ((jPad < 0) || (jPad >= nColMax-1)) {
483 signals[jPad-col+3] = TMath::Nint(TMath::Abs(digitsOut->GetDataUnchecked(row,jPad,time)));
485 cluster->SetSignals(signals);
495 // Compress the arrays
496 track0->Compress(1,0);
497 track1->Compress(1,0);
498 track2->Compress(1,0);
500 // Write the cluster and reset the array
512 //_____________________________________________________________________________
513 Double_t AliTRDclusterizerV1::GetCOG(Double_t signal[5])
517 // Used for clusters with more than 3 pads - where LUT not applicable
520 Double_t sum = signal[0]+signal[1]+signal[2]+signal[3]+signal[4];
521 Double_t res = (0.0*(-signal[0]+signal[4])+(-signal[1]+signal[3]))/sum;
527 //_____________________________________________________________________________
528 Double_t AliTRDclusterizerV1::Unfold(Double_t eps, Int_t plane, Double_t* padSignal)
531 // Method to unfold neighbouring maxima.
532 // The charge ratio on the overlapping pad is calculated
533 // until there is no more change within the range given by eps.
534 // The resulting ratio is then returned to the calling method.
537 AliTRDcalibDB* calibration = AliTRDcalibDB::Instance();
539 AliError("No AliTRDcalibDB instance available\n");
544 Int_t itStep = 0; // Count iteration steps
546 Double_t ratio = 0.5; // Start value for ratio
547 Double_t prevRatio = 0; // Store previous ratio
549 Double_t newLeftSignal[3] = {0}; // Array to store left cluster signal
550 Double_t newRightSignal[3] = {0}; // Array to store right cluster signal
551 Double_t newSignal[3] = {0};
553 // Start the iteration
554 while ((TMath::Abs(prevRatio - ratio) > eps) && (itStep < 10)) {
559 // Cluster position according to charge ratio
560 Double_t maxLeft = (ratio*padSignal[2] - padSignal[0])
561 / (padSignal[0] + padSignal[1] + ratio*padSignal[2]);
562 Double_t maxRight = (padSignal[4] - (1-ratio)*padSignal[2])
563 / ((1-ratio)*padSignal[2] + padSignal[3] + padSignal[4]);
565 // Set cluster charge ratio
566 irc = calibration->PadResponse(1.0,maxLeft ,plane,newSignal);
567 Double_t ampLeft = padSignal[1] / newSignal[1];
568 irc = calibration->PadResponse(1.0,maxRight,plane,newSignal);
569 Double_t ampRight = padSignal[3] / newSignal[1];
571 // Apply pad response to parameters
572 irc = calibration->PadResponse(ampLeft ,maxLeft ,plane,newLeftSignal );
573 irc = calibration->PadResponse(ampRight,maxRight,plane,newRightSignal);
575 // Calculate new overlapping ratio
576 ratio = TMath::Min((Double_t)1.0,newLeftSignal[2] /
577 (newLeftSignal[2] + newRightSignal[0]));
585 //_____________________________________________________________________________
586 void AliTRDclusterizerV1::Transform(AliTRDdataArrayI* digitsIn,
587 AliTRDdataArrayF* digitsOut,
588 Int_t idet, Int_t nRowMax,
589 Int_t nColMax, Int_t nTimeTotal,
590 Float_t ADCthreshold)
594 // Apply tail cancellation: Transform digitsIn to digitsOut
601 AliTRDRecParam* recParam = AliTRDRecParam::Instance();
603 AliError("No AliTRDRecParam instance available\n");
606 AliTRDcalibDB* calibration = AliTRDcalibDB::Instance();
608 AliError("No AliTRDcalibDB instance available\n");
612 Double_t *inADC = new Double_t[nTimeTotal]; // adc data before tail cancellation
613 Double_t *outADC = new Double_t[nTimeTotal]; // adc data after tail cancellation
615 AliDebug(1,Form("Tail cancellation (nExp = %d) for detector %d.\n"
616 ,recParam->GetTCnexp(),idet));
618 for (iRow = 0; iRow < nRowMax; iRow++ ) {
619 for (iCol = 0; iCol < nColMax; iCol++ ) {
620 for (iTime = 0; iTime < nTimeTotal; iTime++) {
625 Double_t gain = calibration->GetGainFactor(idet,iCol,iRow);
627 AliError("Not a valid gain\n");
629 inADC[iTime] = digitsIn->GetDataUnchecked(iRow,iCol,iTime);
630 inADC[iTime] /= gain;
631 outADC[iTime] = inADC[iTime];
635 // Apply the tail cancelation via the digital filter
636 if (recParam->TCOn()) {
637 DeConvExp(inADC,outADC,nTimeTotal,recParam->GetTCnexp());
640 for (iTime = 0; iTime < nTimeTotal; iTime++) {
642 // Store the amplitude of the digit if above threshold
643 if (outADC[iTime] > ADCthreshold) {
644 AliDebug(2,Form(" iRow = %d, iCol = %d, iTime = %d, adc = %f\n"
645 ,iRow,iCol,iTime,outADC[iTime]));
646 digitsOut->SetDataUnchecked(iRow,iCol,iTime,outADC[iTime]);
661 //_____________________________________________________________________________
662 void AliTRDclusterizerV1::DeConvExp(Double_t *source, Double_t *target,
666 // Tail cancellation by deconvolution for PASA v4 TRF
670 Double_t coefficients[2];
672 // Initialization (coefficient = alpha, rates = lambda)
678 if (nexp == 1) { // 1 Exponentials
684 if (nexp == 2) { // 2 Exponentials
691 coefficients[0] = C1;
692 coefficients[1] = C2;
696 rates[0] = TMath::Exp(-Dt/(R1));
697 rates[1] = TMath::Exp(-Dt/(R2));
702 Double_t reminder[2];
706 // Attention: computation order is important
708 for (k = 0; k < nexp; k++) {
711 for (i = 0; i < n; i++) {
712 result = (source[i] - correction); // no rescaling
715 for (k = 0; k < nexp; k++) {
716 reminder[k] = rates[k] * (reminder[k] + coefficients[k] * result);
719 for (k = 0; k < nexp; k++) {
720 correction += reminder[k];