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"
87 ClassImp(AliTPCclustererMI)
91 AliTPCclustererMI::AliTPCclustererMI(const AliTPCParam* par, const AliTPCRecoParam * recoParam):
97 fMaxTime(1006), // 1000>940 so use 1000, add 3 virtual time bins before and 3 after
106 fPedSubtraction(kFALSE),
113 fOutputClonesArray(0),
121 fBDumpSignal(kFALSE),
122 fBClonesArray(kFALSE),
127 fHLTClusterAccess(NULL)
131 // param - tpc parameters for given file
132 // recoparam - reconstruction parameters
137 fRecoParam = recoParam;
139 //set default parameters if not specified
140 fRecoParam = AliTPCReconstructor::GetRecoParam();
141 if (!fRecoParam) fRecoParam = AliTPCRecoParam::GetLowFluxParam();
144 if(AliTPCReconstructor::StreamLevel()>0) {
145 fDebugStreamer = new TTreeSRedirector("TPCsignal.root");
148 // Int_t nPoints = fRecoParam->GetLastBin()-fRecoParam->GetFirstBin();
149 fRowCl= new AliTPCClustersRow("AliTPCclusterMI");
151 // Non-persistent arrays
153 //alocate memory for sector - maximal case
155 AliTPCROC * roc = AliTPCROC::Instance();
156 Int_t nRowsMax = roc->GetNRows(roc->GetNSector()-1);
157 Int_t nPadsMax = roc->GetNPads(roc->GetNSector()-1,nRowsMax-1);
159 fAllBins = new Float_t*[nRowsMax];
160 fAllSigBins = new Int_t*[nRowsMax];
161 fAllNSigBins = new Int_t[nRowsMax];
162 for (Int_t iRow = 0; iRow < nRowsMax; iRow++) {
164 Int_t maxBin = fMaxTime*(nPadsMax+6); // add 3 virtual pads before and 3 after
165 fAllBins[iRow] = new Float_t[maxBin];
166 memset(fAllBins[iRow],0,sizeof(Float_t)*maxBin);
167 fAllSigBins[iRow] = new Int_t[maxBin];
168 fAllNSigBins[iRow]=0;
172 //______________________________________________________________
173 AliTPCclustererMI::~AliTPCclustererMI(){
177 if (fDebugStreamer) delete fDebugStreamer;
179 //fOutputArray->Delete();
182 if (fOutputClonesArray){
183 fOutputClonesArray->Delete();
184 delete fOutputClonesArray;
188 fRowCl->GetArray()->Delete();
192 AliTPCROC * roc = AliTPCROC::Instance();
193 Int_t nRowsMax = roc->GetNRows(roc->GetNSector()-1);
194 for (Int_t iRow = 0; iRow < nRowsMax; iRow++) {
195 delete [] fAllBins[iRow];
196 delete [] fAllSigBins[iRow];
199 delete [] fAllSigBins;
200 delete [] fAllNSigBins;
201 if (fHLTClusterAccess) delete fHLTClusterAccess;
204 void AliTPCclustererMI::SetInput(TTree * tree)
207 // set input tree with digits
210 if (!fInput->GetBranch("Segment")){
211 cerr<<"AliTPC::Digits2Clusters(): no porper input tree !\n";
217 void AliTPCclustererMI::SetOutput(TTree * tree)
220 // Set the output tree
221 // If not set the ObjArray used - Option for HLT
225 AliTPCClustersRow clrow("AliTPCclusterMI");
226 AliTPCClustersRow *pclrow=&clrow;
227 fOutput->Branch("Segment","AliTPCClustersRow",&pclrow,32000,200);
231 void AliTPCclustererMI::FillRow(){
233 // fill the output container -
234 // 2 Options possible
238 if (fOutput) fOutput->Fill();
239 if (!fOutput && !fBClonesArray){
241 if (!fOutputArray) fOutputArray = new TObjArray(fParam->GetNRowsTotal());
242 if (fRowCl && fRowCl->GetArray()->GetEntriesFast()>0) fOutputArray->AddAt(fRowCl->Clone(), fRowCl->GetID());
246 Float_t AliTPCclustererMI::GetSigmaY2(Int_t iz){
247 // sigma y2 = in digits - we don't know the angle
248 Float_t z = iz*fParam->GetZWidth()+fParam->GetNTBinsL1()*fParam->GetZWidth();
249 Float_t sd2 = (z*fParam->GetDiffL()*fParam->GetDiffL())/
250 (fPadWidth*fPadWidth);
252 Float_t res = sd2+sres;
257 Float_t AliTPCclustererMI::GetSigmaZ2(Int_t iz){
258 //sigma z2 = in digits - angle estimated supposing vertex constraint
259 Float_t z = iz*fZWidth+fParam->GetNTBinsL1()*fParam->GetZWidth();
260 Float_t sd2 = (z*fParam->GetDiffL()*fParam->GetDiffL())/(fZWidth*fZWidth);
261 Float_t angular = fPadLength*(fParam->GetZLength(fSector)-z)/(fRx*fZWidth);
264 Float_t sres = fParam->GetZSigma()/fZWidth;
266 Float_t res = angular +sd2+sres;
270 void AliTPCclustererMI::MakeCluster(Int_t k,Int_t max,Float_t *bins, UInt_t /*m*/,
274 // Make cluster: characterized by position ( mean- COG) , shape (RMS) a charge, QMax and Q tot
275 // Additional correction:
276 // a) To correct for charge below threshold, in the +1 neghborhood to the max charge charge
277 // is extrapolated using gaussian approximation assuming given cluster width..
278 // Additional empirical factor is used to account for the charge fluctuation (kVirtualChargeFactor).
279 // Actual value of the kVirtualChargeFactor should obtained minimimizing residuals between the cluster
280 // and track interpolation.
281 // b.) For space points with extended shape (in comparison with expected using parameterization) clusters are
284 // NOTE. Actual/Empirical values for correction are hardwired in the code.
286 // Input paramters for function:
287 // k - Make cluster at position k
288 // bins - 2 D array of signals mapped to 1 dimensional array -
289 // max - the number of time bins er one dimension
290 // c - reference to cluster to be filled
292 Double_t kVirtualChargeFactor=0.5;
293 Int_t i0=k/max; //central pad
294 Int_t j0=k%max; //central time bin
296 // set pointers to data
297 //Int_t dummy[5] ={0,0,0,0,0};
298 Float_t * matrix[5]; //5x5 matrix with digits - indexing i = 0 ..4 j = -2..2
299 for (Int_t di=-2;di<=2;di++){
300 matrix[di+2] = &bins[k+di*max];
302 //build matrix with virtual charge
303 Float_t sigmay2= GetSigmaY2(j0);
304 Float_t sigmaz2= GetSigmaZ2(j0);
306 Float_t vmatrix[5][5];
307 vmatrix[2][2] = matrix[2][0];
309 c.SetMax((UShort_t)(vmatrix[2][2])); // write maximal amplitude
310 for (Int_t di =-1;di <=1;di++)
311 for (Int_t dj =-1;dj <=1;dj++){
312 Float_t amp = matrix[di+2][dj];
313 if ( (amp<2) && (fLoop<2)){
314 // if under threshold - calculate virtual charge
315 Float_t ratio = TMath::Exp(-1.2*TMath::Abs(di)/sigmay2)*TMath::Exp(-1.2*TMath::Abs(dj)/sigmaz2);
316 amp = ((matrix[2][0]-2)*(matrix[2][0]-2)/(matrix[-di+2][-dj]+2))*ratio;
318 vmatrix[2+di][2+dj]= kVirtualChargeFactor*amp;
319 vmatrix[2+2*di][2+2*dj]=0;
322 vmatrix[2+2*di][2+dj] =0;
323 vmatrix[2+di][2+2*dj] =0;
328 //if small amplitude - below 2 x threshold - don't consider other one
329 vmatrix[2+di][2+dj]=amp;
330 vmatrix[2+2*di][2+2*dj]=0; // don't take to the account next bin
333 vmatrix[2+2*di][2+dj] =0;
334 vmatrix[2+di][2+2*dj] =0;
338 //if bigger then take everything
339 vmatrix[2+di][2+dj]=amp;
340 vmatrix[2+2*di][2+2*dj]= matrix[2*di+2][2*dj] ;
343 vmatrix[2+2*di][2+dj] = matrix[2*di+2][dj];
344 vmatrix[2+di][2+2*dj] = matrix[2+di][dj*2];
356 for (Int_t i=-2;i<=2;i++)
357 for (Int_t j=-2;j<=2;j++){
358 Float_t amp = vmatrix[i+2][j+2];
367 Float_t meani = sumiw/sumw;
368 Float_t mi2 = sumi2w/sumw-meani*meani;
369 Float_t meanj = sumjw/sumw;
370 Float_t mj2 = sumj2w/sumw-meanj*meanj;
372 Float_t ry = mi2/sigmay2;
373 Float_t rz = mj2/sigmaz2;
376 if ( ( (ry<0.6) || (rz<0.6) ) && fLoop==2) return;
377 if ( ((ry <1.2) && (rz<1.2)) || (!fRecoParam->GetDoUnfold())) {
379 //if cluster looks like expected or Unfolding not switched on
380 //standard COG is used
381 //+1.2 deviation from expected sigma accepted
382 // c.fMax = FitMax(vmatrix,meani,meanj,TMath::Sqrt(sigmay2),TMath::Sqrt(sigmaz2));
386 //set cluster parameters
389 c.SetTimeBin(meanj-3);
393 AddCluster(c,(Float_t*)vmatrix,k);
397 //unfolding when neccessary
400 Float_t * matrix2[7]; //7x7 matrix with digits - indexing i = 0 ..6 j = -3..3
401 Float_t dummy[7]={0,0,0,0,0,0};
402 for (Int_t di=-3;di<=3;di++){
403 matrix2[di+3] = &bins[k+di*max];
404 if ((k+di*max)<3) matrix2[di+3] = &dummy[3];
405 if ((k+di*max)>fMaxBin-3) matrix2[di+3] = &dummy[3];
407 Float_t vmatrix2[5][5];
410 UnfoldCluster(matrix2,vmatrix2,meani,meanj,sumu,overlap);
412 // c.fMax = FitMax(vmatrix2,meani,meanj,TMath::Sqrt(sigmay2),TMath::Sqrt(sigmaz2));
415 //set cluster parameters
418 c.SetTimeBin(meanj-3);
421 c.SetType(Char_t(overlap)+1);
422 AddCluster(c,(Float_t*)vmatrix,k);
431 void AliTPCclustererMI::UnfoldCluster(Float_t * matrix2[7], Float_t recmatrix[5][5], Float_t & meani, Float_t & meanj,
432 Float_t & sumu, Float_t & overlap )
435 //unfold cluster from input matrix
436 //data corresponding to cluster writen in recmatrix
437 //output meani and meanj
439 //take separatelly y and z
441 Float_t sum3i[7] = {0,0,0,0,0,0,0};
442 Float_t sum3j[7] = {0,0,0,0,0,0,0};
444 for (Int_t k =0;k<7;k++)
445 for (Int_t l = -1; l<=1;l++){
446 sum3i[k]+=matrix2[k][l];
447 sum3j[k]+=matrix2[l+3][k-3];
449 Float_t mratio[3][3]={{1,1,1},{1,1,1},{1,1,1}};
452 Float_t sum3wi = 0; //charge minus overlap
453 Float_t sum3wio = 0; //full charge
454 Float_t sum3iw = 0; //sum for mean value
455 for (Int_t dk=-1;dk<=1;dk++){
456 sum3wio+=sum3i[dk+3];
462 if ( ( ((sum3i[dk+3]+3)/(sum3i[3]-3))+1 < (sum3i[2*dk+3]-3)/(sum3i[dk+3]+3))||
463 (sum3i[dk+3]<=sum3i[2*dk+3] && sum3i[dk+3]>2 )){
464 Float_t xm2 = sum3i[-dk+3];
465 Float_t xm1 = sum3i[+3];
466 Float_t x1 = sum3i[2*dk+3];
467 Float_t x2 = sum3i[3*dk+3];
468 Float_t w11 = TMath::Max((Float_t)(4.*xm1-xm2),(Float_t)0.000001);
469 Float_t w12 = TMath::Max((Float_t)(4 *x1 -x2),(Float_t)0.);
470 ratio = w11/(w11+w12);
471 for (Int_t dl=-1;dl<=1;dl++)
472 mratio[dk+1][dl+1] *= ratio;
474 Float_t amp = sum3i[dk+3]*ratio;
479 meani = sum3iw/sum3wi;
480 Float_t overlapi = (sum3wio-sum3wi)/sum3wio;
485 Float_t sum3wj = 0; //charge minus overlap
486 Float_t sum3wjo = 0; //full charge
487 Float_t sum3jw = 0; //sum for mean value
488 for (Int_t dk=-1;dk<=1;dk++){
489 sum3wjo+=sum3j[dk+3];
495 if ( ( ((sum3j[dk+3]+3)/(sum3j[3]-3))+1 < (sum3j[2*dk+3]-3)/(sum3j[dk+3]+3)) ||
496 (sum3j[dk+3]<=sum3j[2*dk+3] && sum3j[dk+3]>2)){
497 Float_t xm2 = sum3j[-dk+3];
498 Float_t xm1 = sum3j[+3];
499 Float_t x1 = sum3j[2*dk+3];
500 Float_t x2 = sum3j[3*dk+3];
501 Float_t w11 = TMath::Max((Float_t)(4.*xm1-xm2),(Float_t)0.000001);
502 Float_t w12 = TMath::Max((Float_t)(4 *x1 -x2),(Float_t)0.);
503 ratio = w11/(w11+w12);
504 for (Int_t dl=-1;dl<=1;dl++)
505 mratio[dl+1][dk+1] *= ratio;
507 Float_t amp = sum3j[dk+3]*ratio;
512 meanj = sum3jw/sum3wj;
513 Float_t overlapj = (sum3wjo-sum3wj)/sum3wjo;
514 overlap = Int_t(100*TMath::Max(overlapi,overlapj)+3);
515 sumu = (sum3wj+sum3wi)/2.;
518 //if not overlap detected remove everything
519 for (Int_t di =-2; di<=2;di++)
520 for (Int_t dj =-2; dj<=2;dj++){
521 recmatrix[di+2][dj+2] = matrix2[3+di][dj];
525 for (Int_t di =-1; di<=1;di++)
526 for (Int_t dj =-1; dj<=1;dj++){
528 if (mratio[di+1][dj+1]==1){
529 recmatrix[di+2][dj+2] = matrix2[3+di][dj];
530 if (TMath::Abs(di)+TMath::Abs(dj)>1){
531 recmatrix[2*di+2][dj+2] = matrix2[3+2*di][dj];
532 recmatrix[di+2][2*dj+2] = matrix2[3+di][2*dj];
534 recmatrix[2*di+2][2*dj+2] = matrix2[3+2*di][2*dj];
538 //if we have overlap in direction
539 recmatrix[di+2][dj+2] = mratio[di+1][dj+1]* matrix2[3+di][dj];
540 if (TMath::Abs(di)+TMath::Abs(dj)>1){
541 ratio = TMath::Min((Float_t)(recmatrix[di+2][dj+2]/(matrix2[3+0*di][1*dj]+1)),(Float_t)1.);
542 recmatrix[2*di+2][dj+2] = ratio*recmatrix[di+2][dj+2];
544 ratio = TMath::Min((Float_t)(recmatrix[di+2][dj+2]/(matrix2[3+1*di][0*dj]+1)),(Float_t)1.);
545 recmatrix[di+2][2*dj+2] = ratio*recmatrix[di+2][dj+2];
548 ratio = recmatrix[di+2][dj+2]/matrix2[3][0];
549 recmatrix[2*di+2][2*dj+2] = ratio*recmatrix[di+2][dj+2];
557 Float_t AliTPCclustererMI::FitMax(Float_t vmatrix[5][5], Float_t y, Float_t z, Float_t sigmay, Float_t sigmaz)
564 for (Int_t di = -1;di<=1;di++)
565 for (Int_t dj = -1;dj<=1;dj++){
566 if (vmatrix[2+di][2+dj]>2){
567 Float_t teor = TMath::Gaus(di,y,sigmay*1.2)*TMath::Gaus(dj,z,sigmaz*1.2);
568 sumteor += teor*vmatrix[2+di][2+dj];
569 sumamp += vmatrix[2+di][2+dj]*vmatrix[2+di][2+dj];
572 Float_t max = sumamp/sumteor;
576 void AliTPCclustererMI::AddCluster(AliTPCclusterMI &c, Float_t * /*matrix*/, Int_t /*pos*/){
579 // Transform cluster to the rotated global coordinata
580 // Assign labels to the cluster
581 // add the cluster to the array
582 // for more details - See AliTPCTranform::Transform(x,i,0,1)
583 Float_t meani = c.GetPad();
584 Float_t meanj = c.GetTimeBin();
586 Int_t ki = TMath::Nint(meani);
588 if (ki>=fMaxPad) ki = fMaxPad-1;
589 Int_t kj = TMath::Nint(meanj);
591 if (kj>=fMaxTime-3) kj=fMaxTime-4;
592 // ki and kj shifted as integers coordinata
594 c.SetLabel(fRowDig->GetTrackIDFast(kj,ki,0)-2,0);
595 c.SetLabel(fRowDig->GetTrackIDFast(kj,ki,1)-2,1);
596 c.SetLabel(fRowDig->GetTrackIDFast(kj,ki,2)-2,2);
600 c.SetDetector(fSector);
601 Float_t s2 = c.GetSigmaY2();
602 Float_t w=fParam->GetPadPitchWidth(fSector);
603 c.SetSigmaY2(s2*w*w);
605 c.SetSigmaZ2(s2*fZWidth*fZWidth);
610 AliTPCTransform *transform = AliTPCcalibDB::Instance()->GetTransform() ;
612 AliFatal("Tranformations not in calibDB");
615 transform->SetCurrentRecoParam((AliTPCRecoParam*)fRecoParam);
616 Double_t x[3]={c.GetRow(),c.GetPad(),c.GetTimeBin()};
617 Int_t i[1]={fSector};
618 transform->Transform(x,i,0,1);
624 if (ki<=1 || ki>=fMaxPad-1 || kj==1 || kj==fMaxTime-2) {
625 c.SetType(-(c.GetType()+3)); //edge clusters
627 if (fLoop==2) c.SetType(100);
630 TClonesArray * arr = 0;
631 AliTPCclusterMI * cl = 0;
633 if(fBClonesArray==kFALSE) {
634 arr = fRowCl->GetArray();
635 cl = new ((*arr)[fNcluster]) AliTPCclusterMI(c);
637 cl = new ((*fOutputClonesArray)[fNclusters+fNcluster]) AliTPCclusterMI(c);
640 // if (fRecoParam->DumpSignal() &&matrix ) {
642 // Float_t *graph =0;
643 // if (fRecoParam->GetCalcPedestal() && cl->GetMax()>fRecoParam->GetDumpAmplitudeMin() &&fBDumpSignal){
645 // graph = &(fBins[fMaxTime*(pos/fMaxTime)]);
647 // AliTPCclusterInfo * info = new AliTPCclusterInfo(matrix,nbins,graph);
648 // cl->SetInfo(info);
650 if (!fRecoParam->DumpSignal()) {
654 if (AliTPCReconstructor::StreamLevel()>1) {
656 cl->GetGlobalXYZ(xyz);
657 (*fDebugStreamer)<<"Clusters"<<
669 //_____________________________________________________________________________
670 void AliTPCclustererMI::Digits2Clusters()
672 //-----------------------------------------------------------------
673 // This is a simple cluster finder.
674 //-----------------------------------------------------------------
677 Error("Digits2Clusters", "input tree not initialised");
680 fRecoParam = AliTPCReconstructor::GetRecoParam();
682 AliFatal("Can not get the reconstruction parameters");
684 if(AliTPCReconstructor::StreamLevel()>5) {
685 AliInfo("Parameter Dumps");
691 //-----------------------------------------------------------------
693 //-----------------------------------------------------------------
694 if (fUseHLTClusters == 3 || fUseHLTClusters == 4) {
695 AliInfo("Using HLT clusters for TPC off-line reconstruction");
696 fZWidth = fParam->GetZWidth();
697 Int_t iResult = ReadHLTClusters();
699 // HLT clusters present
700 if (iResult >= 0 && fNclusters > 0)
703 // HLT clusters not present
704 if (iResult < 0 || fNclusters == 0) {
705 if (fUseHLTClusters == 3) {
706 AliError("No HLT clusters present, but requested.");
710 AliInfo("Now trying to read from TPC RAW");
713 // Some other problem during cluster reading
715 AliWarning("Some problem while unpacking of HLT clusters.");
718 } // if (fUseHLTClusters == 3 || fUseHLTClusters == 4) {
720 //-----------------------------------------------------------------
721 // Run TPC off-line clusterer
722 //-----------------------------------------------------------------
723 AliTPCCalPad * gainTPC = AliTPCcalibDB::Instance()->GetPadGainFactor();
724 AliTPCCalPad * noiseTPC = AliTPCcalibDB::Instance()->GetPadNoise();
725 AliSimDigits digarr, *dummy=&digarr;
727 fInput->GetBranch("Segment")->SetAddress(&dummy);
728 Stat_t nentries = fInput->GetEntries();
730 fMaxTime=fRecoParam->GetLastBin()+6; // add 3 virtual time bins before and 3 after
734 for (Int_t n=0; n<nentries; n++) {
736 if (!fParam->AdjustSectorRow(digarr.GetID(),fSector,fRow)) {
737 cerr<<"AliTPC warning: invalid segment ID ! "<<digarr.GetID()<<endl;
741 AliTPCCalROC * gainROC = gainTPC->GetCalROC(fSector); // pad gains per given sector
742 AliTPCCalROC * noiseROC = noiseTPC->GetCalROC(fSector); // noise per given sector
745 fRowCl->SetID(digarr.GetID());
746 if (fOutput) fOutput->GetBranch("Segment")->SetAddress(&fRowCl);
747 fRx=fParam->GetPadRowRadii(fSector,row);
750 const Int_t kNIS=fParam->GetNInnerSector(), kNOS=fParam->GetNOuterSector();
751 fZWidth = fParam->GetZWidth();
752 if (fSector < kNIS) {
753 fMaxPad = fParam->GetNPadsLow(row);
754 fSign = (fSector < kNIS/2) ? 1 : -1;
755 fPadLength = fParam->GetPadPitchLength(fSector,row);
756 fPadWidth = fParam->GetPadPitchWidth();
758 fMaxPad = fParam->GetNPadsUp(row);
759 fSign = ((fSector-kNIS) < kNOS/2) ? 1 : -1;
760 fPadLength = fParam->GetPadPitchLength(fSector,row);
761 fPadWidth = fParam->GetPadPitchWidth();
765 fMaxBin=fMaxTime*(fMaxPad+6); // add 3 virtual pads before and 3 after
766 fBins =new Float_t[fMaxBin];
767 fSigBins =new Int_t[fMaxBin];
769 memset(fBins,0,sizeof(Float_t)*fMaxBin);
771 if (digarr.First()) //MI change
773 Float_t dig=digarr.CurrentDigit();
774 if (dig<=fParam->GetZeroSup()) continue;
775 Int_t j=digarr.CurrentRow()+3, i=digarr.CurrentColumn()+3;
776 Float_t gain = gainROC->GetValue(row,digarr.CurrentColumn());
777 Int_t bin = i*fMaxTime+j;
783 fSigBins[fNSigBins++]=bin;
784 } while (digarr.Next());
785 digarr.ExpandTrackBuffer();
787 FindClusters(noiseROC);
789 fRowCl->GetArray()->Clear("C");
790 nclusters+=fNcluster;
796 Info("Digits2Clusters", "Number of found clusters : %d", nclusters);
798 if (fUseHLTClusters == 2 && nclusters == 0) {
799 AliInfo("No clusters from TPC Raw data, now trying to read HLT clusters.");
801 fZWidth = fParam->GetZWidth();
806 void AliTPCclustererMI::ProcessSectorData(){
808 // Process the data for the current sector
811 AliTPCCalPad * pedestalTPC = AliTPCcalibDB::Instance()->GetPedestals();
812 AliTPCCalPad * noiseTPC = AliTPCcalibDB::Instance()->GetPadNoise();
813 AliTPCCalROC * pedestalROC = pedestalTPC->GetCalROC(fSector); // pedestal per given sector
814 AliTPCCalROC * noiseROC = noiseTPC->GetCalROC(fSector); // noise per given sector
815 //check the presence of the calibration
816 if (!noiseROC ||!pedestalROC ) {
817 AliError(Form("Missing calibration per sector\t%d\n",fSector));
820 Int_t nRows=fParam->GetNRow(fSector);
821 Bool_t calcPedestal = fRecoParam->GetCalcPedestal();
822 Int_t zeroSup = fParam->GetZeroSup();
823 // if (calcPedestal) {
825 for (Int_t iRow = 0; iRow < nRows; iRow++) {
826 Int_t maxPad = fParam->GetNPads(fSector, iRow);
828 for (Int_t iPad = 3; iPad < maxPad + 3; iPad++) {
830 // Temporary fix for data production - !!!! MARIAN
831 // The noise calibration should take mean and RMS - currently the Gaussian fit used
832 // In case of double peak - the pad should be rejected
834 // Line mean - if more than given digits over threshold - make a noise calculation
835 // and pedestal substration
836 if (!calcPedestal && fAllBins[iRow][iPad*fMaxTime+0]<50) continue;
838 if (fAllBins[iRow][iPad*fMaxTime+0] <1 ) continue; // no data
839 Float_t *p = &fAllBins[iRow][iPad*fMaxTime+3];
840 //Float_t pedestal = TMath::Median(fMaxTime, p);
841 Int_t id[3] = {fSector, iRow, iPad-3};
843 Double_t rmsCalib= noiseROC->GetValue(iRow,iPad-3);
844 Double_t pedestalCalib = pedestalROC->GetValue(iRow,iPad-3);
845 Double_t rmsEvent = rmsCalib;
846 Double_t pedestalEvent = pedestalCalib;
847 ProcesSignal(p, fMaxTime, id, rmsEvent, pedestalEvent);
848 if (rmsEvent<rmsCalib) rmsEvent = rmsCalib; // take worst scenario
849 if (TMath::Abs(pedestalEvent-pedestalCalib)<1.0) pedestalEvent = pedestalCalib;
852 for (Int_t iTimeBin = 0; iTimeBin < fMaxTime; iTimeBin++) {
853 Int_t bin = iPad*fMaxTime+iTimeBin;
854 fAllBins[iRow][bin] -= pedestalEvent;
855 if (iTimeBin < fRecoParam->GetFirstBin())
856 fAllBins[iRow][bin] = 0;
857 if (iTimeBin > fRecoParam->GetLastBin())
858 fAllBins[iRow][bin] = 0;
859 if (fAllBins[iRow][iPad*fMaxTime+iTimeBin] < zeroSup)
860 fAllBins[iRow][bin] = 0;
861 if (fAllBins[iRow][bin] < 3.0*rmsEvent) // 3 sigma cut on RMS
862 fAllBins[iRow][bin] = 0;
863 if (fAllBins[iRow][bin]) fAllSigBins[iRow][fAllNSigBins[iRow]++] = bin;
869 if (AliTPCReconstructor::StreamLevel()>5) {
870 for (Int_t iRow = 0; iRow < nRows; iRow++) {
871 Int_t maxPad = fParam->GetNPads(fSector,iRow);
873 for (Int_t iPad = 3; iPad < maxPad + 3; iPad++) {
874 for (Int_t iTimeBin = 0; iTimeBin < fMaxTime; iTimeBin++) {
875 Int_t bin = iPad*fMaxTime+iTimeBin;
876 Float_t signal = fAllBins[iRow][bin];
877 if (AliTPCReconstructor::StreamLevel()>3 && signal>3) {
878 Double_t x[]={iRow,iPad-3,iTimeBin-3};
880 AliTPCTransform trafo;
881 trafo.Transform(x,i,0,1);
882 Double_t gx[3]={x[0],x[1],x[2]};
883 trafo.RotatedGlobal2Global(fSector,gx);
884 // fAllSigBins[iRow][fAllNSigBins[iRow]++]
885 Int_t rowsigBins = fAllNSigBins[iRow];
886 Int_t first=fAllSigBins[iRow][0];
888 // if (rowsigBins>0) fAllSigBins[iRow][fAllNSigBins[iRow]-1];
890 if (AliTPCReconstructor::StreamLevel()>5) {
891 (*fDebugStreamer)<<"Digits"<<
904 "rowsigBins="<<rowsigBins<<
915 // Now loop over rows and find clusters
916 for (fRow = 0; fRow < nRows; fRow++) {
917 fRowCl->SetID(fParam->GetIndex(fSector, fRow));
918 if (fOutput) fOutput->GetBranch("Segment")->SetAddress(&fRowCl);
920 fRx = fParam->GetPadRowRadii(fSector, fRow);
921 fPadLength = fParam->GetPadPitchLength(fSector, fRow);
922 fPadWidth = fParam->GetPadPitchWidth();
923 fMaxPad = fParam->GetNPads(fSector,fRow);
924 fMaxBin = fMaxTime*(fMaxPad+6); // add 3 virtual pads before and 3 after
926 fBins = fAllBins[fRow];
927 fSigBins = fAllSigBins[fRow];
928 fNSigBins = fAllNSigBins[fRow];
930 FindClusters(noiseROC);
933 if(fBClonesArray == kFALSE) fRowCl->GetArray()->Clear("C");
934 fNclusters += fNcluster;
936 } // End of loop to find clusters
940 void AliTPCclustererMI::Digits2Clusters(AliRawReader* rawReader)
942 //-----------------------------------------------------------------
943 // This is a cluster finder for the TPC raw data.
944 // The method assumes NO ordering of the altro channels.
945 // The pedestal subtraction can be switched on and off
946 // using an option of the TPC reconstructor
947 //-----------------------------------------------------------------
948 fRecoParam = AliTPCReconstructor::GetRecoParam();
950 AliFatal("Can not get the reconstruction parameters");
952 if(AliTPCReconstructor::StreamLevel()>5) {
953 AliInfo("Parameter Dumps");
959 //-----------------------------------------------------------------
961 //-----------------------------------------------------------------
962 if (fUseHLTClusters == 3 || fUseHLTClusters == 4) {
963 AliInfo("Using HLT clusters for TPC off-line reconstruction");
964 fZWidth = fParam->GetZWidth();
965 Int_t iResult = ReadHLTClusters();
967 // HLT clusters present
968 if (iResult >= 0 && fNclusters > 0)
971 // HLT clusters not present
972 if (iResult < 0 || fNclusters == 0) {
973 if (fUseHLTClusters == 3) {
974 AliError("No HLT clusters present, but requested.");
978 AliInfo("Now trying to read TPC RAW");
981 // Some other problem during cluster reading
983 AliWarning("Some problem while unpacking of HLT clusters.");
986 } // if (fUseHLTClusters == 3 || fUseHLTClusters == 4) {
988 //-----------------------------------------------------------------
989 // Run TPC off-line clusterer
990 //-----------------------------------------------------------------
991 AliTPCCalPad * gainTPC = AliTPCcalibDB::Instance()->GetPadGainFactor();
992 AliTPCAltroMapping** mapping =AliTPCcalibDB::Instance()->GetMapping();
994 AliTPCRawStreamV3 input(rawReader,(AliAltroMapping**)mapping);
995 fEventHeader = (AliRawEventHeaderBase*)rawReader->GetEventHeader();
997 fTimeStamp = fEventHeader->Get("Timestamp");
998 fEventType = fEventHeader->Get("Type");
999 AliTPCTransform *transform = AliTPCcalibDB::Instance()->GetTransform() ;
1000 transform->SetCurrentTimeStamp(fTimeStamp);
1001 transform->SetCurrentRun(rawReader->GetRunNumber());
1004 // creaate one TClonesArray for all clusters
1005 if(fBClonesArray && !fOutputClonesArray) fOutputClonesArray = new TClonesArray("AliTPCclusterMI",1000);
1009 fMaxTime = fRecoParam->GetLastBin() + 6; // add 3 virtual time bins before and 3 after
1010 // const Int_t kNIS = fParam->GetNInnerSector();
1011 // const Int_t kNOS = fParam->GetNOuterSector();
1012 // const Int_t kNS = kNIS + kNOS;
1013 fZWidth = fParam->GetZWidth();
1014 Int_t zeroSup = fParam->GetZeroSup();
1018 AliTPCROC * roc = AliTPCROC::Instance();
1019 Int_t nRowsMax = roc->GetNRows(roc->GetNSector()-1);
1020 Int_t nPadsMax = roc->GetNPads(roc->GetNSector()-1,nRowsMax-1);
1021 for (Int_t iRow = 0; iRow < nRowsMax; iRow++) {
1023 Int_t maxBin = fMaxTime*(nPadsMax+6); // add 3 virtual pads before and 3 after
1024 memset(fAllBins[iRow],0,sizeof(Float_t)*maxBin);
1025 fAllNSigBins[iRow]=0;
1028 Int_t prevSector=-1;
1034 const Int_t kNIS = fParam->GetNInnerSector();
1035 const Int_t kNOS = fParam->GetNOuterSector();
1036 const Int_t kNS = kNIS + kNOS;
1038 for(fSector = 0; fSector < kNS; fSector++) {
1041 Int_t nDDLs = 0, indexDDL = 0;
1042 if (fSector < kNIS) {
1043 nRows = fParam->GetNRowLow();
1044 fSign = (fSector < kNIS/2) ? 1 : -1;
1046 indexDDL = fSector * 2;
1049 nRows = fParam->GetNRowUp();
1050 fSign = ((fSector-kNIS) < kNOS/2) ? 1 : -1;
1052 indexDDL = (fSector-kNIS) * 4 + kNIS * 2;
1055 // load the raw data for corresponding DDLs
1057 rawReader->Select("TPC",indexDDL,indexDDL+nDDLs-1);
1059 while (input.NextDDL()){
1060 if (input.GetSector() != fSector)
1061 AliFatal(Form("Sector index mismatch ! Expected (%d), but got (%d) !",fSector,input.GetSector()));
1063 //Int_t nRows = fParam->GetNRow(fSector);
1065 AliTPCCalROC * gainROC = gainTPC->GetCalROC(fSector); // pad gains per given sector
1066 // Begin loop over altro data
1067 Bool_t calcPedestal = fRecoParam->GetCalcPedestal();
1071 while ( input.NextChannel() ) {
1072 Int_t iRow = input.GetRow();
1077 AliError(Form("Pad-row index (%d) outside the range (%d -> %d) !",
1078 iRow, 0, nRows -1));
1082 Int_t iPad = input.GetPad();
1083 if (iPad < 0 || iPad >= nPadsMax) {
1084 AliError(Form("Pad index (%d) outside the range (%d -> %d) !",
1085 iPad, 0, nPadsMax-1));
1088 gain = gainROC->GetValue(iRow,iPad);
1092 while ( input.NextBunch() ){
1093 Int_t startTbin = (Int_t)input.GetStartTimeBin();
1094 Int_t bunchlength = (Int_t)input.GetBunchLength();
1095 const UShort_t *sig = input.GetSignals();
1096 for (Int_t iTime = 0; iTime<bunchlength; iTime++){
1097 Int_t iTimeBin=startTbin-iTime;
1098 if ( iTimeBin < fRecoParam->GetFirstBin() || iTimeBin >= fRecoParam->GetLastBin()){
1100 AliFatal(Form("Timebin index (%d) outside the range (%d -> %d) !",
1101 iTimeBin, 0, iTimeBin -1));
1105 Float_t signal=(Float_t)sig[iTime];
1106 if (!calcPedestal && signal <= zeroSup) continue;
1108 if (!calcPedestal) {
1109 Int_t bin = iPad*fMaxTime+iTimeBin;
1111 fAllBins[iRow][bin] = signal/gain;
1113 fAllBins[iRow][bin] =0;
1115 fAllSigBins[iRow][fAllNSigBins[iRow]++] = bin;
1117 fAllBins[iRow][iPad*fMaxTime+iTimeBin] = signal;
1119 fAllBins[iRow][iPad*fMaxTime+0]+=1.; // pad with signal
1123 }// end loop signals in bunch
1124 }// end loop bunches
1130 // Now loop over rows and perform pedestal subtraction
1131 if (digCounter==0) continue;
1132 } // End of loop over sectors
1133 //process last sector
1134 if ( digCounter>0 ){
1135 ProcessSectorData();
1136 for (Int_t iRow = 0; iRow < fParam->GetNRow(fSector); iRow++) {
1137 Int_t maxPad = fParam->GetNPads(fSector,iRow);
1138 Int_t maxBin = fMaxTime*(maxPad+6); // add 3 virtual pads before and 3 after
1139 memset(fAllBins[iRow],0,sizeof(Float_t)*maxBin);
1140 fAllNSigBins[iRow] = 0;
1147 if (rawReader->GetEventId() && fOutput ){
1148 Info("Digits2Clusters", "File %s Event\t%u\tNumber of found clusters : %d\n", fOutput->GetName(),*(rawReader->GetEventId()), fNclusters);
1151 if(rawReader->GetEventId()) {
1152 Info("Digits2Clusters", "Event\t%u\tNumber of found clusters : %d\n",*(rawReader->GetEventId()), fNclusters);
1156 //Info("Digits2Clusters", "Number of found clusters : %d\n",fOutputClonesArray->GetEntriesFast());
1159 if (fUseHLTClusters == 2 && fNclusters == 0) {
1160 AliInfo("No clusters from TPC Raw data, now trying to read HLT clusters.");
1162 fZWidth = fParam->GetZWidth();
1171 void AliTPCclustererMI::Digits2ClustersOld
1172 (AliRawReader* rawReader)
1174 //-----------------------------------------------------------------
1175 // This is a cluster finder for the TPC raw data.
1176 // The method assumes NO ordering of the altro channels.
1177 // The pedestal subtraction can be switched on and off
1178 // using an option of the TPC reconstructor
1179 //-----------------------------------------------------------------
1180 fRecoParam = AliTPCReconstructor::GetRecoParam();
1182 AliFatal("Can not get the reconstruction parameters");
1184 if(AliTPCReconstructor::StreamLevel()>5) {
1185 AliInfo("Parameter Dumps");
1191 AliTPCCalPad * gainTPC = AliTPCcalibDB::Instance()->GetPadGainFactor();
1192 AliTPCAltroMapping** mapping =AliTPCcalibDB::Instance()->GetMapping();
1194 AliTPCRawStream input(rawReader,(AliAltroMapping**)mapping);
1195 fEventHeader = (AliRawEventHeaderBase*)rawReader->GetEventHeader();
1197 fTimeStamp = fEventHeader->Get("Timestamp");
1198 fEventType = fEventHeader->Get("Type");
1201 // creaate one TClonesArray for all clusters
1202 if(fBClonesArray && !fOutputClonesArray) fOutputClonesArray = new TClonesArray("AliTPCclusterMI",1000);
1206 fMaxTime = fRecoParam->GetLastBin() + 6; // add 3 virtual time bins before and 3 after
1207 const Int_t kNIS = fParam->GetNInnerSector();
1208 const Int_t kNOS = fParam->GetNOuterSector();
1209 const Int_t kNS = kNIS + kNOS;
1210 fZWidth = fParam->GetZWidth();
1211 Int_t zeroSup = fParam->GetZeroSup();
1216 AliTPCROC * roc = AliTPCROC::Instance();
1217 Int_t nRowsMax = roc->GetNRows(roc->GetNSector()-1);
1218 Int_t nPadsMax = roc->GetNPads(roc->GetNSector()-1,nRowsMax-1);
1219 for (Int_t iRow = 0; iRow < nRowsMax; iRow++) {
1221 Int_t maxBin = fMaxTime*(nPadsMax+6); // add 3 virtual pads before and 3 after
1222 memset(fAllBins[iRow],0,sizeof(Float_t)*maxBin);
1223 fAllNSigBins[iRow]=0;
1226 // Loop over sectors
1228 for(fSector = 0; fSector < kNS; fSector++) {
1231 Int_t nDDLs = 0, indexDDL = 0;
1232 if (fSector < kNIS) {
1233 nRows = fParam->GetNRowLow();
1234 fSign = (fSector < kNIS/2) ? 1 : -1;
1236 indexDDL = fSector * 2;
1239 nRows = fParam->GetNRowUp();
1240 fSign = ((fSector-kNIS) < kNOS/2) ? 1 : -1;
1242 indexDDL = (fSector-kNIS) * 4 + kNIS * 2;
1245 // load the raw data for corresponding DDLs
1247 rawReader->Select("TPC",indexDDL,indexDDL+nDDLs-1);
1249 // select only good sector
1250 if (!input.Next()) continue;
1251 if(input.GetSector() != fSector) continue;
1253 AliTPCCalROC * gainROC = gainTPC->GetCalROC(fSector); // pad gains per given sector
1255 for (Int_t iRow = 0; iRow < nRows; iRow++) {
1258 maxPad = fParam->GetNPadsLow(iRow);
1260 maxPad = fParam->GetNPadsUp(iRow);
1262 Int_t maxBin = fMaxTime*(maxPad+6); // add 3 virtual pads before and 3 after
1263 memset(fAllBins[iRow],0,sizeof(Float_t)*maxBin);
1264 fAllNSigBins[iRow] = 0;
1268 // Begin loop over altro data
1269 Bool_t calcPedestal = fRecoParam->GetCalcPedestal();
1274 while (input.Next()) {
1275 if (input.GetSector() != fSector)
1276 AliFatal(Form("Sector index mismatch ! Expected (%d), but got (%d) !",fSector,input.GetSector()));
1279 Int_t iRow = input.GetRow();
1284 if (iRow < 0 || iRow >= nRows){
1285 AliError(Form("Pad-row index (%d) outside the range (%d -> %d) !",
1286 iRow, 0, nRows -1));
1290 Int_t iPad = input.GetPad();
1291 if (iPad < 0 || iPad >= nPadsMax) {
1292 AliError(Form("Pad index (%d) outside the range (%d -> %d) !",
1293 iPad, 0, nPadsMax-1));
1297 gain = gainROC->GetValue(iRow,iPad);
1302 Int_t iTimeBin = input.GetTime();
1303 if ( iTimeBin < fRecoParam->GetFirstBin() || iTimeBin >= fRecoParam->GetLastBin()){
1305 AliFatal(Form("Timebin index (%d) outside the range (%d -> %d) !",
1306 iTimeBin, 0, iTimeBin -1));
1311 Float_t signal = input.GetSignal();
1312 if (!calcPedestal && signal <= zeroSup) continue;
1314 if (!calcPedestal) {
1315 Int_t bin = iPad*fMaxTime+iTimeBin;
1317 fAllBins[iRow][bin] = signal/gain;
1319 fAllBins[iRow][bin] =0;
1321 fAllSigBins[iRow][fAllNSigBins[iRow]++] = bin;
1323 fAllBins[iRow][iPad*fMaxTime+iTimeBin] = signal;
1325 fAllBins[iRow][iPad*fMaxTime+0]+=1.; // pad with signal
1329 } // End of the loop over altro data
1334 // Now loop over rows and perform pedestal subtraction
1335 if (digCounter==0) continue;
1336 ProcessSectorData();
1337 } // End of loop over sectors
1339 if (rawReader->GetEventId() && fOutput ){
1340 Info("Digits2Clusters", "File %s Event\t%d\tNumber of found clusters : %d\n", fOutput->GetName(),*(rawReader->GetEventId()), fNclusters);
1343 if(rawReader->GetEventId()) {
1344 Info("Digits2Clusters", "Event\t%d\tNumber of found clusters : %d\n",*(rawReader->GetEventId()), fNclusters);
1348 //Info("Digits2Clusters", "Number of found clusters : %d\n",fOutputClonesArray->GetEntriesFast());
1352 void AliTPCclustererMI::FindClusters(AliTPCCalROC * noiseROC)
1356 // add virtual charge at the edge
1358 Double_t kMaxDumpSize = 500000;
1360 fBDumpSignal =kFALSE;
1362 if (fRecoParam->GetCalcPedestal() && fOutput->GetZipBytes()< kMaxDumpSize) fBDumpSignal =kTRUE; //dump signal flag
1367 Int_t crtime = Int_t((fParam->GetZLength(fSector)-fRecoParam->GetCtgRange()*fRx)/fZWidth-fParam->GetNTBinsL1()-5);
1368 Float_t minMaxCutAbs = fRecoParam->GetMinMaxCutAbs();
1369 Float_t minLeftRightCutAbs = fRecoParam->GetMinLeftRightCutAbs();
1370 Float_t minUpDownCutAbs = fRecoParam->GetMinUpDownCutAbs();
1371 Float_t minMaxCutSigma = fRecoParam->GetMinMaxCutSigma();
1372 Float_t minLeftRightCutSigma = fRecoParam->GetMinLeftRightCutSigma();
1373 Float_t minUpDownCutSigma = fRecoParam->GetMinUpDownCutSigma();
1374 Int_t useOnePadCluster = fRecoParam->GetUseOnePadCluster();
1375 for (Int_t iSig = 0; iSig < fNSigBins; iSig++) {
1376 Int_t i = fSigBins[iSig];
1377 if (i%fMaxTime<=crtime) continue;
1378 Float_t *b = &fBins[i];
1380 if (b[0]<minMaxCutAbs) continue; //threshold for maxima
1382 if (useOnePadCluster==0){
1383 if (b[-1]+b[1]+b[-fMaxTime]+b[fMaxTime]<=0) continue; // cut on isolated clusters
1384 if (b[-1]+b[1]<=0) continue; // cut on isolated clusters
1385 if (b[-fMaxTime]+b[fMaxTime]<=0) continue; // cut on isolated clusters
1388 if ((b[0]+b[-1]+b[1])<minUpDownCutAbs) continue; //threshold for up down (TRF)
1389 if ((b[0]+b[-fMaxTime]+b[fMaxTime])<minLeftRightCutAbs) continue; //threshold for left right (PRF)
1390 if (!IsMaximum(*b,fMaxTime,b)) continue;
1392 Float_t noise = noiseROC->GetValue(fRow, i/fMaxTime);
1393 if (noise>fRecoParam->GetMaxNoise()) continue;
1395 if (b[0]<minMaxCutSigma*noise) continue; //threshold form maxima
1396 if ((b[0]+b[-1]+b[1])<minUpDownCutSigma*noise) continue; //threshold for up town TRF
1397 if ((b[0]+b[-fMaxTime]+b[fMaxTime])<minLeftRightCutSigma*noise) continue; //threshold for left right (PRF)
1399 AliTPCclusterMI c; // default cosntruction without info
1401 MakeCluster(i, fMaxTime, fBins, dummy,c);
1407 Bool_t AliTPCclustererMI::AcceptCluster(AliTPCclusterMI *cl){
1409 // Currently hack to filter digital noise (15.06.2008)
1410 // To be parameterized in the AliTPCrecoParam
1411 // More inteligent way to be used in future
1412 // Acces to the proper pedestal file needed
1414 if (cl->GetMax()<400) return kTRUE;
1415 Double_t ratio = cl->GetQ()/cl->GetMax();
1416 if (cl->GetMax()>700){
1417 if ((ratio - int(ratio)>0.8)) return kFALSE;
1419 if ((ratio - int(ratio)<0.95)) return kTRUE;
1424 Double_t AliTPCclustererMI::ProcesSignal(Float_t *signal, Int_t nchannels, Int_t id[3], Double_t &rmsEvent, Double_t &pedestalEvent){
1426 // process signal on given pad - + streaming of additional information in special mode
1433 // ESTIMATE pedestal and the noise
1435 const Int_t kPedMax = 100;
1441 Float_t rmsCalib = rmsEvent; // backup initial value ( from calib)
1442 Float_t pedestalCalib = pedestalEvent;// backup initial value ( from calib)
1443 Int_t firstBin = fRecoParam->GetFirstBin();
1445 UShort_t histo[kPedMax];
1446 //memset(histo,0,kPedMax*sizeof(UShort_t));
1447 for (Int_t i=0; i<kPedMax; i++) histo[i]=0;
1448 for (Int_t i=0; i<fMaxTime; i++){
1449 if (signal[i]<=0) continue;
1450 if (signal[i]>max && i>firstBin) {
1454 if (signal[i]>kPedMax-1) continue;
1455 histo[int(signal[i]+0.5)]++;
1459 for (Int_t i=1; i<kPedMax; i++){
1460 if (count1<count0*0.5) median=i;
1465 Float_t count10=histo[median] ,mean=histo[median]*median, rms=histo[median]*median*median ;
1466 Float_t count06=histo[median] ,mean06=histo[median]*median, rms06=histo[median]*median*median ;
1467 Float_t count09=histo[median] ,mean09=histo[median]*median, rms09=histo[median]*median*median ;
1469 for (Int_t idelta=1; idelta<10; idelta++){
1470 if (median-idelta<=0) continue;
1471 if (median+idelta>kPedMax) continue;
1472 if (count06<0.6*count1){
1473 count06+=histo[median-idelta];
1474 mean06 +=histo[median-idelta]*(median-idelta);
1475 rms06 +=histo[median-idelta]*(median-idelta)*(median-idelta);
1476 count06+=histo[median+idelta];
1477 mean06 +=histo[median+idelta]*(median+idelta);
1478 rms06 +=histo[median+idelta]*(median+idelta)*(median+idelta);
1480 if (count09<0.9*count1){
1481 count09+=histo[median-idelta];
1482 mean09 +=histo[median-idelta]*(median-idelta);
1483 rms09 +=histo[median-idelta]*(median-idelta)*(median-idelta);
1484 count09+=histo[median+idelta];
1485 mean09 +=histo[median+idelta]*(median+idelta);
1486 rms09 +=histo[median+idelta]*(median+idelta)*(median+idelta);
1488 if (count10<0.95*count1){
1489 count10+=histo[median-idelta];
1490 mean +=histo[median-idelta]*(median-idelta);
1491 rms +=histo[median-idelta]*(median-idelta)*(median-idelta);
1492 count10+=histo[median+idelta];
1493 mean +=histo[median+idelta]*(median+idelta);
1494 rms +=histo[median+idelta]*(median+idelta)*(median+idelta);
1499 rms = TMath::Sqrt(TMath::Abs(rms/count10-mean*mean));
1503 rms06 = TMath::Sqrt(TMath::Abs(rms06/count06-mean06*mean06));
1507 rms09 = TMath::Sqrt(TMath::Abs(rms09/count09-mean09*mean09));
1511 pedestalEvent = median;
1512 if (AliLog::GetDebugLevel("","AliTPCclustererMI")==0) return median;
1514 UInt_t uid[3] = {UInt_t(id[0]),UInt_t(id[1]),UInt_t(id[2])};
1516 // Dump mean signal info
1518 if (AliTPCReconstructor::StreamLevel()>0) {
1519 (*fDebugStreamer)<<"Signal"<<
1520 "TimeStamp="<<fTimeStamp<<
1521 "EventType="<<fEventType<<
1535 "RMSCalib="<<rmsCalib<<
1536 "PedCalib="<<pedestalCalib<<
1540 // fill pedestal histogram
1545 Float_t kMin =fRecoParam->GetDumpAmplitudeMin(); // minimal signal to be dumped
1546 Float_t *dsignal = new Float_t[nchannels];
1547 Float_t *dtime = new Float_t[nchannels];
1548 for (Int_t i=0; i<nchannels; i++){
1550 dsignal[i] = signal[i];
1554 // Big signals dumping
1556 if (AliTPCReconstructor::StreamLevel()>0) {
1557 if (max-median>kMin &&maxPos>fRecoParam->GetFirstBin())
1558 (*fDebugStreamer)<<"SignalB"<< // pads with signal
1559 "TimeStamp="<<fTimeStamp<<
1560 "EventType="<<fEventType<<
1564 // "Graph="<<graph<<
1580 if (rms06>fRecoParam->GetMaxNoise()) {
1581 pedestalEvent+=1024.;
1582 return 1024+median; // sign noisy channel in debug mode
1587 Int_t AliTPCclustererMI::ReadHLTClusters()
1590 // read HLT clusters instead of off line custers,
1591 // used in Digits2Clusters
1594 if (!fHLTClusterAccess) {
1596 ROOT::NewFunc_t pNewFunc=NULL;
1598 pCl=TClass::GetClass("AliHLTTPCClusterAccessHLTOUT");
1599 } while (!pCl && gSystem->Load("libAliHLTTPC.so")==0);
1600 if (!pCl || (pNewFunc=pCl->GetNew())==NULL) {
1601 AliError("can not load class description of AliHLTTPCClusterAccessHLTOUT, aborting ...");
1605 void* p=(*pNewFunc)(NULL);
1607 AliError("unable to create instance of AliHLTTPCClusterAccessHLTOUT");
1610 fHLTClusterAccess=reinterpret_cast<TObject*>(p);
1613 TObject* pClusterAccess=fHLTClusterAccess;
1615 const Int_t kNIS = fParam->GetNInnerSector();
1616 const Int_t kNOS = fParam->GetNOuterSector();
1617 const Int_t kNS = kNIS + kNOS;
1620 // make sure that all clusters from the previous event are cleared
1621 pClusterAccess->Clear("event");
1622 for(fSector = 0; fSector < kNS; fSector++) {
1625 TString param("sector="); param+=fSector;
1626 // prepare for next sector
1627 pClusterAccess->Clear("sector");
1628 pClusterAccess->Execute("read", param, &iResult);
1631 AliError("HLT Clusters can not be found");
1634 TObject* pObj=pClusterAccess->FindObject("clusterarray");
1636 AliError("HLT clusters requested, but not cluster array not present");
1640 TObjArray* clusterArray=dynamic_cast<TClonesArray*>(pObj);
1641 if (!clusterArray) {
1642 AliError("HLT cluster array is not of class type TClonesArray");
1646 AliDebug(4,Form("Reading %d clusters from HLT for sector %d", clusterArray->GetEntriesFast(), fSector));
1648 Int_t nClusterSector=0;
1649 Int_t nRows=fParam->GetNRow(fSector);
1651 for (fRow = 0; fRow < nRows; fRow++) {
1652 fRowCl->SetID(fParam->GetIndex(fSector, fRow));
1653 if (fOutput) fOutput->GetBranch("Segment")->SetAddress(&fRowCl);
1654 fNcluster=0; // reset clusters per row
1656 fRx = fParam->GetPadRowRadii(fSector, fRow);
1657 fPadLength = fParam->GetPadPitchLength(fSector, fRow);
1658 fPadWidth = fParam->GetPadPitchWidth();
1659 fMaxPad = fParam->GetNPads(fSector,fRow);
1660 fMaxBin = fMaxTime*(fMaxPad+6); // add 3 virtual pads before and 3 after
1662 fBins = fAllBins[fRow];
1663 fSigBins = fAllSigBins[fRow];
1664 fNSigBins = fAllNSigBins[fRow];
1666 for (Int_t i=0; i<clusterArray->GetEntriesFast(); i++) {
1667 if (!clusterArray->At(i))
1670 AliTPCclusterMI* cluster=dynamic_cast<AliTPCclusterMI*>(clusterArray->At(i));
1671 if (!cluster) continue;
1672 if (cluster->GetRow()!=fRow) continue;
1674 AddCluster(*cluster, NULL, 0);
1678 fRowCl->GetArray()->Clear("c");
1680 } // for (fRow = 0; fRow < nRows; fRow++) {
1681 if (nClusterSector!=clusterArray->GetEntriesFast()) {
1682 AliError(Form("Failed to read %d out of %d HLT clusters",
1683 clusterArray->GetEntriesFast()-nClusterSector,
1684 clusterArray->GetEntriesFast()));
1686 fNclusters+=nClusterSector;
1687 } // for(fSector = 0; fSector < kNS; fSector++) {
1689 pClusterAccess->Clear("event");
1691 Info("Digits2Clusters", "Number of converted HLT clusters : %d", fNclusters);