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 //-------------------------------------------------------
19 // Implementation of the TPC clusterer
21 // Origin: Marian Ivanov
22 //-------------------------------------------------------
24 #include "AliTPCReconstructor.h"
25 #include "AliTPCclustererMI.h"
26 #include "AliTPCclusterMI.h"
27 #include <TObjArray.h>
32 #include "AliMathBase.h"
34 #include "AliTPCClustersArray.h"
35 #include "AliTPCClustersRow.h"
36 #include "AliDigits.h"
37 #include "AliSimDigits.h"
38 #include "AliTPCParam.h"
39 #include "AliTPCRecoParam.h"
40 #include "AliRawReader.h"
41 #include "AliTPCRawStream.h"
42 #include "AliRunLoader.h"
43 #include "AliLoader.h"
44 #include "Riostream.h"
46 #include "AliTPCcalibDB.h"
47 #include "AliTPCCalPad.h"
48 #include "AliTPCCalROC.h"
49 #include "TTreeStream.h"
53 ClassImp(AliTPCclustererMI)
57 AliTPCclustererMI::AliTPCclustererMI(const AliTPCParam* par, const AliTPCRecoParam * recoParam)
59 fIsOldRCUFormat = kFALSE;
64 fRecoParam = recoParam;
66 //set default parameters if not specified
67 fRecoParam = AliTPCReconstructor::GetRecoParam();
68 if (!fRecoParam) fRecoParam = AliTPCRecoParam::GetLowFluxParam();
70 fDebugStreamer = new TTreeSRedirector("TPCsignal.root");
74 AliTPCclustererMI::~AliTPCclustererMI(){
76 if (fAmplitudeHisto) delete fAmplitudeHisto;
77 if (fDebugStreamer) delete fDebugStreamer;
80 void AliTPCclustererMI::SetInput(TTree * tree)
83 // set input tree with digits
86 if (!fInput->GetBranch("Segment")){
87 cerr<<"AliTPC::Digits2Clusters(): no porper input tree !\n";
93 void AliTPCclustererMI::SetOutput(TTree * tree)
98 AliTPCClustersRow clrow;
99 AliTPCClustersRow *pclrow=&clrow;
100 clrow.SetClass("AliTPCclusterMI");
101 clrow.SetArray(1); // to make Clones array
102 fOutput->Branch("Segment","AliTPCClustersRow",&pclrow,32000,200);
106 Float_t AliTPCclustererMI::GetSigmaY2(Int_t iz){
107 // sigma y2 = in digits - we don't know the angle
108 Float_t z = iz*fParam->GetZWidth()+fParam->GetNTBinsL1()*fParam->GetZWidth();
109 Float_t sd2 = (z*fParam->GetDiffL()*fParam->GetDiffL())/
110 (fPadWidth*fPadWidth);
112 Float_t res = sd2+sres;
117 Float_t AliTPCclustererMI::GetSigmaZ2(Int_t iz){
118 //sigma z2 = in digits - angle estimated supposing vertex constraint
119 Float_t z = iz*fZWidth+fParam->GetNTBinsL1()*fParam->GetZWidth();
120 Float_t sd2 = (z*fParam->GetDiffL()*fParam->GetDiffL())/(fZWidth*fZWidth);
121 Float_t angular = fPadLength*(fParam->GetZLength()-z)/(fRx*fZWidth);
124 Float_t sres = fParam->GetZSigma()/fZWidth;
126 Float_t res = angular +sd2+sres;
130 void AliTPCclustererMI::MakeCluster(Int_t k,Int_t max,Float_t *bins, UInt_t /*m*/,
133 Int_t i0=k/max; //central pad
134 Int_t j0=k%max; //central time bin
136 // set pointers to data
137 //Int_t dummy[5] ={0,0,0,0,0};
138 Float_t * matrix[5]; //5x5 matrix with digits - indexing i = 0 ..4 j = -2..2
139 Float_t * resmatrix[5];
140 for (Int_t di=-2;di<=2;di++){
141 matrix[di+2] = &bins[k+di*max];
142 resmatrix[di+2] = &fResBins[k+di*max];
144 //build matrix with virtual charge
145 Float_t sigmay2= GetSigmaY2(j0);
146 Float_t sigmaz2= GetSigmaZ2(j0);
148 Float_t vmatrix[5][5];
149 vmatrix[2][2] = matrix[2][0];
151 c.SetMax((UShort_t)(vmatrix[2][2])); // write maximal amplitude
152 for (Int_t di =-1;di <=1;di++)
153 for (Int_t dj =-1;dj <=1;dj++){
154 Float_t amp = matrix[di+2][dj];
155 if ( (amp<2) && (fLoop<2)){
156 // if under threshold - calculate virtual charge
157 Float_t ratio = TMath::Exp(-1.2*TMath::Abs(di)/sigmay2)*TMath::Exp(-1.2*TMath::Abs(dj)/sigmaz2);
158 amp = ((matrix[2][0]-2)*(matrix[2][0]-2)/(matrix[-di+2][-dj]+2))*ratio;
160 vmatrix[2+di][2+dj]=amp;
161 vmatrix[2+2*di][2+2*dj]=0;
164 vmatrix[2+2*di][2+dj] =0;
165 vmatrix[2+di][2+2*dj] =0;
170 //if small amplitude - below 2 x threshold - don't consider other one
171 vmatrix[2+di][2+dj]=amp;
172 vmatrix[2+2*di][2+2*dj]=0; // don't take to the account next bin
175 vmatrix[2+2*di][2+dj] =0;
176 vmatrix[2+di][2+2*dj] =0;
180 //if bigger then take everything
181 vmatrix[2+di][2+dj]=amp;
182 vmatrix[2+2*di][2+2*dj]= matrix[2*di+2][2*dj] ;
185 vmatrix[2+2*di][2+dj] = matrix[2*di+2][dj];
186 vmatrix[2+di][2+2*dj] = matrix[2+di][dj*2];
198 for (Int_t i=-2;i<=2;i++)
199 for (Int_t j=-2;j<=2;j++){
200 Float_t amp = vmatrix[i+2][j+2];
209 Float_t meani = sumiw/sumw;
210 Float_t mi2 = sumi2w/sumw-meani*meani;
211 Float_t meanj = sumjw/sumw;
212 Float_t mj2 = sumj2w/sumw-meanj*meanj;
214 Float_t ry = mi2/sigmay2;
215 Float_t rz = mj2/sigmaz2;
218 if ( ( (ry<0.6) || (rz<0.6) ) && fLoop==2) return;
219 if ( (ry <1.2) && (rz<1.2) || (!fRecoParam->GetDoUnfold())) {
221 //if cluster looks like expected or Unfolding not switched on
222 //standard COG is used
223 //+1.2 deviation from expected sigma accepted
224 // c.fMax = FitMax(vmatrix,meani,meanj,TMath::Sqrt(sigmay2),TMath::Sqrt(sigmaz2));
228 //set cluster parameters
230 c.SetY(meani*fPadWidth);
231 c.SetZ(meanj*fZWidth);
235 //remove cluster data from data
236 for (Int_t di=-2;di<=2;di++)
237 for (Int_t dj=-2;dj<=2;dj++){
238 resmatrix[di+2][dj] -= vmatrix[di+2][dj+2];
239 if (resmatrix[di+2][dj]<0) resmatrix[di+2][dj]=0;
245 //unfolding when neccessary
248 Float_t * matrix2[7]; //7x7 matrix with digits - indexing i = 0 ..6 j = -3..3
249 Float_t dummy[7]={0,0,0,0,0,0};
250 for (Int_t di=-3;di<=3;di++){
251 matrix2[di+3] = &bins[k+di*max];
252 if ((k+di*max)<3) matrix2[di+3] = &dummy[3];
253 if ((k+di*max)>fMaxBin-3) matrix2[di+3] = &dummy[3];
255 Float_t vmatrix2[5][5];
258 UnfoldCluster(matrix2,vmatrix2,meani,meanj,sumu,overlap);
260 // c.fMax = FitMax(vmatrix2,meani,meanj,TMath::Sqrt(sigmay2),TMath::Sqrt(sigmaz2));
263 //set cluster parameters
265 c.SetY(meani*fPadWidth);
266 c.SetZ(meanj*fZWidth);
269 c.SetType(Char_t(overlap)+1);
276 printf("%f\t%f\n", vmatrix2[2][2], vmatrix[2][2]);
281 void AliTPCclustererMI::UnfoldCluster(Float_t * matrix2[7], Float_t recmatrix[5][5], Float_t & meani, Float_t & meanj,
282 Float_t & sumu, Float_t & overlap )
285 //unfold cluster from input matrix
286 //data corresponding to cluster writen in recmatrix
287 //output meani and meanj
289 //take separatelly y and z
291 Float_t sum3i[7] = {0,0,0,0,0,0,0};
292 Float_t sum3j[7] = {0,0,0,0,0,0,0};
294 for (Int_t k =0;k<7;k++)
295 for (Int_t l = -1; l<=1;l++){
296 sum3i[k]+=matrix2[k][l];
297 sum3j[k]+=matrix2[l+3][k-3];
299 Float_t mratio[3][3]={{1,1,1},{1,1,1},{1,1,1}};
302 Float_t sum3wi = 0; //charge minus overlap
303 Float_t sum3wio = 0; //full charge
304 Float_t sum3iw = 0; //sum for mean value
305 for (Int_t dk=-1;dk<=1;dk++){
306 sum3wio+=sum3i[dk+3];
312 if ( ( ((sum3i[dk+3]+3)/(sum3i[3]-3))+1 < (sum3i[2*dk+3]-3)/(sum3i[dk+3]+3))||
313 sum3i[dk+3]<=sum3i[2*dk+3] && sum3i[dk+3]>2 ){
314 Float_t xm2 = sum3i[-dk+3];
315 Float_t xm1 = sum3i[+3];
316 Float_t x1 = sum3i[2*dk+3];
317 Float_t x2 = sum3i[3*dk+3];
318 Float_t w11 = TMath::Max((Float_t)(4.*xm1-xm2),(Float_t)0.000001);
319 Float_t w12 = TMath::Max((Float_t)(4 *x1 -x2),(Float_t)0.);
320 ratio = w11/(w11+w12);
321 for (Int_t dl=-1;dl<=1;dl++)
322 mratio[dk+1][dl+1] *= ratio;
324 Float_t amp = sum3i[dk+3]*ratio;
329 meani = sum3iw/sum3wi;
330 Float_t overlapi = (sum3wio-sum3wi)/sum3wio;
335 Float_t sum3wj = 0; //charge minus overlap
336 Float_t sum3wjo = 0; //full charge
337 Float_t sum3jw = 0; //sum for mean value
338 for (Int_t dk=-1;dk<=1;dk++){
339 sum3wjo+=sum3j[dk+3];
345 if ( ( ((sum3j[dk+3]+3)/(sum3j[3]-3))+1 < (sum3j[2*dk+3]-3)/(sum3j[dk+3]+3)) ||
346 (sum3j[dk+3]<=sum3j[2*dk+3] && sum3j[dk+3]>2)){
347 Float_t xm2 = sum3j[-dk+3];
348 Float_t xm1 = sum3j[+3];
349 Float_t x1 = sum3j[2*dk+3];
350 Float_t x2 = sum3j[3*dk+3];
351 Float_t w11 = TMath::Max((Float_t)(4.*xm1-xm2),(Float_t)0.000001);
352 Float_t w12 = TMath::Max((Float_t)(4 *x1 -x2),(Float_t)0.);
353 ratio = w11/(w11+w12);
354 for (Int_t dl=-1;dl<=1;dl++)
355 mratio[dl+1][dk+1] *= ratio;
357 Float_t amp = sum3j[dk+3]*ratio;
362 meanj = sum3jw/sum3wj;
363 Float_t overlapj = (sum3wjo-sum3wj)/sum3wjo;
364 overlap = Int_t(100*TMath::Max(overlapi,overlapj)+3);
365 sumu = (sum3wj+sum3wi)/2.;
368 //if not overlap detected remove everything
369 for (Int_t di =-2; di<=2;di++)
370 for (Int_t dj =-2; dj<=2;dj++){
371 recmatrix[di+2][dj+2] = matrix2[3+di][dj];
375 for (Int_t di =-1; di<=1;di++)
376 for (Int_t dj =-1; dj<=1;dj++){
378 if (mratio[di+1][dj+1]==1){
379 recmatrix[di+2][dj+2] = matrix2[3+di][dj];
380 if (TMath::Abs(di)+TMath::Abs(dj)>1){
381 recmatrix[2*di+2][dj+2] = matrix2[3+2*di][dj];
382 recmatrix[di+2][2*dj+2] = matrix2[3+di][2*dj];
384 recmatrix[2*di+2][2*dj+2] = matrix2[3+2*di][2*dj];
388 //if we have overlap in direction
389 recmatrix[di+2][dj+2] = mratio[di+1][dj+1]* matrix2[3+di][dj];
390 if (TMath::Abs(di)+TMath::Abs(dj)>1){
391 ratio = TMath::Min((Float_t)(recmatrix[di+2][dj+2]/(matrix2[3+0*di][1*dj]+1)),(Float_t)1.);
392 recmatrix[2*di+2][dj+2] = ratio*recmatrix[di+2][dj+2];
394 ratio = TMath::Min((Float_t)(recmatrix[di+2][dj+2]/(matrix2[3+1*di][0*dj]+1)),(Float_t)1.);
395 recmatrix[di+2][2*dj+2] = ratio*recmatrix[di+2][dj+2];
398 ratio = recmatrix[di+2][dj+2]/matrix2[3][0];
399 recmatrix[2*di+2][2*dj+2] = ratio*recmatrix[di+2][dj+2];
405 printf("%f\n", recmatrix[2][2]);
409 Float_t AliTPCclustererMI::FitMax(Float_t vmatrix[5][5], Float_t y, Float_t z, Float_t sigmay, Float_t sigmaz)
416 for (Int_t di = -1;di<=1;di++)
417 for (Int_t dj = -1;dj<=1;dj++){
418 if (vmatrix[2+di][2+dj]>2){
419 Float_t teor = TMath::Gaus(di,y,sigmay*1.2)*TMath::Gaus(dj,z,sigmaz*1.2);
420 sumteor += teor*vmatrix[2+di][2+dj];
421 sumamp += vmatrix[2+di][2+dj]*vmatrix[2+di][2+dj];
424 Float_t max = sumamp/sumteor;
428 void AliTPCclustererMI::AddCluster(AliTPCclusterMI &c){
430 // transform cluster to the global coordinata
431 // add the cluster to the array
433 Float_t meani = c.GetY()/fPadWidth;
434 Float_t meanj = c.GetZ()/fZWidth;
436 Int_t ki = TMath::Nint(meani-3);
438 if (ki>=fMaxPad) ki = fMaxPad-1;
439 Int_t kj = TMath::Nint(meanj-3);
441 if (kj>=fMaxTime-3) kj=fMaxTime-4;
442 // ki and kj shifted to "real" coordinata
444 c.SetLabel(fRowDig->GetTrackIDFast(kj,ki,0)-2,0);
445 c.SetLabel(fRowDig->GetTrackIDFast(kj,ki,1)-2,1);
446 c.SetLabel(fRowDig->GetTrackIDFast(kj,ki,2)-2,2);
450 Float_t s2 = c.GetSigmaY2();
451 Float_t w=fParam->GetPadPitchWidth(fSector);
453 c.SetSigmaY2(s2*w*w);
456 c.SetSigmaZ2(s2*w*w);
457 c.SetY((meani - 2.5 - 0.5*fMaxPad)*fParam->GetPadPitchWidth(fSector));
458 c.SetZ(fZWidth*(meanj-3));
459 c.SetZ(c.GetZ() - 3.*fParam->GetZSigma() + fParam->GetNTBinsL1()*fParam->GetZWidth()); // PASA delay + L1 delay
460 c.SetZ(fSign*(fParam->GetZLength() - c.GetZ()));
462 c.SetDetector(fSector);
465 if (ki<=1 || ki>=fMaxPad-1 || kj==1 || kj==fMaxTime-2) {
466 //c.SetSigmaY2(c.GetSigmaY2()*25.);
467 //c.SetSigmaZ2(c.GetSigmaZ2()*4.);
468 c.SetType(-(c.GetType()+3)); //edge clusters
470 if (fLoop==2) c.SetType(100);
472 TClonesArray * arr = fRowCl->GetArray();
473 // AliTPCclusterMI * cl =
474 new ((*arr)[fNcluster]) AliTPCclusterMI(c);
480 //_____________________________________________________________________________
481 void AliTPCclustererMI::Digits2Clusters()
483 //-----------------------------------------------------------------
484 // This is a simple cluster finder.
485 //-----------------------------------------------------------------
488 Error("Digits2Clusters", "input tree not initialised");
493 Error("Digits2Clusters", "output tree not initialised");
497 AliTPCCalPad * gainTPC = AliTPCcalibDB::Instance()->GetPadGainFactor();
499 AliSimDigits digarr, *dummy=&digarr;
501 fInput->GetBranch("Segment")->SetAddress(&dummy);
502 Stat_t nentries = fInput->GetEntries();
504 fMaxTime=fParam->GetMaxTBin()+6; // add 3 virtual time bins before and 3 after
508 for (Int_t n=0; n<nentries; n++) {
510 if (!fParam->AdjustSectorRow(digarr.GetID(),fSector,fRow)) {
511 cerr<<"AliTPC warning: invalid segment ID ! "<<digarr.GetID()<<endl;
515 AliTPCCalROC * gainROC = gainTPC->GetCalROC(fSector); // pad gains per given sector
518 AliTPCClustersRow *clrow= new AliTPCClustersRow();
520 clrow->SetClass("AliTPCclusterMI");
523 clrow->SetID(digarr.GetID());
524 fOutput->GetBranch("Segment")->SetAddress(&clrow);
525 fRx=fParam->GetPadRowRadii(fSector,row);
528 const Int_t kNIS=fParam->GetNInnerSector(), kNOS=fParam->GetNOuterSector();
529 fZWidth = fParam->GetZWidth();
530 if (fSector < kNIS) {
531 fMaxPad = fParam->GetNPadsLow(row);
532 fSign = (fSector < kNIS/2) ? 1 : -1;
533 fPadLength = fParam->GetPadPitchLength(fSector,row);
534 fPadWidth = fParam->GetPadPitchWidth();
536 fMaxPad = fParam->GetNPadsUp(row);
537 fSign = ((fSector-kNIS) < kNOS/2) ? 1 : -1;
538 fPadLength = fParam->GetPadPitchLength(fSector,row);
539 fPadWidth = fParam->GetPadPitchWidth();
543 fMaxBin=fMaxTime*(fMaxPad+6); // add 3 virtual pads before and 3 after
544 fBins =new Float_t[fMaxBin];
545 fResBins =new Float_t[fMaxBin]; //fBins with residuals after 1 finder loop
546 memset(fBins,0,sizeof(Float_t)*fMaxBin);
547 memset(fResBins,0,sizeof(Float_t)*fMaxBin);
549 if (digarr.First()) //MI change
551 Float_t dig=digarr.CurrentDigit();
552 if (dig<=fParam->GetZeroSup()) continue;
553 Int_t j=digarr.CurrentRow()+3, i=digarr.CurrentColumn()+3;
554 Float_t gain = gainROC->GetValue(row,digarr.CurrentColumn());
555 fBins[i*fMaxTime+j]=dig/gain;
556 } while (digarr.Next());
557 digarr.ExpandTrackBuffer();
563 nclusters+=fNcluster;
568 Info("Digits2Clusters", "Number of found clusters : %d", nclusters);
571 void AliTPCclustererMI::Digits2Clusters(AliRawReader* rawReader)
573 //-----------------------------------------------------------------
574 // This is a cluster finder for the TPC raw data.
575 // The method assumes NO ordering of the altro channels.
576 // The pedestal subtraction can be switched on and off
577 // using an option of the TPC reconstructor
578 //-----------------------------------------------------------------
581 Error("Digits2Clusters", "output tree not initialised");
587 AliTPCCalPad * gainTPC = AliTPCcalibDB::Instance()->GetPadGainFactor();
591 fMaxTime = fParam->GetMaxTBin() + 6; // add 3 virtual time bins before and 3 after
592 const Int_t kNIS = fParam->GetNInnerSector();
593 const Int_t kNOS = fParam->GetNOuterSector();
594 const Int_t kNS = kNIS + kNOS;
595 fZWidth = fParam->GetZWidth();
596 Int_t zeroSup = fParam->GetZeroSup();
598 Float_t** allBins = NULL;
599 Float_t** allBinsRes = NULL;
602 for(fSector = 0; fSector < kNS; fSector++) {
604 AliTPCCalROC * gainROC = gainTPC->GetCalROC(fSector); // pad gains per given sector
607 Int_t nDDLs = 0, indexDDL = 0;
608 if (fSector < kNIS) {
609 nRows = fParam->GetNRowLow();
610 fSign = (fSector < kNIS/2) ? 1 : -1;
612 indexDDL = fSector * 2;
615 nRows = fParam->GetNRowUp();
616 fSign = ((fSector-kNIS) < kNOS/2) ? 1 : -1;
618 indexDDL = (fSector-kNIS) * 4 + kNIS * 2;
621 allBins = new Float_t*[nRows];
622 allBinsRes = new Float_t*[nRows];
624 for (Int_t iRow = 0; iRow < nRows; iRow++) {
627 maxPad = fParam->GetNPadsLow(iRow);
629 maxPad = fParam->GetNPadsUp(iRow);
631 Int_t maxBin = fMaxTime*(maxPad+6); // add 3 virtual pads before and 3 after
632 allBins[iRow] = new Float_t[maxBin];
633 allBinsRes[iRow] = new Float_t[maxBin];
634 memset(allBins[iRow],0,sizeof(Float_t)*maxBin);
635 memset(allBinsRes[iRow],0,sizeof(Float_t)*maxBin);
638 // Loas the raw data for corresponding DDLs
640 AliTPCRawStream input(rawReader);
641 input.SetOldRCUFormat(fIsOldRCUFormat);
642 rawReader->Select("TPC",indexDDL,indexDDL+nDDLs-1);
644 // Begin loop over altro data
645 while (input.Next()) {
647 if (input.GetSector() != fSector)
648 AliFatal(Form("Sector index mismatch ! Expected (%d), but got (%d) !",fSector,input.GetSector()));
650 if (input.GetTime() < AliTPCReconstructor::GetRecoParam()->GetFirstBin()) continue;
651 if (input.GetTime() > AliTPCReconstructor::GetRecoParam()->GetLastBin()) continue;
653 Int_t iRow = input.GetRow();
654 if (iRow < 0 || iRow >= nRows)
655 AliFatal(Form("Pad-row index (%d) outside the range (%d -> %d) !",
658 Int_t iPad = input.GetPad() + 3;
662 maxPad = fParam->GetNPadsLow(iRow);
664 maxPad = fParam->GetNPadsUp(iRow);
666 if (input.GetPad() < 0 || input.GetPad() >= maxPad)
667 AliFatal(Form("Pad index (%d) outside the range (%d -> %d) !",
668 input.GetPad(), 0, maxPad -1));
670 Int_t iTimeBin = input.GetTime() + 3;
671 if ( input.GetTime() < 0 || input.GetTime() >= fParam->GetMaxTBin())
672 AliFatal(Form("Timebin index (%d) outside the range (%d -> %d) !",
673 input.GetTime(), 0, fParam->GetMaxTBin() -1));
675 Int_t maxBin = fMaxTime*(maxPad+6); // add 3 virtual pads before and 3 after
677 if (((iPad*fMaxTime+iTimeBin) >= maxBin) ||
678 ((iPad*fMaxTime+iTimeBin) < 0))
679 AliFatal(Form("Index outside the allowed range"
680 " Sector=%d Row=%d Pad=%d Timebin=%d"
681 " (Max.index=%d)",fSector,iRow,iPad,iTimeBin,maxBin));
683 Float_t signal = input.GetSignal();
684 // if (!fPedSubtraction && signal <= zeroSup) continue;
685 if (!fRecoParam->GetCalcPedestal() && signal <= zeroSup) continue;
687 Float_t gain = gainROC->GetValue(iRow,input.GetPad());
688 allBins[iRow][iPad*fMaxTime+iTimeBin] = signal/gain;
689 allBins[iRow][iPad*fMaxTime+0] = 1; // pad with signal
690 } // End of the loop over altro data
692 // Now loop over rows and perform pedestal subtraction
693 // if (fPedSubtraction) {
694 if (fRecoParam->GetCalcPedestal()) {
695 for (Int_t iRow = 0; iRow < nRows; iRow++) {
698 maxPad = fParam->GetNPadsLow(iRow);
700 maxPad = fParam->GetNPadsUp(iRow);
702 for (Int_t iPad = 0; iPad < maxPad + 6; iPad++) {
703 if (allBins[iRow][iPad*fMaxTime+0] !=1) continue; // no data
704 Float_t *p = &allBins[iRow][iPad*fMaxTime+3];
705 //Float_t pedestal = TMath::Median(fMaxTime, p);
706 Int_t id[3] = {fSector, iRow, iPad-3};
707 Float_t pedestal = ProcesSignal(p, fMaxTime, id);
708 for (Int_t iTimeBin = 0; iTimeBin < fMaxTime; iTimeBin++) {
709 allBins[iRow][iPad*fMaxTime+iTimeBin] -= pedestal;
710 if (allBins[iRow][iPad*fMaxTime+iTimeBin] < zeroSup)
711 allBins[iRow][iPad*fMaxTime+iTimeBin] = 0;
717 // Now loop over rows and find clusters
718 for (fRow = 0; fRow < nRows; fRow++) {
719 fRowCl = new AliTPCClustersRow;
720 fRowCl->SetClass("AliTPCclusterMI");
722 fRowCl->SetID(fParam->GetIndex(fSector, fRow));
723 fOutput->GetBranch("Segment")->SetAddress(&fRowCl);
725 fRx = fParam->GetPadRowRadii(fSector, fRow);
726 fPadLength = fParam->GetPadPitchLength(fSector, fRow);
727 fPadWidth = fParam->GetPadPitchWidth();
729 fMaxPad = fParam->GetNPadsLow(fRow);
731 fMaxPad = fParam->GetNPadsUp(fRow);
732 fMaxBin = fMaxTime*(fMaxPad+6); // add 3 virtual pads before and 3 after
734 fBins = allBins[fRow];
735 fResBins = allBinsRes[fRow];
741 nclusters += fNcluster;
742 } // End of loop to find clusters
745 for (Int_t iRow = 0; iRow < nRows; iRow++) {
746 delete [] allBins[iRow];
747 delete [] allBinsRes[iRow];
751 delete [] allBinsRes;
753 } // End of loop over sectors
755 Info("Digits2Clusters", "Number of found clusters : %d\n", nclusters);
759 void AliTPCclustererMI::FindClusters()
761 //add virtual charge at the edge
762 for (Int_t i=0; i<fMaxTime; i++){
763 Float_t amp1 = fBins[i+3*fMaxTime];
766 Float_t amp2 = fBins[i+4*fMaxTime];
767 if (amp2==0) amp2=0.5;
768 Float_t sigma2 = GetSigmaY2(i);
769 amp0 = (amp1*amp1/amp2)*TMath::Exp(-1./sigma2);
770 if (gDebug>4) printf("\n%f\n",amp0);
772 fBins[i+2*fMaxTime] = amp0;
774 amp1 = fBins[(fMaxPad+2)*fMaxTime+i];
776 Float_t amp2 = fBins[i+(fMaxPad+1)*fMaxTime];
777 if (amp2==0) amp2=0.5;
778 Float_t sigma2 = GetSigmaY2(i);
779 amp0 = (amp1*amp1/amp2)*TMath::Exp(-1./sigma2);
780 if (gDebug>4) printf("\n%f\n",amp0);
782 fBins[(fMaxPad+3)*fMaxTime+i] = amp0;
785 // memcpy(fResBins,fBins, fMaxBin*2);
786 memcpy(fResBins,fBins, fMaxBin);
789 //first loop - for "gold cluster"
791 Float_t *b=&fBins[-1]+2*fMaxTime;
792 Int_t crtime = Int_t((fParam->GetZLength()-fRecoParam->GetCtgRange()*fRx)/fZWidth-fParam->GetNTBinsL1()-5);
794 for (Int_t i=2*fMaxTime; i<fMaxBin-2*fMaxTime; i++) {
796 if (*b<8) continue; //threshold form maxima
797 if (i%fMaxTime<crtime) {
798 Int_t delta = -(i%fMaxTime)+crtime;
804 if (!IsMaximum(*b,fMaxTime,b)) continue;
807 MakeCluster(i, fMaxTime, fBins, dummy,c);
810 //memcpy(fBins,fResBins, fMaxBin*2);
811 //second loop - for rest cluster
814 b=&fResBins[-1]+2*fMaxTime;
815 for (Int_t i=2*fMaxTime; i<fMaxBin-2*fMaxTime; i++) {
817 if (*b<25) continue; // bigger threshold for maxima
818 if (!IsMaximum(*b,fMaxTime,b)) continue;
821 MakeCluster(i, fMaxTime, fResBins, dummy,c);
828 Double_t AliTPCclustererMI::ProcesSignal(Float_t *signal, Int_t nchannels, Int_t id[3]){
830 // process signal on given pad - + streaming of additional information in special mode
836 Float_t kMin =fRecoParam->GetDumpAmplitudeMin(); // minimal signal to be dumped
837 Double_t median = TMath::Median(nchannels-offset, &(signal[offset]));
838 if (AliLog::GetDebugLevel("","AliTPCclustererMI")==0) return median;
841 Double_t mean = TMath::Mean(nchannels-offset, &(signal[offset]));
842 Double_t rms = TMath::RMS(nchannels-offset, &(signal[offset]));
843 Double_t *dsignal = new Double_t[nchannels];
844 Double_t *dtime = new Double_t[nchannels];
847 for (Int_t i=0; i<nchannels; i++){
849 dsignal[i] = signal[i];
850 if (signal[i]>max && i <fMaxTime-100) { // temporary remove spike signals at the end
856 Double_t mean06=0, mean09=0;
857 Double_t rms06=0, rms09=0;
858 AliMathBase::EvaluateUni(nchannels-offset, &(dsignal[offset]), mean06, rms06, int(0.6*(nchannels-offset)));
859 AliMathBase::EvaluateUni(nchannels-offset, &(dsignal[offset]), mean09, rms09, int(0.9*(nchannels-offset)));
862 UInt_t uid[3] = {UInt_t(id[0]),UInt_t(id[1]),UInt_t(id[2])};
863 if (uid[0]< AliTPCROC::Instance()->GetNSectors()
864 && uid[1]< AliTPCROC::Instance()->GetNRows(uid[0]) &&
865 uid[2] < AliTPCROC::Instance()->GetNPads(uid[0], uid[1])){
866 if (!fAmplitudeHisto){
867 fAmplitudeHisto = new TObjArray(72);
869 TObjArray * sectorArray = (TObjArray*)fAmplitudeHisto->UncheckedAt(uid[0]);
871 Int_t npads = AliTPCROC::Instance()->GetNChannels(uid[0]);
872 sectorArray = new TObjArray(npads);
873 fAmplitudeHisto->AddAt(sectorArray, uid[0]);
875 Int_t position = uid[2]+ AliTPCROC::Instance()->GetRowIndexes(uid[0])[uid[1]];
876 TH1F * histo = (TH1F*)sectorArray->UncheckedAt(position);
879 sprintf(hname,"Amp_%d_%d_%d",uid[0],uid[1],uid[2]);
880 TFile * backup = gFile;
881 fDebugStreamer->GetFile()->cd();
882 histo = new TH1F(hname, hname, 100, 5,100);
883 //histo->SetDirectory(0); // histogram not connected to directory -(File)
884 sectorArray->AddAt(histo, position);
885 if (backup) backup->cd();
887 for (Int_t i=0; i<nchannels; i++){
888 if (signal[i]>0) histo->Fill(signal[i]);
893 Bool_t random = (gRandom->Rndm()<0.0001);
894 if (max-median>kMin || rms06>2.*fParam->GetZeroSup() || random){
895 graph =new TGraph(nchannels, dtime, dsignal);
896 if (rms06>2.*fParam->GetZeroSup() || random)
897 (*fDebugStreamer)<<"SignalN"<< //noise pads - or random sample of pads
914 (*fDebugStreamer)<<"SignalB"<< // pads with signal
933 (*fDebugStreamer)<<"Signal"<<
950 if (rms06>fRecoParam->GetMaxNoise()) return 1024+median; // sign noisy channel in debug mode
956 void AliTPCclustererMI::DumpHistos(){
957 if (!fAmplitudeHisto) return;
958 for (UInt_t isector=0; isector<AliTPCROC::Instance()->GetNSectors(); isector++){
959 TObjArray * array = (TObjArray*)fAmplitudeHisto->UncheckedAt(isector);
960 if (!array) continue;
961 for (UInt_t ipad = 0; ipad <(UInt_t)array->GetEntriesFast(); ipad++){
962 TH1F * histo = (TH1F*) array->UncheckedAt(ipad);
963 if (!histo) continue;
964 if (histo->GetEntries()<100) continue;
965 histo->Fit("gaus","q");
966 Float_t mean = histo->GetMean();
967 Float_t rms = histo->GetRMS();
968 Float_t gmean = histo->GetFunction("gaus")->GetParameter(1);
969 Float_t gsigma = histo->GetFunction("gaus")->GetParameter(2);
970 Float_t max = histo->GetFunction("gaus")->GetParameter(0);
973 UInt_t row=0, pad =0;
974 const UInt_t *indexes = AliTPCROC::Instance()->GetRowIndexes(isector);
975 for (UInt_t irow=0; irow< AliTPCROC::Instance()->GetNRows(isector); irow++){
976 if (indexes[irow]<ipad){
978 pad = ipad-indexes[irow];
981 Int_t rpad = pad - (AliTPCROC::Instance()->GetNPads(isector,row))/2;
983 (*fDebugStreamer)<<"Fit"<<
994 if (array->UncheckedAt(ipad)) fDebugStreamer->StoreObject(array->UncheckedAt(ipad));