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 // 1. The Input data for reconstruction - Options
22 // 1.a Simulated data - TTree - invoked Digits2Clusters()
23 // 1.b Raw data - Digits2Clusters(AliRawReader* rawReader);
24 // 1.c HLT clusters - Digits2Clusters and Digits2Clusters(AliRawReader* rawReader)
25 // invoke ReadHLTClusters()
27 // fUseHLTClusters - switches between different inputs
28 // 1 -> only TPC raw/sim data
29 // 2 -> if present TPC raw/sim data, otherwise HLT clusters
30 // 3 -> only HLT clusters
31 // 4 -> if present HLT clusters, otherwise TPC raw/sim data
34 // 2.a TTree with clusters - if SetOutput(TTree * tree) invoked
35 // 2.b TObjArray - Faster option for HLT
36 // 2.c TClonesArray - Faster option for HLT (smaller memory consumption), activate with fBClonesArray flag
38 // 3. Reconstruction setup
39 // see AliTPCRecoParam for list of parameters
40 // The reconstruction parameterization taken from the
41 // AliTPCReconstructor::GetRecoParam()
42 // Possible to setup it in reconstruction macro AliTPCReconstructor::SetRecoParam(...)
46 // Origin: Marian Ivanov
47 //-------------------------------------------------------
49 #include "Riostream.h"
54 #include <TObjArray.h>
55 #include <TClonesArray.h>
58 #include <TTreeStream.h>
62 #include "AliDigits.h"
63 #include "AliLoader.h"
65 #include "AliMathBase.h"
66 #include "AliRawEventHeaderBase.h"
67 #include "AliRawReader.h"
68 #include "AliRunLoader.h"
69 #include "AliSimDigits.h"
70 #include "AliTPCCalPad.h"
71 #include "AliTPCCalROC.h"
72 #include "AliTPCClustersArray.h"
73 #include "AliTPCClustersRow.h"
74 #include "AliTPCParam.h"
75 #include "AliTPCRawStream.h"
76 #include "AliTPCRawStreamV3.h"
77 #include "AliTPCRecoParam.h"
78 #include "AliTPCReconstructor.h"
79 #include "AliTPCcalibDB.h"
80 #include "AliTPCclusterInfo.h"
81 #include "AliTPCclusterMI.h"
82 #include "AliTPCTransform.h"
83 #include "AliTPCclustererMI.h"
85 ClassImp(AliTPCclustererMI)
89 AliTPCclustererMI::AliTPCclustererMI(const AliTPCParam* par, const AliTPCRecoParam * recoParam):
95 fMaxTime(1006), // 1000>940 so use 1000, add 3 virtual time bins before and 3 after
104 fPedSubtraction(kFALSE),
111 fOutputClonesArray(0),
119 fBDumpSignal(kFALSE),
120 fBClonesArray(kFALSE),
125 fHLTClusterAccess(NULL)
129 // param - tpc parameters for given file
130 // recoparam - reconstruction parameters
135 fRecoParam = recoParam;
137 //set default parameters if not specified
138 fRecoParam = AliTPCReconstructor::GetRecoParam();
139 if (!fRecoParam) fRecoParam = AliTPCRecoParam::GetLowFluxParam();
142 if(AliTPCReconstructor::StreamLevel()>0) {
143 fDebugStreamer = new TTreeSRedirector("TPCsignal.root");
146 // Int_t nPoints = fRecoParam->GetLastBin()-fRecoParam->GetFirstBin();
147 fRowCl= new AliTPCClustersRow("AliTPCclusterMI");
149 // Non-persistent arrays
151 //alocate memory for sector - maximal case
153 AliTPCROC * roc = AliTPCROC::Instance();
154 Int_t nRowsMax = roc->GetNRows(roc->GetNSector()-1);
155 Int_t nPadsMax = roc->GetNPads(roc->GetNSector()-1,nRowsMax-1);
157 fAllBins = new Float_t*[nRowsMax];
158 fAllSigBins = new Int_t*[nRowsMax];
159 fAllNSigBins = new Int_t[nRowsMax];
160 for (Int_t iRow = 0; iRow < nRowsMax; iRow++) {
162 Int_t maxBin = fMaxTime*(nPadsMax+6); // add 3 virtual pads before and 3 after
163 fAllBins[iRow] = new Float_t[maxBin];
164 memset(fAllBins[iRow],0,sizeof(Float_t)*maxBin);
165 fAllSigBins[iRow] = new Int_t[maxBin];
166 fAllNSigBins[iRow]=0;
170 //______________________________________________________________
171 AliTPCclustererMI::~AliTPCclustererMI(){
175 if (fDebugStreamer) delete fDebugStreamer;
177 //fOutputArray->Delete();
180 if (fOutputClonesArray){
181 fOutputClonesArray->Delete();
182 delete fOutputClonesArray;
186 fRowCl->GetArray()->Delete();
190 AliTPCROC * roc = AliTPCROC::Instance();
191 Int_t nRowsMax = roc->GetNRows(roc->GetNSector()-1);
192 for (Int_t iRow = 0; iRow < nRowsMax; iRow++) {
193 delete [] fAllBins[iRow];
194 delete [] fAllSigBins[iRow];
197 delete [] fAllSigBins;
198 delete [] fAllNSigBins;
199 if (fHLTClusterAccess) delete fHLTClusterAccess;
202 void AliTPCclustererMI::SetInput(TTree * tree)
205 // set input tree with digits
208 if (!fInput->GetBranch("Segment")){
209 cerr<<"AliTPC::Digits2Clusters(): no porper input tree !\n";
215 void AliTPCclustererMI::SetOutput(TTree * tree)
218 // Set the output tree
219 // If not set the ObjArray used - Option for HLT
223 AliTPCClustersRow clrow("AliTPCclusterMI");
224 AliTPCClustersRow *pclrow=&clrow;
225 fOutput->Branch("Segment","AliTPCClustersRow",&pclrow,32000,200);
229 void AliTPCclustererMI::FillRow(){
231 // fill the output container -
232 // 2 Options possible
236 if (fOutput) fOutput->Fill();
237 if (!fOutput && !fBClonesArray){
239 if (!fOutputArray) fOutputArray = new TObjArray(fParam->GetNRowsTotal());
240 if (fRowCl && fRowCl->GetArray()->GetEntriesFast()>0) fOutputArray->AddAt(fRowCl->Clone(), fRowCl->GetID());
244 Float_t AliTPCclustererMI::GetSigmaY2(Int_t iz){
245 // sigma y2 = in digits - we don't know the angle
246 Float_t z = iz*fParam->GetZWidth()+fParam->GetNTBinsL1()*fParam->GetZWidth();
247 Float_t sd2 = (z*fParam->GetDiffL()*fParam->GetDiffL())/
248 (fPadWidth*fPadWidth);
250 Float_t res = sd2+sres;
255 Float_t AliTPCclustererMI::GetSigmaZ2(Int_t iz){
256 //sigma z2 = in digits - angle estimated supposing vertex constraint
257 Float_t z = iz*fZWidth+fParam->GetNTBinsL1()*fParam->GetZWidth();
258 Float_t sd2 = (z*fParam->GetDiffL()*fParam->GetDiffL())/(fZWidth*fZWidth);
259 Float_t angular = fPadLength*(fParam->GetZLength(fSector)-z)/(fRx*fZWidth);
262 Float_t sres = fParam->GetZSigma()/fZWidth;
264 Float_t res = angular +sd2+sres;
268 void AliTPCclustererMI::MakeCluster(Int_t k,Int_t max,Float_t *bins, UInt_t /*m*/,
272 // k - Make cluster at position k
273 // bins - 2 D array of signals mapped to 1 dimensional array -
274 // max - the number of time bins er one dimension
275 // c - refernce to cluster to be filled
277 Int_t i0=k/max; //central pad
278 Int_t j0=k%max; //central time bin
280 // set pointers to data
281 //Int_t dummy[5] ={0,0,0,0,0};
282 Float_t * matrix[5]; //5x5 matrix with digits - indexing i = 0 ..4 j = -2..2
283 for (Int_t di=-2;di<=2;di++){
284 matrix[di+2] = &bins[k+di*max];
286 //build matrix with virtual charge
287 Float_t sigmay2= GetSigmaY2(j0);
288 Float_t sigmaz2= GetSigmaZ2(j0);
290 Float_t vmatrix[5][5];
291 vmatrix[2][2] = matrix[2][0];
293 c.SetMax((UShort_t)(vmatrix[2][2])); // write maximal amplitude
294 for (Int_t di =-1;di <=1;di++)
295 for (Int_t dj =-1;dj <=1;dj++){
296 Float_t amp = matrix[di+2][dj];
297 if ( (amp<2) && (fLoop<2)){
298 // if under threshold - calculate virtual charge
299 Float_t ratio = TMath::Exp(-1.2*TMath::Abs(di)/sigmay2)*TMath::Exp(-1.2*TMath::Abs(dj)/sigmaz2);
300 amp = ((matrix[2][0]-2)*(matrix[2][0]-2)/(matrix[-di+2][-dj]+2))*ratio;
302 vmatrix[2+di][2+dj]=amp;
303 vmatrix[2+2*di][2+2*dj]=0;
306 vmatrix[2+2*di][2+dj] =0;
307 vmatrix[2+di][2+2*dj] =0;
312 //if small amplitude - below 2 x threshold - don't consider other one
313 vmatrix[2+di][2+dj]=amp;
314 vmatrix[2+2*di][2+2*dj]=0; // don't take to the account next bin
317 vmatrix[2+2*di][2+dj] =0;
318 vmatrix[2+di][2+2*dj] =0;
322 //if bigger then take everything
323 vmatrix[2+di][2+dj]=amp;
324 vmatrix[2+2*di][2+2*dj]= matrix[2*di+2][2*dj] ;
327 vmatrix[2+2*di][2+dj] = matrix[2*di+2][dj];
328 vmatrix[2+di][2+2*dj] = matrix[2+di][dj*2];
340 for (Int_t i=-2;i<=2;i++)
341 for (Int_t j=-2;j<=2;j++){
342 Float_t amp = vmatrix[i+2][j+2];
351 Float_t meani = sumiw/sumw;
352 Float_t mi2 = sumi2w/sumw-meani*meani;
353 Float_t meanj = sumjw/sumw;
354 Float_t mj2 = sumj2w/sumw-meanj*meanj;
356 Float_t ry = mi2/sigmay2;
357 Float_t rz = mj2/sigmaz2;
360 if ( ( (ry<0.6) || (rz<0.6) ) && fLoop==2) return;
361 if ( ((ry <1.2) && (rz<1.2)) || (!fRecoParam->GetDoUnfold())) {
363 //if cluster looks like expected or Unfolding not switched on
364 //standard COG is used
365 //+1.2 deviation from expected sigma accepted
366 // c.fMax = FitMax(vmatrix,meani,meanj,TMath::Sqrt(sigmay2),TMath::Sqrt(sigmaz2));
370 //set cluster parameters
373 c.SetTimeBin(meanj-3);
377 AddCluster(c,(Float_t*)vmatrix,k);
381 //unfolding when neccessary
384 Float_t * matrix2[7]; //7x7 matrix with digits - indexing i = 0 ..6 j = -3..3
385 Float_t dummy[7]={0,0,0,0,0,0};
386 for (Int_t di=-3;di<=3;di++){
387 matrix2[di+3] = &bins[k+di*max];
388 if ((k+di*max)<3) matrix2[di+3] = &dummy[3];
389 if ((k+di*max)>fMaxBin-3) matrix2[di+3] = &dummy[3];
391 Float_t vmatrix2[5][5];
394 UnfoldCluster(matrix2,vmatrix2,meani,meanj,sumu,overlap);
396 // c.fMax = FitMax(vmatrix2,meani,meanj,TMath::Sqrt(sigmay2),TMath::Sqrt(sigmaz2));
399 //set cluster parameters
402 c.SetTimeBin(meanj-3);
405 c.SetType(Char_t(overlap)+1);
406 AddCluster(c,(Float_t*)vmatrix,k);
415 void AliTPCclustererMI::UnfoldCluster(Float_t * matrix2[7], Float_t recmatrix[5][5], Float_t & meani, Float_t & meanj,
416 Float_t & sumu, Float_t & overlap )
419 //unfold cluster from input matrix
420 //data corresponding to cluster writen in recmatrix
421 //output meani and meanj
423 //take separatelly y and z
425 Float_t sum3i[7] = {0,0,0,0,0,0,0};
426 Float_t sum3j[7] = {0,0,0,0,0,0,0};
428 for (Int_t k =0;k<7;k++)
429 for (Int_t l = -1; l<=1;l++){
430 sum3i[k]+=matrix2[k][l];
431 sum3j[k]+=matrix2[l+3][k-3];
433 Float_t mratio[3][3]={{1,1,1},{1,1,1},{1,1,1}};
436 Float_t sum3wi = 0; //charge minus overlap
437 Float_t sum3wio = 0; //full charge
438 Float_t sum3iw = 0; //sum for mean value
439 for (Int_t dk=-1;dk<=1;dk++){
440 sum3wio+=sum3i[dk+3];
446 if ( ( ((sum3i[dk+3]+3)/(sum3i[3]-3))+1 < (sum3i[2*dk+3]-3)/(sum3i[dk+3]+3))||
447 (sum3i[dk+3]<=sum3i[2*dk+3] && sum3i[dk+3]>2 )){
448 Float_t xm2 = sum3i[-dk+3];
449 Float_t xm1 = sum3i[+3];
450 Float_t x1 = sum3i[2*dk+3];
451 Float_t x2 = sum3i[3*dk+3];
452 Float_t w11 = TMath::Max((Float_t)(4.*xm1-xm2),(Float_t)0.000001);
453 Float_t w12 = TMath::Max((Float_t)(4 *x1 -x2),(Float_t)0.);
454 ratio = w11/(w11+w12);
455 for (Int_t dl=-1;dl<=1;dl++)
456 mratio[dk+1][dl+1] *= ratio;
458 Float_t amp = sum3i[dk+3]*ratio;
463 meani = sum3iw/sum3wi;
464 Float_t overlapi = (sum3wio-sum3wi)/sum3wio;
469 Float_t sum3wj = 0; //charge minus overlap
470 Float_t sum3wjo = 0; //full charge
471 Float_t sum3jw = 0; //sum for mean value
472 for (Int_t dk=-1;dk<=1;dk++){
473 sum3wjo+=sum3j[dk+3];
479 if ( ( ((sum3j[dk+3]+3)/(sum3j[3]-3))+1 < (sum3j[2*dk+3]-3)/(sum3j[dk+3]+3)) ||
480 (sum3j[dk+3]<=sum3j[2*dk+3] && sum3j[dk+3]>2)){
481 Float_t xm2 = sum3j[-dk+3];
482 Float_t xm1 = sum3j[+3];
483 Float_t x1 = sum3j[2*dk+3];
484 Float_t x2 = sum3j[3*dk+3];
485 Float_t w11 = TMath::Max((Float_t)(4.*xm1-xm2),(Float_t)0.000001);
486 Float_t w12 = TMath::Max((Float_t)(4 *x1 -x2),(Float_t)0.);
487 ratio = w11/(w11+w12);
488 for (Int_t dl=-1;dl<=1;dl++)
489 mratio[dl+1][dk+1] *= ratio;
491 Float_t amp = sum3j[dk+3]*ratio;
496 meanj = sum3jw/sum3wj;
497 Float_t overlapj = (sum3wjo-sum3wj)/sum3wjo;
498 overlap = Int_t(100*TMath::Max(overlapi,overlapj)+3);
499 sumu = (sum3wj+sum3wi)/2.;
502 //if not overlap detected remove everything
503 for (Int_t di =-2; di<=2;di++)
504 for (Int_t dj =-2; dj<=2;dj++){
505 recmatrix[di+2][dj+2] = matrix2[3+di][dj];
509 for (Int_t di =-1; di<=1;di++)
510 for (Int_t dj =-1; dj<=1;dj++){
512 if (mratio[di+1][dj+1]==1){
513 recmatrix[di+2][dj+2] = matrix2[3+di][dj];
514 if (TMath::Abs(di)+TMath::Abs(dj)>1){
515 recmatrix[2*di+2][dj+2] = matrix2[3+2*di][dj];
516 recmatrix[di+2][2*dj+2] = matrix2[3+di][2*dj];
518 recmatrix[2*di+2][2*dj+2] = matrix2[3+2*di][2*dj];
522 //if we have overlap in direction
523 recmatrix[di+2][dj+2] = mratio[di+1][dj+1]* matrix2[3+di][dj];
524 if (TMath::Abs(di)+TMath::Abs(dj)>1){
525 ratio = TMath::Min((Float_t)(recmatrix[di+2][dj+2]/(matrix2[3+0*di][1*dj]+1)),(Float_t)1.);
526 recmatrix[2*di+2][dj+2] = ratio*recmatrix[di+2][dj+2];
528 ratio = TMath::Min((Float_t)(recmatrix[di+2][dj+2]/(matrix2[3+1*di][0*dj]+1)),(Float_t)1.);
529 recmatrix[di+2][2*dj+2] = ratio*recmatrix[di+2][dj+2];
532 ratio = recmatrix[di+2][dj+2]/matrix2[3][0];
533 recmatrix[2*di+2][2*dj+2] = ratio*recmatrix[di+2][dj+2];
541 Float_t AliTPCclustererMI::FitMax(Float_t vmatrix[5][5], Float_t y, Float_t z, Float_t sigmay, Float_t sigmaz)
548 for (Int_t di = -1;di<=1;di++)
549 for (Int_t dj = -1;dj<=1;dj++){
550 if (vmatrix[2+di][2+dj]>2){
551 Float_t teor = TMath::Gaus(di,y,sigmay*1.2)*TMath::Gaus(dj,z,sigmaz*1.2);
552 sumteor += teor*vmatrix[2+di][2+dj];
553 sumamp += vmatrix[2+di][2+dj]*vmatrix[2+di][2+dj];
556 Float_t max = sumamp/sumteor;
560 void AliTPCclustererMI::AddCluster(AliTPCclusterMI &c, Float_t * /*matrix*/, Int_t /*pos*/){
563 // Transform cluster to the rotated global coordinata
564 // Assign labels to the cluster
565 // add the cluster to the array
566 // for more details - See AliTPCTranform::Transform(x,i,0,1)
567 Float_t meani = c.GetPad();
568 Float_t meanj = c.GetTimeBin();
570 Int_t ki = TMath::Nint(meani);
572 if (ki>=fMaxPad) ki = fMaxPad-1;
573 Int_t kj = TMath::Nint(meanj);
575 if (kj>=fMaxTime-3) kj=fMaxTime-4;
576 // ki and kj shifted as integers coordinata
578 c.SetLabel(fRowDig->GetTrackIDFast(kj,ki,0)-2,0);
579 c.SetLabel(fRowDig->GetTrackIDFast(kj,ki,1)-2,1);
580 c.SetLabel(fRowDig->GetTrackIDFast(kj,ki,2)-2,2);
584 c.SetDetector(fSector);
585 Float_t s2 = c.GetSigmaY2();
586 Float_t w=fParam->GetPadPitchWidth(fSector);
587 c.SetSigmaY2(s2*w*w);
589 c.SetSigmaZ2(s2*fZWidth*fZWidth);
594 AliTPCTransform *transform = AliTPCcalibDB::Instance()->GetTransform() ;
596 AliFatal("Tranformations not in calibDB");
599 transform->SetCurrentRecoParam((AliTPCRecoParam*)fRecoParam);
600 Double_t x[3]={c.GetRow(),c.GetPad(),c.GetTimeBin()};
601 Int_t i[1]={fSector};
602 transform->Transform(x,i,0,1);
608 if (ki<=1 || ki>=fMaxPad-1 || kj==1 || kj==fMaxTime-2) {
609 c.SetType(-(c.GetType()+3)); //edge clusters
611 if (fLoop==2) c.SetType(100);
614 TClonesArray * arr = 0;
615 AliTPCclusterMI * cl = 0;
617 if(fBClonesArray==kFALSE) {
618 arr = fRowCl->GetArray();
619 cl = new ((*arr)[fNcluster]) AliTPCclusterMI(c);
621 cl = new ((*fOutputClonesArray)[fNclusters+fNcluster]) AliTPCclusterMI(c);
624 // if (fRecoParam->DumpSignal() &&matrix ) {
626 // Float_t *graph =0;
627 // if (fRecoParam->GetCalcPedestal() && cl->GetMax()>fRecoParam->GetDumpAmplitudeMin() &&fBDumpSignal){
629 // graph = &(fBins[fMaxTime*(pos/fMaxTime)]);
631 // AliTPCclusterInfo * info = new AliTPCclusterInfo(matrix,nbins,graph);
632 // cl->SetInfo(info);
634 if (!fRecoParam->DumpSignal()) {
638 if (AliTPCReconstructor::StreamLevel()>1) {
640 cl->GetGlobalXYZ(xyz);
641 (*fDebugStreamer)<<"Clusters"<<
653 //_____________________________________________________________________________
654 void AliTPCclustererMI::Digits2Clusters()
656 //-----------------------------------------------------------------
657 // This is a simple cluster finder.
658 //-----------------------------------------------------------------
661 Error("Digits2Clusters", "input tree not initialised");
664 fRecoParam = AliTPCReconstructor::GetRecoParam();
666 AliFatal("Can not get the reconstruction parameters");
668 if(AliTPCReconstructor::StreamLevel()>5) {
669 AliInfo("Parameter Dumps");
675 //-----------------------------------------------------------------
677 //-----------------------------------------------------------------
678 if (fUseHLTClusters == 3 || fUseHLTClusters == 4) {
679 AliInfo("Using HLT clusters for TPC off-line reconstruction");
680 fZWidth = fParam->GetZWidth();
681 Int_t iResult = ReadHLTClusters();
683 // HLT clusters present
684 if (iResult >= 0 && fNclusters > 0)
687 // HLT clusters not present
688 if (iResult < 0 || fNclusters == 0) {
689 if (fUseHLTClusters == 3) {
690 AliError("No HLT clusters present, but requested.");
694 AliInfo("Now trying to read from TPC RAW");
697 // Some other problem during cluster reading
699 AliWarning("Some problem while unpacking of HLT clusters.");
702 } // if (fUseHLTClusters == 3 || fUseHLTClusters == 4) {
704 //-----------------------------------------------------------------
705 // Run TPC off-line clusterer
706 //-----------------------------------------------------------------
707 AliTPCCalPad * gainTPC = AliTPCcalibDB::Instance()->GetPadGainFactor();
708 AliTPCCalPad * noiseTPC = AliTPCcalibDB::Instance()->GetPadNoise();
709 AliSimDigits digarr, *dummy=&digarr;
711 fInput->GetBranch("Segment")->SetAddress(&dummy);
712 Stat_t nentries = fInput->GetEntries();
714 fMaxTime=fRecoParam->GetLastBin()+6; // add 3 virtual time bins before and 3 after
718 for (Int_t n=0; n<nentries; n++) {
720 if (!fParam->AdjustSectorRow(digarr.GetID(),fSector,fRow)) {
721 cerr<<"AliTPC warning: invalid segment ID ! "<<digarr.GetID()<<endl;
725 AliTPCCalROC * gainROC = gainTPC->GetCalROC(fSector); // pad gains per given sector
726 AliTPCCalROC * noiseROC = noiseTPC->GetCalROC(fSector); // noise per given sector
729 fRowCl->SetID(digarr.GetID());
730 if (fOutput) fOutput->GetBranch("Segment")->SetAddress(&fRowCl);
731 fRx=fParam->GetPadRowRadii(fSector,row);
734 const Int_t kNIS=fParam->GetNInnerSector(), kNOS=fParam->GetNOuterSector();
735 fZWidth = fParam->GetZWidth();
736 if (fSector < kNIS) {
737 fMaxPad = fParam->GetNPadsLow(row);
738 fSign = (fSector < kNIS/2) ? 1 : -1;
739 fPadLength = fParam->GetPadPitchLength(fSector,row);
740 fPadWidth = fParam->GetPadPitchWidth();
742 fMaxPad = fParam->GetNPadsUp(row);
743 fSign = ((fSector-kNIS) < kNOS/2) ? 1 : -1;
744 fPadLength = fParam->GetPadPitchLength(fSector,row);
745 fPadWidth = fParam->GetPadPitchWidth();
749 fMaxBin=fMaxTime*(fMaxPad+6); // add 3 virtual pads before and 3 after
750 fBins =new Float_t[fMaxBin];
751 fSigBins =new Int_t[fMaxBin];
753 memset(fBins,0,sizeof(Float_t)*fMaxBin);
755 if (digarr.First()) //MI change
757 Float_t dig=digarr.CurrentDigit();
758 if (dig<=fParam->GetZeroSup()) continue;
759 Int_t j=digarr.CurrentRow()+3, i=digarr.CurrentColumn()+3;
760 Float_t gain = gainROC->GetValue(row,digarr.CurrentColumn());
761 Int_t bin = i*fMaxTime+j;
767 fSigBins[fNSigBins++]=bin;
768 } while (digarr.Next());
769 digarr.ExpandTrackBuffer();
771 FindClusters(noiseROC);
773 fRowCl->GetArray()->Clear("C");
774 nclusters+=fNcluster;
780 Info("Digits2Clusters", "Number of found clusters : %d", nclusters);
782 if (fUseHLTClusters == 2 && nclusters == 0) {
783 AliInfo("No clusters from TPC Raw data, now trying to read HLT clusters.");
785 fZWidth = fParam->GetZWidth();
790 void AliTPCclustererMI::ProcessSectorData(){
792 // Process the data for the current sector
795 AliTPCCalPad * pedestalTPC = AliTPCcalibDB::Instance()->GetPedestals();
796 AliTPCCalPad * noiseTPC = AliTPCcalibDB::Instance()->GetPadNoise();
797 AliTPCCalROC * pedestalROC = pedestalTPC->GetCalROC(fSector); // pedestal per given sector
798 AliTPCCalROC * noiseROC = noiseTPC->GetCalROC(fSector); // noise per given sector
799 //check the presence of the calibration
800 if (!noiseROC ||!pedestalROC ) {
801 AliError(Form("Missing calibration per sector\t%d\n",fSector));
804 Int_t nRows=fParam->GetNRow(fSector);
805 Bool_t calcPedestal = fRecoParam->GetCalcPedestal();
806 Int_t zeroSup = fParam->GetZeroSup();
807 // if (calcPedestal) {
809 for (Int_t iRow = 0; iRow < nRows; iRow++) {
810 Int_t maxPad = fParam->GetNPads(fSector, iRow);
812 for (Int_t iPad = 3; iPad < maxPad + 3; iPad++) {
814 // Temporary fix for data production - !!!! MARIAN
815 // The noise calibration should take mean and RMS - currently the Gaussian fit used
816 // In case of double peak - the pad should be rejected
818 // Line mean - if more than given digits over threshold - make a noise calculation
819 // and pedestal substration
820 if (!calcPedestal && fAllBins[iRow][iPad*fMaxTime+0]<50) continue;
822 if (fAllBins[iRow][iPad*fMaxTime+0] <1 ) continue; // no data
823 Float_t *p = &fAllBins[iRow][iPad*fMaxTime+3];
824 //Float_t pedestal = TMath::Median(fMaxTime, p);
825 Int_t id[3] = {fSector, iRow, iPad-3};
827 Double_t rmsCalib= noiseROC->GetValue(iRow,iPad-3);
828 Double_t pedestalCalib = pedestalROC->GetValue(iRow,iPad-3);
829 Double_t rmsEvent = rmsCalib;
830 Double_t pedestalEvent = pedestalCalib;
831 ProcesSignal(p, fMaxTime, id, rmsEvent, pedestalEvent);
832 if (rmsEvent<rmsCalib) rmsEvent = rmsCalib; // take worst scenario
833 if (TMath::Abs(pedestalEvent-pedestalCalib)<1.0) pedestalEvent = pedestalCalib;
836 for (Int_t iTimeBin = 0; iTimeBin < fMaxTime; iTimeBin++) {
837 Int_t bin = iPad*fMaxTime+iTimeBin;
838 fAllBins[iRow][bin] -= pedestalEvent;
839 if (iTimeBin < fRecoParam->GetFirstBin())
840 fAllBins[iRow][bin] = 0;
841 if (iTimeBin > fRecoParam->GetLastBin())
842 fAllBins[iRow][bin] = 0;
843 if (fAllBins[iRow][iPad*fMaxTime+iTimeBin] < zeroSup)
844 fAllBins[iRow][bin] = 0;
845 if (fAllBins[iRow][bin] < 3.0*rmsEvent) // 3 sigma cut on RMS
846 fAllBins[iRow][bin] = 0;
847 if (fAllBins[iRow][bin]) fAllSigBins[iRow][fAllNSigBins[iRow]++] = bin;
853 if (AliTPCReconstructor::StreamLevel()>5) {
854 for (Int_t iRow = 0; iRow < nRows; iRow++) {
855 Int_t maxPad = fParam->GetNPads(fSector,iRow);
857 for (Int_t iPad = 3; iPad < maxPad + 3; iPad++) {
858 for (Int_t iTimeBin = 0; iTimeBin < fMaxTime; iTimeBin++) {
859 Int_t bin = iPad*fMaxTime+iTimeBin;
860 Float_t signal = fAllBins[iRow][bin];
861 if (AliTPCReconstructor::StreamLevel()>3 && signal>3) {
862 Double_t x[]={iRow,iPad-3,iTimeBin-3};
864 AliTPCTransform trafo;
865 trafo.Transform(x,i,0,1);
866 Double_t gx[3]={x[0],x[1],x[2]};
867 trafo.RotatedGlobal2Global(fSector,gx);
868 // fAllSigBins[iRow][fAllNSigBins[iRow]++]
869 Int_t rowsigBins = fAllNSigBins[iRow];
870 Int_t first=fAllSigBins[iRow][0];
872 // if (rowsigBins>0) fAllSigBins[iRow][fAllNSigBins[iRow]-1];
874 if (AliTPCReconstructor::StreamLevel()>5) {
875 (*fDebugStreamer)<<"Digits"<<
888 "rowsigBins="<<rowsigBins<<
899 // Now loop over rows and find clusters
900 for (fRow = 0; fRow < nRows; fRow++) {
901 fRowCl->SetID(fParam->GetIndex(fSector, fRow));
902 if (fOutput) fOutput->GetBranch("Segment")->SetAddress(&fRowCl);
904 fRx = fParam->GetPadRowRadii(fSector, fRow);
905 fPadLength = fParam->GetPadPitchLength(fSector, fRow);
906 fPadWidth = fParam->GetPadPitchWidth();
907 fMaxPad = fParam->GetNPads(fSector,fRow);
908 fMaxBin = fMaxTime*(fMaxPad+6); // add 3 virtual pads before and 3 after
910 fBins = fAllBins[fRow];
911 fSigBins = fAllSigBins[fRow];
912 fNSigBins = fAllNSigBins[fRow];
914 FindClusters(noiseROC);
917 if(fBClonesArray == kFALSE) fRowCl->GetArray()->Clear("C");
918 fNclusters += fNcluster;
920 } // End of loop to find clusters
924 void AliTPCclustererMI::Digits2Clusters(AliRawReader* rawReader)
926 //-----------------------------------------------------------------
927 // This is a cluster finder for the TPC raw data.
928 // The method assumes NO ordering of the altro channels.
929 // The pedestal subtraction can be switched on and off
930 // using an option of the TPC reconstructor
931 //-----------------------------------------------------------------
932 fRecoParam = AliTPCReconstructor::GetRecoParam();
934 AliFatal("Can not get the reconstruction parameters");
936 if(AliTPCReconstructor::StreamLevel()>5) {
937 AliInfo("Parameter Dumps");
943 //-----------------------------------------------------------------
945 //-----------------------------------------------------------------
946 if (fUseHLTClusters == 3 || fUseHLTClusters == 4) {
947 AliInfo("Using HLT clusters for TPC off-line reconstruction");
948 fZWidth = fParam->GetZWidth();
949 Int_t iResult = ReadHLTClusters();
951 // HLT clusters present
952 if (iResult >= 0 && fNclusters > 0)
955 // HLT clusters not present
956 if (iResult < 0 || fNclusters == 0) {
957 if (fUseHLTClusters == 3) {
958 AliError("No HLT clusters present, but requested.");
962 AliInfo("Now trying to read TPC RAW");
965 // Some other problem during cluster reading
967 AliWarning("Some problem while unpacking of HLT clusters.");
970 } // if (fUseHLTClusters == 3 || fUseHLTClusters == 4) {
972 //-----------------------------------------------------------------
973 // Run TPC off-line clusterer
974 //-----------------------------------------------------------------
975 AliTPCCalPad * gainTPC = AliTPCcalibDB::Instance()->GetPadGainFactor();
976 AliTPCAltroMapping** mapping =AliTPCcalibDB::Instance()->GetMapping();
978 AliTPCRawStreamV3 input(rawReader,(AliAltroMapping**)mapping);
979 fEventHeader = (AliRawEventHeaderBase*)rawReader->GetEventHeader();
981 fTimeStamp = fEventHeader->Get("Timestamp");
982 fEventType = fEventHeader->Get("Type");
983 AliTPCTransform *transform = AliTPCcalibDB::Instance()->GetTransform() ;
984 transform->SetCurrentTimeStamp(fTimeStamp);
985 transform->SetCurrentRun(rawReader->GetRunNumber());
988 // creaate one TClonesArray for all clusters
989 if(fBClonesArray && !fOutputClonesArray) fOutputClonesArray = new TClonesArray("AliTPCclusterMI",1000);
993 fMaxTime = fRecoParam->GetLastBin() + 6; // add 3 virtual time bins before and 3 after
994 // const Int_t kNIS = fParam->GetNInnerSector();
995 // const Int_t kNOS = fParam->GetNOuterSector();
996 // const Int_t kNS = kNIS + kNOS;
997 fZWidth = fParam->GetZWidth();
998 Int_t zeroSup = fParam->GetZeroSup();
1002 AliTPCROC * roc = AliTPCROC::Instance();
1003 Int_t nRowsMax = roc->GetNRows(roc->GetNSector()-1);
1004 Int_t nPadsMax = roc->GetNPads(roc->GetNSector()-1,nRowsMax-1);
1005 for (Int_t iRow = 0; iRow < nRowsMax; iRow++) {
1007 Int_t maxBin = fMaxTime*(nPadsMax+6); // add 3 virtual pads before and 3 after
1008 memset(fAllBins[iRow],0,sizeof(Float_t)*maxBin);
1009 fAllNSigBins[iRow]=0;
1012 Int_t prevSector=-1;
1018 const Int_t kNIS = fParam->GetNInnerSector();
1019 const Int_t kNOS = fParam->GetNOuterSector();
1020 const Int_t kNS = kNIS + kNOS;
1022 for(fSector = 0; fSector < kNS; fSector++) {
1025 Int_t nDDLs = 0, indexDDL = 0;
1026 if (fSector < kNIS) {
1027 nRows = fParam->GetNRowLow();
1028 fSign = (fSector < kNIS/2) ? 1 : -1;
1030 indexDDL = fSector * 2;
1033 nRows = fParam->GetNRowUp();
1034 fSign = ((fSector-kNIS) < kNOS/2) ? 1 : -1;
1036 indexDDL = (fSector-kNIS) * 4 + kNIS * 2;
1039 // load the raw data for corresponding DDLs
1041 rawReader->Select("TPC",indexDDL,indexDDL+nDDLs-1);
1043 while (input.NextDDL()){
1044 if (input.GetSector() != fSector)
1045 AliFatal(Form("Sector index mismatch ! Expected (%d), but got (%d) !",fSector,input.GetSector()));
1047 //Int_t nRows = fParam->GetNRow(fSector);
1049 AliTPCCalROC * gainROC = gainTPC->GetCalROC(fSector); // pad gains per given sector
1050 // Begin loop over altro data
1051 Bool_t calcPedestal = fRecoParam->GetCalcPedestal();
1055 while ( input.NextChannel() ) {
1056 Int_t iRow = input.GetRow();
1061 AliError(Form("Pad-row index (%d) outside the range (%d -> %d) !",
1062 iRow, 0, nRows -1));
1066 Int_t iPad = input.GetPad();
1067 if (iPad < 0 || iPad >= nPadsMax) {
1068 AliError(Form("Pad index (%d) outside the range (%d -> %d) !",
1069 iPad, 0, nPadsMax-1));
1072 gain = gainROC->GetValue(iRow,iPad);
1076 while ( input.NextBunch() ){
1077 Int_t startTbin = (Int_t)input.GetStartTimeBin();
1078 Int_t bunchlength = (Int_t)input.GetBunchLength();
1079 const UShort_t *sig = input.GetSignals();
1080 for (Int_t iTime = 0; iTime<bunchlength; iTime++){
1081 Int_t iTimeBin=startTbin-iTime;
1082 if ( iTimeBin < fRecoParam->GetFirstBin() || iTimeBin >= fRecoParam->GetLastBin()){
1084 AliFatal(Form("Timebin index (%d) outside the range (%d -> %d) !",
1085 iTimeBin, 0, iTimeBin -1));
1089 Float_t signal=(Float_t)sig[iTime];
1090 if (!calcPedestal && signal <= zeroSup) continue;
1092 if (!calcPedestal) {
1093 Int_t bin = iPad*fMaxTime+iTimeBin;
1095 fAllBins[iRow][bin] = signal/gain;
1097 fAllBins[iRow][bin] =0;
1099 fAllSigBins[iRow][fAllNSigBins[iRow]++] = bin;
1101 fAllBins[iRow][iPad*fMaxTime+iTimeBin] = signal;
1103 fAllBins[iRow][iPad*fMaxTime+0]+=1.; // pad with signal
1107 }// end loop signals in bunch
1108 }// end loop bunches
1114 // Now loop over rows and perform pedestal subtraction
1115 if (digCounter==0) continue;
1116 } // End of loop over sectors
1117 //process last sector
1118 if ( digCounter>0 ){
1119 ProcessSectorData();
1120 for (Int_t iRow = 0; iRow < fParam->GetNRow(fSector); iRow++) {
1121 Int_t maxPad = fParam->GetNPads(fSector,iRow);
1122 Int_t maxBin = fMaxTime*(maxPad+6); // add 3 virtual pads before and 3 after
1123 memset(fAllBins[iRow],0,sizeof(Float_t)*maxBin);
1124 fAllNSigBins[iRow] = 0;
1131 if (rawReader->GetEventId() && fOutput ){
1132 Info("Digits2Clusters", "File %s Event\t%u\tNumber of found clusters : %d\n", fOutput->GetName(),*(rawReader->GetEventId()), fNclusters);
1135 if(rawReader->GetEventId()) {
1136 Info("Digits2Clusters", "Event\t%u\tNumber of found clusters : %d\n",*(rawReader->GetEventId()), fNclusters);
1140 //Info("Digits2Clusters", "Number of found clusters : %d\n",fOutputClonesArray->GetEntriesFast());
1143 if (fUseHLTClusters == 2 && fNclusters == 0) {
1144 AliInfo("No clusters from TPC Raw data, now trying to read HLT clusters.");
1146 fZWidth = fParam->GetZWidth();
1155 void AliTPCclustererMI::Digits2ClustersOld
1156 (AliRawReader* rawReader)
1158 //-----------------------------------------------------------------
1159 // This is a cluster finder for the TPC raw data.
1160 // The method assumes NO ordering of the altro channels.
1161 // The pedestal subtraction can be switched on and off
1162 // using an option of the TPC reconstructor
1163 //-----------------------------------------------------------------
1164 fRecoParam = AliTPCReconstructor::GetRecoParam();
1166 AliFatal("Can not get the reconstruction parameters");
1168 if(AliTPCReconstructor::StreamLevel()>5) {
1169 AliInfo("Parameter Dumps");
1175 AliTPCCalPad * gainTPC = AliTPCcalibDB::Instance()->GetPadGainFactor();
1176 AliTPCAltroMapping** mapping =AliTPCcalibDB::Instance()->GetMapping();
1178 AliTPCRawStream input(rawReader,(AliAltroMapping**)mapping);
1179 fEventHeader = (AliRawEventHeaderBase*)rawReader->GetEventHeader();
1181 fTimeStamp = fEventHeader->Get("Timestamp");
1182 fEventType = fEventHeader->Get("Type");
1185 // creaate one TClonesArray for all clusters
1186 if(fBClonesArray && !fOutputClonesArray) fOutputClonesArray = new TClonesArray("AliTPCclusterMI",1000);
1190 fMaxTime = fRecoParam->GetLastBin() + 6; // add 3 virtual time bins before and 3 after
1191 const Int_t kNIS = fParam->GetNInnerSector();
1192 const Int_t kNOS = fParam->GetNOuterSector();
1193 const Int_t kNS = kNIS + kNOS;
1194 fZWidth = fParam->GetZWidth();
1195 Int_t zeroSup = fParam->GetZeroSup();
1200 AliTPCROC * roc = AliTPCROC::Instance();
1201 Int_t nRowsMax = roc->GetNRows(roc->GetNSector()-1);
1202 Int_t nPadsMax = roc->GetNPads(roc->GetNSector()-1,nRowsMax-1);
1203 for (Int_t iRow = 0; iRow < nRowsMax; iRow++) {
1205 Int_t maxBin = fMaxTime*(nPadsMax+6); // add 3 virtual pads before and 3 after
1206 memset(fAllBins[iRow],0,sizeof(Float_t)*maxBin);
1207 fAllNSigBins[iRow]=0;
1210 // Loop over sectors
1212 for(fSector = 0; fSector < kNS; fSector++) {
1215 Int_t nDDLs = 0, indexDDL = 0;
1216 if (fSector < kNIS) {
1217 nRows = fParam->GetNRowLow();
1218 fSign = (fSector < kNIS/2) ? 1 : -1;
1220 indexDDL = fSector * 2;
1223 nRows = fParam->GetNRowUp();
1224 fSign = ((fSector-kNIS) < kNOS/2) ? 1 : -1;
1226 indexDDL = (fSector-kNIS) * 4 + kNIS * 2;
1229 // load the raw data for corresponding DDLs
1231 rawReader->Select("TPC",indexDDL,indexDDL+nDDLs-1);
1233 // select only good sector
1234 if (!input.Next()) continue;
1235 if(input.GetSector() != fSector) continue;
1237 AliTPCCalROC * gainROC = gainTPC->GetCalROC(fSector); // pad gains per given sector
1239 for (Int_t iRow = 0; iRow < nRows; iRow++) {
1242 maxPad = fParam->GetNPadsLow(iRow);
1244 maxPad = fParam->GetNPadsUp(iRow);
1246 Int_t maxBin = fMaxTime*(maxPad+6); // add 3 virtual pads before and 3 after
1247 memset(fAllBins[iRow],0,sizeof(Float_t)*maxBin);
1248 fAllNSigBins[iRow] = 0;
1252 // Begin loop over altro data
1253 Bool_t calcPedestal = fRecoParam->GetCalcPedestal();
1258 while (input.Next()) {
1259 if (input.GetSector() != fSector)
1260 AliFatal(Form("Sector index mismatch ! Expected (%d), but got (%d) !",fSector,input.GetSector()));
1263 Int_t iRow = input.GetRow();
1268 if (iRow < 0 || iRow >= nRows){
1269 AliError(Form("Pad-row index (%d) outside the range (%d -> %d) !",
1270 iRow, 0, nRows -1));
1274 Int_t iPad = input.GetPad();
1275 if (iPad < 0 || iPad >= nPadsMax) {
1276 AliError(Form("Pad index (%d) outside the range (%d -> %d) !",
1277 iPad, 0, nPadsMax-1));
1281 gain = gainROC->GetValue(iRow,iPad);
1286 Int_t iTimeBin = input.GetTime();
1287 if ( iTimeBin < fRecoParam->GetFirstBin() || iTimeBin >= fRecoParam->GetLastBin()){
1289 AliFatal(Form("Timebin index (%d) outside the range (%d -> %d) !",
1290 iTimeBin, 0, iTimeBin -1));
1295 Float_t signal = input.GetSignal();
1296 if (!calcPedestal && signal <= zeroSup) continue;
1298 if (!calcPedestal) {
1299 Int_t bin = iPad*fMaxTime+iTimeBin;
1301 fAllBins[iRow][bin] = signal/gain;
1303 fAllBins[iRow][bin] =0;
1305 fAllSigBins[iRow][fAllNSigBins[iRow]++] = bin;
1307 fAllBins[iRow][iPad*fMaxTime+iTimeBin] = signal;
1309 fAllBins[iRow][iPad*fMaxTime+0]+=1.; // pad with signal
1313 } // End of the loop over altro data
1318 // Now loop over rows and perform pedestal subtraction
1319 if (digCounter==0) continue;
1320 ProcessSectorData();
1321 } // End of loop over sectors
1323 if (rawReader->GetEventId() && fOutput ){
1324 Info("Digits2Clusters", "File %s Event\t%d\tNumber of found clusters : %d\n", fOutput->GetName(),*(rawReader->GetEventId()), fNclusters);
1327 if(rawReader->GetEventId()) {
1328 Info("Digits2Clusters", "Event\t%d\tNumber of found clusters : %d\n",*(rawReader->GetEventId()), fNclusters);
1332 //Info("Digits2Clusters", "Number of found clusters : %d\n",fOutputClonesArray->GetEntriesFast());
1336 void AliTPCclustererMI::FindClusters(AliTPCCalROC * noiseROC)
1340 // add virtual charge at the edge
1342 Double_t kMaxDumpSize = 500000;
1344 fBDumpSignal =kFALSE;
1346 if (fRecoParam->GetCalcPedestal() && fOutput->GetZipBytes()< kMaxDumpSize) fBDumpSignal =kTRUE; //dump signal flag
1351 Int_t crtime = Int_t((fParam->GetZLength(fSector)-fRecoParam->GetCtgRange()*fRx)/fZWidth-fParam->GetNTBinsL1()-5);
1352 Float_t minMaxCutAbs = fRecoParam->GetMinMaxCutAbs();
1353 Float_t minLeftRightCutAbs = fRecoParam->GetMinLeftRightCutAbs();
1354 Float_t minUpDownCutAbs = fRecoParam->GetMinUpDownCutAbs();
1355 Float_t minMaxCutSigma = fRecoParam->GetMinMaxCutSigma();
1356 Float_t minLeftRightCutSigma = fRecoParam->GetMinLeftRightCutSigma();
1357 Float_t minUpDownCutSigma = fRecoParam->GetMinUpDownCutSigma();
1358 Int_t useOnePadCluster = fRecoParam->GetUseOnePadCluster();
1359 for (Int_t iSig = 0; iSig < fNSigBins; iSig++) {
1360 Int_t i = fSigBins[iSig];
1361 if (i%fMaxTime<=crtime) continue;
1362 Float_t *b = &fBins[i];
1364 if (b[0]<minMaxCutAbs) continue; //threshold for maxima
1366 if (useOnePadCluster==0){
1367 if (b[-1]+b[1]+b[-fMaxTime]+b[fMaxTime]<=0) continue; // cut on isolated clusters
1368 if (b[-1]+b[1]<=0) continue; // cut on isolated clusters
1369 if (b[-fMaxTime]+b[fMaxTime]<=0) continue; // cut on isolated clusters
1372 if ((b[0]+b[-1]+b[1])<minUpDownCutAbs) continue; //threshold for up down (TRF)
1373 if ((b[0]+b[-fMaxTime]+b[fMaxTime])<minLeftRightCutAbs) continue; //threshold for left right (PRF)
1374 if (!IsMaximum(*b,fMaxTime,b)) continue;
1376 Float_t noise = noiseROC->GetValue(fRow, i/fMaxTime);
1377 if (noise>fRecoParam->GetMaxNoise()) continue;
1379 if (b[0]<minMaxCutSigma*noise) continue; //threshold form maxima
1380 if ((b[0]+b[-1]+b[1])<minUpDownCutSigma*noise) continue; //threshold for up town TRF
1381 if ((b[0]+b[-fMaxTime]+b[fMaxTime])<minLeftRightCutSigma*noise) continue; //threshold for left right (PRF)
1383 AliTPCclusterMI c; // default cosntruction without info
1385 MakeCluster(i, fMaxTime, fBins, dummy,c);
1391 Bool_t AliTPCclustererMI::AcceptCluster(AliTPCclusterMI *cl){
1393 // Currently hack to filter digital noise (15.06.2008)
1394 // To be parameterized in the AliTPCrecoParam
1395 // More inteligent way to be used in future
1396 // Acces to the proper pedestal file needed
1398 if (cl->GetMax()<400) return kTRUE;
1399 Double_t ratio = cl->GetQ()/cl->GetMax();
1400 if (cl->GetMax()>700){
1401 if ((ratio - int(ratio)>0.8)) return kFALSE;
1403 if ((ratio - int(ratio)<0.95)) return kTRUE;
1408 Double_t AliTPCclustererMI::ProcesSignal(Float_t *signal, Int_t nchannels, Int_t id[3], Double_t &rmsEvent, Double_t &pedestalEvent){
1410 // process signal on given pad - + streaming of additional information in special mode
1417 // ESTIMATE pedestal and the noise
1419 const Int_t kPedMax = 100;
1425 Float_t rmsCalib = rmsEvent; // backup initial value ( from calib)
1426 Float_t pedestalCalib = pedestalEvent;// backup initial value ( from calib)
1427 Int_t firstBin = fRecoParam->GetFirstBin();
1429 UShort_t histo[kPedMax];
1430 //memset(histo,0,kPedMax*sizeof(UShort_t));
1431 for (Int_t i=0; i<kPedMax; i++) histo[i]=0;
1432 for (Int_t i=0; i<fMaxTime; i++){
1433 if (signal[i]<=0) continue;
1434 if (signal[i]>max && i>firstBin) {
1438 if (signal[i]>kPedMax-1) continue;
1439 histo[int(signal[i]+0.5)]++;
1443 for (Int_t i=1; i<kPedMax; i++){
1444 if (count1<count0*0.5) median=i;
1449 Float_t count10=histo[median] ,mean=histo[median]*median, rms=histo[median]*median*median ;
1450 Float_t count06=histo[median] ,mean06=histo[median]*median, rms06=histo[median]*median*median ;
1451 Float_t count09=histo[median] ,mean09=histo[median]*median, rms09=histo[median]*median*median ;
1453 for (Int_t idelta=1; idelta<10; idelta++){
1454 if (median-idelta<=0) continue;
1455 if (median+idelta>kPedMax) continue;
1456 if (count06<0.6*count1){
1457 count06+=histo[median-idelta];
1458 mean06 +=histo[median-idelta]*(median-idelta);
1459 rms06 +=histo[median-idelta]*(median-idelta)*(median-idelta);
1460 count06+=histo[median+idelta];
1461 mean06 +=histo[median+idelta]*(median+idelta);
1462 rms06 +=histo[median+idelta]*(median+idelta)*(median+idelta);
1464 if (count09<0.9*count1){
1465 count09+=histo[median-idelta];
1466 mean09 +=histo[median-idelta]*(median-idelta);
1467 rms09 +=histo[median-idelta]*(median-idelta)*(median-idelta);
1468 count09+=histo[median+idelta];
1469 mean09 +=histo[median+idelta]*(median+idelta);
1470 rms09 +=histo[median+idelta]*(median+idelta)*(median+idelta);
1472 if (count10<0.95*count1){
1473 count10+=histo[median-idelta];
1474 mean +=histo[median-idelta]*(median-idelta);
1475 rms +=histo[median-idelta]*(median-idelta)*(median-idelta);
1476 count10+=histo[median+idelta];
1477 mean +=histo[median+idelta]*(median+idelta);
1478 rms +=histo[median+idelta]*(median+idelta)*(median+idelta);
1483 rms = TMath::Sqrt(TMath::Abs(rms/count10-mean*mean));
1487 rms06 = TMath::Sqrt(TMath::Abs(rms06/count06-mean06*mean06));
1491 rms09 = TMath::Sqrt(TMath::Abs(rms09/count09-mean09*mean09));
1495 pedestalEvent = median;
1496 if (AliLog::GetDebugLevel("","AliTPCclustererMI")==0) return median;
1498 UInt_t uid[3] = {UInt_t(id[0]),UInt_t(id[1]),UInt_t(id[2])};
1500 // Dump mean signal info
1502 if (AliTPCReconstructor::StreamLevel()>0) {
1503 (*fDebugStreamer)<<"Signal"<<
1504 "TimeStamp="<<fTimeStamp<<
1505 "EventType="<<fEventType<<
1519 "RMSCalib="<<rmsCalib<<
1520 "PedCalib="<<pedestalCalib<<
1524 // fill pedestal histogram
1529 Float_t kMin =fRecoParam->GetDumpAmplitudeMin(); // minimal signal to be dumped
1530 Float_t *dsignal = new Float_t[nchannels];
1531 Float_t *dtime = new Float_t[nchannels];
1532 for (Int_t i=0; i<nchannels; i++){
1534 dsignal[i] = signal[i];
1538 // Big signals dumping
1540 if (AliTPCReconstructor::StreamLevel()>0) {
1541 if (max-median>kMin &&maxPos>fRecoParam->GetFirstBin())
1542 (*fDebugStreamer)<<"SignalB"<< // pads with signal
1543 "TimeStamp="<<fTimeStamp<<
1544 "EventType="<<fEventType<<
1548 // "Graph="<<graph<<
1564 if (rms06>fRecoParam->GetMaxNoise()) {
1565 pedestalEvent+=1024.;
1566 return 1024+median; // sign noisy channel in debug mode
1571 Int_t AliTPCclustererMI::ReadHLTClusters()
1574 // read HLT clusters instead of off line custers,
1575 // used in Digits2Clusters
1578 if (!fHLTClusterAccess) {
1580 ROOT::NewFunc_t pNewFunc=NULL;
1582 pCl=TClass::GetClass("AliHLTTPCClusterAccessHLTOUT");
1583 } while (!pCl && gSystem->Load("libAliHLTTPC.so")==0);
1584 if (!pCl || (pNewFunc=pCl->GetNew())==NULL) {
1585 AliError("can not load class description of AliHLTTPCClusterAccessHLTOUT, aborting ...");
1589 void* p=(*pNewFunc)(NULL);
1591 AliError("unable to create instance of AliHLTTPCClusterAccessHLTOUT");
1594 fHLTClusterAccess=reinterpret_cast<TObject*>(p);
1597 TObject* pClusterAccess=fHLTClusterAccess;
1599 const Int_t kNIS = fParam->GetNInnerSector();
1600 const Int_t kNOS = fParam->GetNOuterSector();
1601 const Int_t kNS = kNIS + kNOS;
1604 // make sure that all clusters from the previous event are cleared
1605 pClusterAccess->Clear("event");
1606 for(fSector = 0; fSector < kNS; fSector++) {
1609 TString param("sector="); param+=fSector;
1610 // prepare for next sector
1611 pClusterAccess->Clear("sector");
1612 pClusterAccess->Execute("read", param, &iResult);
1615 AliError("HLT Clusters can not be found");
1618 if (pClusterAccess->FindObject("clusterarray")==NULL) {
1619 AliError("HLT clusters requested, but not cluster array not present");
1623 TClonesArray* clusterArray=dynamic_cast<TClonesArray*>(pClusterAccess->FindObject("clusterarray"));
1624 if (!clusterArray) {
1625 AliError("HLT cluster array is not of class type TClonesArray");
1629 AliDebug(4,Form("Reading %d clusters from HLT for sector %d", clusterArray->GetEntriesFast(), fSector));
1631 Int_t nClusterSector=0;
1632 Int_t nRows=fParam->GetNRow(fSector);
1634 for (fRow = 0; fRow < nRows; fRow++) {
1635 fRowCl->SetID(fParam->GetIndex(fSector, fRow));
1636 if (fOutput) fOutput->GetBranch("Segment")->SetAddress(&fRowCl);
1637 fNcluster=0; // reset clusters per row
1639 fRx = fParam->GetPadRowRadii(fSector, fRow);
1640 fPadLength = fParam->GetPadPitchLength(fSector, fRow);
1641 fPadWidth = fParam->GetPadPitchWidth();
1642 fMaxPad = fParam->GetNPads(fSector,fRow);
1643 fMaxBin = fMaxTime*(fMaxPad+6); // add 3 virtual pads before and 3 after
1645 fBins = fAllBins[fRow];
1646 fSigBins = fAllSigBins[fRow];
1647 fNSigBins = fAllNSigBins[fRow];
1649 for (Int_t i=0; i<clusterArray->GetEntriesFast(); i++) {
1650 if (!clusterArray->At(i))
1653 AliTPCclusterMI* cluster=dynamic_cast<AliTPCclusterMI*>(clusterArray->At(i));
1654 if (!cluster) continue;
1655 if (cluster->GetRow()!=fRow) continue;
1657 AddCluster(*cluster, NULL, 0);
1661 fRowCl->GetArray()->Clear("c");
1663 } // for (fRow = 0; fRow < nRows; fRow++) {
1664 if (nClusterSector!=clusterArray->GetEntriesFast()) {
1665 AliError(Form("Failed to read %d out of %d HLT clusters",
1666 clusterArray->GetEntriesFast()-nClusterSector,
1667 clusterArray->GetEntriesFast()));
1669 fNclusters+=nClusterSector;
1670 } // for(fSector = 0; fSector < kNS; fSector++) {
1672 pClusterAccess->Clear("event");
1674 Info("Digits2Clusters", "Number of converted HLT clusters : %d", fNclusters);