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 "Riostream.h"
29 #include <TObjArray.h>
32 #include <TTreeStream.h>
33 #include <TVirtualFFT.h>
35 #include "AliDigits.h"
36 #include "AliLoader.h"
38 #include "AliMathBase.h"
39 #include "AliRawEventHeaderBase.h"
40 #include "AliRawReader.h"
41 #include "AliRunLoader.h"
42 #include "AliSimDigits.h"
43 #include "AliTPCCalPad.h"
44 #include "AliTPCCalROC.h"
45 #include "AliTPCClustersArray.h"
46 #include "AliTPCClustersRow.h"
47 #include "AliTPCParam.h"
48 #include "AliTPCRawStream.h"
49 #include "AliTPCRecoParam.h"
50 #include "AliTPCReconstructor.h"
51 #include "AliTPCcalibDB.h"
52 #include "AliTPCclusterInfo.h"
53 #include "AliTPCclusterMI.h"
54 #include "AliTPCTransform.h"
55 #include "AliTPCclustererMI.h"
57 ClassImp(AliTPCclustererMI)
61 AliTPCclustererMI::AliTPCclustererMI(const AliTPCParam* par, const AliTPCRecoParam * recoParam):
76 fPedSubtraction(kFALSE),
77 fIsOldRCUFormat(kFALSE),
95 // param - tpc parameters for given file
96 // recoparam - reconstruction parameters
98 fIsOldRCUFormat = kFALSE;
103 fRecoParam = recoParam;
105 //set default parameters if not specified
106 fRecoParam = AliTPCReconstructor::GetRecoParam();
107 if (!fRecoParam) fRecoParam = AliTPCRecoParam::GetLowFluxParam();
109 fDebugStreamer = new TTreeSRedirector("TPCsignal.root");
111 Int_t nPoints = fRecoParam->GetLastBin()-fRecoParam->GetFirstBin();
112 fFFTr2c = TVirtualFFT::FFT(1, &nPoints, "R2C K");
114 //______________________________________________________________
115 AliTPCclustererMI::AliTPCclustererMI(const AliTPCclustererMI ¶m)
131 fPedSubtraction(kFALSE),
132 fIsOldRCUFormat(kFALSE),
145 fBDumpSignal(kFALSE),
151 fMaxBin = param.fMaxBin;
153 //______________________________________________________________
154 AliTPCclustererMI & AliTPCclustererMI::operator =(const AliTPCclustererMI & param)
157 // assignment operator - dummy
159 fMaxBin=param.fMaxBin;
162 //______________________________________________________________
163 AliTPCclustererMI::~AliTPCclustererMI(){
165 if (fAmplitudeHisto) delete fAmplitudeHisto;
166 if (fDebugStreamer) delete fDebugStreamer;
169 void AliTPCclustererMI::SetInput(TTree * tree)
172 // set input tree with digits
175 if (!fInput->GetBranch("Segment")){
176 cerr<<"AliTPC::Digits2Clusters(): no porper input tree !\n";
182 void AliTPCclustererMI::SetOutput(TTree * tree)
187 AliTPCClustersRow clrow;
188 AliTPCClustersRow *pclrow=&clrow;
189 clrow.SetClass("AliTPCclusterMI");
190 clrow.SetArray(1); // to make Clones array
191 fOutput->Branch("Segment","AliTPCClustersRow",&pclrow,32000,200);
195 Float_t AliTPCclustererMI::GetSigmaY2(Int_t iz){
196 // sigma y2 = in digits - we don't know the angle
197 Float_t z = iz*fParam->GetZWidth()+fParam->GetNTBinsL1()*fParam->GetZWidth();
198 Float_t sd2 = (z*fParam->GetDiffL()*fParam->GetDiffL())/
199 (fPadWidth*fPadWidth);
201 Float_t res = sd2+sres;
206 Float_t AliTPCclustererMI::GetSigmaZ2(Int_t iz){
207 //sigma z2 = in digits - angle estimated supposing vertex constraint
208 Float_t z = iz*fZWidth+fParam->GetNTBinsL1()*fParam->GetZWidth();
209 Float_t sd2 = (z*fParam->GetDiffL()*fParam->GetDiffL())/(fZWidth*fZWidth);
210 Float_t angular = fPadLength*(fParam->GetZLength(fSector)-z)/(fRx*fZWidth);
213 Float_t sres = fParam->GetZSigma()/fZWidth;
215 Float_t res = angular +sd2+sres;
219 void AliTPCclustererMI::MakeCluster(Int_t k,Int_t max,Float_t *bins, UInt_t /*m*/,
223 // k - Make cluster at position k
224 // bins - 2 D array of signals mapped to 1 dimensional array -
225 // max - the number of time bins er one dimension
226 // c - refernce to cluster to be filled
228 Int_t i0=k/max; //central pad
229 Int_t j0=k%max; //central time bin
231 // set pointers to data
232 //Int_t dummy[5] ={0,0,0,0,0};
233 Float_t * matrix[5]; //5x5 matrix with digits - indexing i = 0 ..4 j = -2..2
234 for (Int_t di=-2;di<=2;di++){
235 matrix[di+2] = &bins[k+di*max];
237 //build matrix with virtual charge
238 Float_t sigmay2= GetSigmaY2(j0);
239 Float_t sigmaz2= GetSigmaZ2(j0);
241 Float_t vmatrix[5][5];
242 vmatrix[2][2] = matrix[2][0];
244 c.SetMax((UShort_t)(vmatrix[2][2])); // write maximal amplitude
245 for (Int_t di =-1;di <=1;di++)
246 for (Int_t dj =-1;dj <=1;dj++){
247 Float_t amp = matrix[di+2][dj];
248 if ( (amp<2) && (fLoop<2)){
249 // if under threshold - calculate virtual charge
250 Float_t ratio = TMath::Exp(-1.2*TMath::Abs(di)/sigmay2)*TMath::Exp(-1.2*TMath::Abs(dj)/sigmaz2);
251 amp = ((matrix[2][0]-2)*(matrix[2][0]-2)/(matrix[-di+2][-dj]+2))*ratio;
253 vmatrix[2+di][2+dj]=amp;
254 vmatrix[2+2*di][2+2*dj]=0;
257 vmatrix[2+2*di][2+dj] =0;
258 vmatrix[2+di][2+2*dj] =0;
263 //if small amplitude - below 2 x threshold - don't consider other one
264 vmatrix[2+di][2+dj]=amp;
265 vmatrix[2+2*di][2+2*dj]=0; // don't take to the account next bin
268 vmatrix[2+2*di][2+dj] =0;
269 vmatrix[2+di][2+2*dj] =0;
273 //if bigger then take everything
274 vmatrix[2+di][2+dj]=amp;
275 vmatrix[2+2*di][2+2*dj]= matrix[2*di+2][2*dj] ;
278 vmatrix[2+2*di][2+dj] = matrix[2*di+2][dj];
279 vmatrix[2+di][2+2*dj] = matrix[2+di][dj*2];
291 for (Int_t i=-2;i<=2;i++)
292 for (Int_t j=-2;j<=2;j++){
293 Float_t amp = vmatrix[i+2][j+2];
302 Float_t meani = sumiw/sumw;
303 Float_t mi2 = sumi2w/sumw-meani*meani;
304 Float_t meanj = sumjw/sumw;
305 Float_t mj2 = sumj2w/sumw-meanj*meanj;
307 Float_t ry = mi2/sigmay2;
308 Float_t rz = mj2/sigmaz2;
311 if ( ( (ry<0.6) || (rz<0.6) ) && fLoop==2) return;
312 if ( (ry <1.2) && (rz<1.2) || (!fRecoParam->GetDoUnfold())) {
314 //if cluster looks like expected or Unfolding not switched on
315 //standard COG is used
316 //+1.2 deviation from expected sigma accepted
317 // c.fMax = FitMax(vmatrix,meani,meanj,TMath::Sqrt(sigmay2),TMath::Sqrt(sigmaz2));
321 //set cluster parameters
324 c.SetTimeBin(meanj-2.5);
328 AddCluster(c,(Float_t*)vmatrix,k);
332 //unfolding when neccessary
335 Float_t * matrix2[7]; //7x7 matrix with digits - indexing i = 0 ..6 j = -3..3
336 Float_t dummy[7]={0,0,0,0,0,0};
337 for (Int_t di=-3;di<=3;di++){
338 matrix2[di+3] = &bins[k+di*max];
339 if ((k+di*max)<3) matrix2[di+3] = &dummy[3];
340 if ((k+di*max)>fMaxBin-3) matrix2[di+3] = &dummy[3];
342 Float_t vmatrix2[5][5];
345 UnfoldCluster(matrix2,vmatrix2,meani,meanj,sumu,overlap);
347 // c.fMax = FitMax(vmatrix2,meani,meanj,TMath::Sqrt(sigmay2),TMath::Sqrt(sigmaz2));
350 //set cluster parameters
353 c.SetTimeBin(meanj-3);
356 c.SetType(Char_t(overlap)+1);
357 AddCluster(c,(Float_t*)vmatrix,k);
363 printf("%f\t%f\n", vmatrix2[2][2], vmatrix[2][2]);
368 void AliTPCclustererMI::UnfoldCluster(Float_t * matrix2[7], Float_t recmatrix[5][5], Float_t & meani, Float_t & meanj,
369 Float_t & sumu, Float_t & overlap )
372 //unfold cluster from input matrix
373 //data corresponding to cluster writen in recmatrix
374 //output meani and meanj
376 //take separatelly y and z
378 Float_t sum3i[7] = {0,0,0,0,0,0,0};
379 Float_t sum3j[7] = {0,0,0,0,0,0,0};
381 for (Int_t k =0;k<7;k++)
382 for (Int_t l = -1; l<=1;l++){
383 sum3i[k]+=matrix2[k][l];
384 sum3j[k]+=matrix2[l+3][k-3];
386 Float_t mratio[3][3]={{1,1,1},{1,1,1},{1,1,1}};
389 Float_t sum3wi = 0; //charge minus overlap
390 Float_t sum3wio = 0; //full charge
391 Float_t sum3iw = 0; //sum for mean value
392 for (Int_t dk=-1;dk<=1;dk++){
393 sum3wio+=sum3i[dk+3];
399 if ( ( ((sum3i[dk+3]+3)/(sum3i[3]-3))+1 < (sum3i[2*dk+3]-3)/(sum3i[dk+3]+3))||
400 sum3i[dk+3]<=sum3i[2*dk+3] && sum3i[dk+3]>2 ){
401 Float_t xm2 = sum3i[-dk+3];
402 Float_t xm1 = sum3i[+3];
403 Float_t x1 = sum3i[2*dk+3];
404 Float_t x2 = sum3i[3*dk+3];
405 Float_t w11 = TMath::Max((Float_t)(4.*xm1-xm2),(Float_t)0.000001);
406 Float_t w12 = TMath::Max((Float_t)(4 *x1 -x2),(Float_t)0.);
407 ratio = w11/(w11+w12);
408 for (Int_t dl=-1;dl<=1;dl++)
409 mratio[dk+1][dl+1] *= ratio;
411 Float_t amp = sum3i[dk+3]*ratio;
416 meani = sum3iw/sum3wi;
417 Float_t overlapi = (sum3wio-sum3wi)/sum3wio;
422 Float_t sum3wj = 0; //charge minus overlap
423 Float_t sum3wjo = 0; //full charge
424 Float_t sum3jw = 0; //sum for mean value
425 for (Int_t dk=-1;dk<=1;dk++){
426 sum3wjo+=sum3j[dk+3];
432 if ( ( ((sum3j[dk+3]+3)/(sum3j[3]-3))+1 < (sum3j[2*dk+3]-3)/(sum3j[dk+3]+3)) ||
433 (sum3j[dk+3]<=sum3j[2*dk+3] && sum3j[dk+3]>2)){
434 Float_t xm2 = sum3j[-dk+3];
435 Float_t xm1 = sum3j[+3];
436 Float_t x1 = sum3j[2*dk+3];
437 Float_t x2 = sum3j[3*dk+3];
438 Float_t w11 = TMath::Max((Float_t)(4.*xm1-xm2),(Float_t)0.000001);
439 Float_t w12 = TMath::Max((Float_t)(4 *x1 -x2),(Float_t)0.);
440 ratio = w11/(w11+w12);
441 for (Int_t dl=-1;dl<=1;dl++)
442 mratio[dl+1][dk+1] *= ratio;
444 Float_t amp = sum3j[dk+3]*ratio;
449 meanj = sum3jw/sum3wj;
450 Float_t overlapj = (sum3wjo-sum3wj)/sum3wjo;
451 overlap = Int_t(100*TMath::Max(overlapi,overlapj)+3);
452 sumu = (sum3wj+sum3wi)/2.;
455 //if not overlap detected remove everything
456 for (Int_t di =-2; di<=2;di++)
457 for (Int_t dj =-2; dj<=2;dj++){
458 recmatrix[di+2][dj+2] = matrix2[3+di][dj];
462 for (Int_t di =-1; di<=1;di++)
463 for (Int_t dj =-1; dj<=1;dj++){
465 if (mratio[di+1][dj+1]==1){
466 recmatrix[di+2][dj+2] = matrix2[3+di][dj];
467 if (TMath::Abs(di)+TMath::Abs(dj)>1){
468 recmatrix[2*di+2][dj+2] = matrix2[3+2*di][dj];
469 recmatrix[di+2][2*dj+2] = matrix2[3+di][2*dj];
471 recmatrix[2*di+2][2*dj+2] = matrix2[3+2*di][2*dj];
475 //if we have overlap in direction
476 recmatrix[di+2][dj+2] = mratio[di+1][dj+1]* matrix2[3+di][dj];
477 if (TMath::Abs(di)+TMath::Abs(dj)>1){
478 ratio = TMath::Min((Float_t)(recmatrix[di+2][dj+2]/(matrix2[3+0*di][1*dj]+1)),(Float_t)1.);
479 recmatrix[2*di+2][dj+2] = ratio*recmatrix[di+2][dj+2];
481 ratio = TMath::Min((Float_t)(recmatrix[di+2][dj+2]/(matrix2[3+1*di][0*dj]+1)),(Float_t)1.);
482 recmatrix[di+2][2*dj+2] = ratio*recmatrix[di+2][dj+2];
485 ratio = recmatrix[di+2][dj+2]/matrix2[3][0];
486 recmatrix[2*di+2][2*dj+2] = ratio*recmatrix[di+2][dj+2];
492 printf("%f\n", recmatrix[2][2]);
496 Float_t AliTPCclustererMI::FitMax(Float_t vmatrix[5][5], Float_t y, Float_t z, Float_t sigmay, Float_t sigmaz)
503 for (Int_t di = -1;di<=1;di++)
504 for (Int_t dj = -1;dj<=1;dj++){
505 if (vmatrix[2+di][2+dj]>2){
506 Float_t teor = TMath::Gaus(di,y,sigmay*1.2)*TMath::Gaus(dj,z,sigmaz*1.2);
507 sumteor += teor*vmatrix[2+di][2+dj];
508 sumamp += vmatrix[2+di][2+dj]*vmatrix[2+di][2+dj];
511 Float_t max = sumamp/sumteor;
515 void AliTPCclustererMI::AddCluster(AliTPCclusterMI &c, Float_t * matrix, Int_t pos){
518 // Transform cluster to the rotated global coordinata
519 // Assign labels to the cluster
520 // add the cluster to the array
521 // for more details - See AliTPCTranform::Transform(x,i,0,1)
522 Float_t meani = c.GetPad();
523 Float_t meanj = c.GetTimeBin();
525 Int_t ki = TMath::Nint(meani);
527 if (ki>=fMaxPad) ki = fMaxPad-1;
528 Int_t kj = TMath::Nint(meanj);
530 if (kj>=fMaxTime-3) kj=fMaxTime-4;
531 // ki and kj shifted as integers coordinata
533 c.SetLabel(fRowDig->GetTrackIDFast(kj,ki,0)-2,0);
534 c.SetLabel(fRowDig->GetTrackIDFast(kj,ki,1)-2,1);
535 c.SetLabel(fRowDig->GetTrackIDFast(kj,ki,2)-2,2);
539 c.SetDetector(fSector);
540 Float_t s2 = c.GetSigmaY2();
541 Float_t w=fParam->GetPadPitchWidth(fSector);
542 c.SetSigmaY2(s2*w*w);
544 c.SetSigmaZ2(s2*fZWidth*fZWidth);
548 AliTPCTransform *transform = AliTPCcalibDB::Instance()->GetTransform() ;
550 AliFatal("Tranformations not in calibDB");
552 Double_t x[3]={c.GetRow(),c.GetPad(),c.GetTimeBin()};
553 Int_t i[1]={fSector};
554 transform->Transform(x,i,0,1);
560 if (!fRecoParam->GetBYMirror()){
566 if (ki<=1 || ki>=fMaxPad-1 || kj==1 || kj==fMaxTime-2) {
567 c.SetType(-(c.GetType()+3)); //edge clusters
569 if (fLoop==2) c.SetType(100);
571 TClonesArray * arr = fRowCl->GetArray();
572 AliTPCclusterMI * cl = new ((*arr)[fNcluster]) AliTPCclusterMI(c);
573 if (fRecoParam->DumpSignal() &&matrix ) {
576 if (fRecoParam->GetCalcPedestal() && cl->GetMax()>fRecoParam->GetDumpAmplitudeMin() &&fBDumpSignal){
578 graph = &(fBins[fMaxTime*(pos/fMaxTime)]);
580 AliTPCclusterInfo * info = new AliTPCclusterInfo(matrix,nbins,graph);
583 if (!fRecoParam->DumpSignal()) {
591 //_____________________________________________________________________________
592 void AliTPCclustererMI::Digits2Clusters()
594 //-----------------------------------------------------------------
595 // This is a simple cluster finder.
596 //-----------------------------------------------------------------
599 Error("Digits2Clusters", "input tree not initialised");
604 Error("Digits2Clusters", "output tree not initialised");
608 AliTPCCalPad * gainTPC = AliTPCcalibDB::Instance()->GetPadGainFactor();
609 AliTPCCalPad * noiseTPC = AliTPCcalibDB::Instance()->GetPadNoise();
610 AliSimDigits digarr, *dummy=&digarr;
612 fInput->GetBranch("Segment")->SetAddress(&dummy);
613 Stat_t nentries = fInput->GetEntries();
615 fMaxTime=fParam->GetMaxTBin()+6; // add 3 virtual time bins before and 3 after
619 for (Int_t n=0; n<nentries; n++) {
621 if (!fParam->AdjustSectorRow(digarr.GetID(),fSector,fRow)) {
622 cerr<<"AliTPC warning: invalid segment ID ! "<<digarr.GetID()<<endl;
626 AliTPCCalROC * gainROC = gainTPC->GetCalROC(fSector); // pad gains per given sector
627 AliTPCCalROC * noiseROC = noiseTPC->GetCalROC(fSector); // noise per given sector
629 AliTPCClustersRow *clrow= new AliTPCClustersRow();
631 clrow->SetClass("AliTPCclusterMI");
634 clrow->SetID(digarr.GetID());
635 fOutput->GetBranch("Segment")->SetAddress(&clrow);
636 fRx=fParam->GetPadRowRadii(fSector,row);
639 const Int_t kNIS=fParam->GetNInnerSector(), kNOS=fParam->GetNOuterSector();
640 fZWidth = fParam->GetZWidth();
641 if (fSector < kNIS) {
642 fMaxPad = fParam->GetNPadsLow(row);
643 fSign = (fSector < kNIS/2) ? 1 : -1;
644 fPadLength = fParam->GetPadPitchLength(fSector,row);
645 fPadWidth = fParam->GetPadPitchWidth();
647 fMaxPad = fParam->GetNPadsUp(row);
648 fSign = ((fSector-kNIS) < kNOS/2) ? 1 : -1;
649 fPadLength = fParam->GetPadPitchLength(fSector,row);
650 fPadWidth = fParam->GetPadPitchWidth();
654 fMaxBin=fMaxTime*(fMaxPad+6); // add 3 virtual pads before and 3 after
655 fBins =new Float_t[fMaxBin];
656 fSigBins =new Int_t[fMaxBin];
658 memset(fBins,0,sizeof(Float_t)*fMaxBin);
660 if (digarr.First()) //MI change
662 Float_t dig=digarr.CurrentDigit();
663 if (dig<=fParam->GetZeroSup()) continue;
664 Int_t j=digarr.CurrentRow()+3, i=digarr.CurrentColumn()+3;
665 Float_t gain = gainROC->GetValue(row,digarr.CurrentColumn());
666 Int_t bin = i*fMaxTime+j;
668 fSigBins[fNSigBins++]=bin;
669 } while (digarr.Next());
670 digarr.ExpandTrackBuffer();
672 FindClusters(noiseROC);
676 nclusters+=fNcluster;
681 Info("Digits2Clusters", "Number of found clusters : %d", nclusters);
684 void AliTPCclustererMI::Digits2Clusters(AliRawReader* rawReader)
686 //-----------------------------------------------------------------
687 // This is a cluster finder for the TPC raw data.
688 // The method assumes NO ordering of the altro channels.
689 // The pedestal subtraction can be switched on and off
690 // using an option of the TPC reconstructor
691 //-----------------------------------------------------------------
694 Error("Digits2Clusters", "output tree not initialised");
699 AliTPCROC * roc = AliTPCROC::Instance();
700 AliTPCCalPad * gainTPC = AliTPCcalibDB::Instance()->GetPadGainFactor();
701 AliTPCCalPad * pedestalTPC = AliTPCcalibDB::Instance()->GetPedestals();
702 AliTPCCalPad * noiseTPC = AliTPCcalibDB::Instance()->GetPadNoise();
703 AliTPCAltroMapping** mapping =AliTPCcalibDB::Instance()->GetMapping();
705 AliTPCRawStream input(rawReader,(AliAltroMapping**)mapping);
706 fEventHeader = (AliRawEventHeaderBase*)rawReader->GetEventHeader();
708 fTimeStamp = fEventHeader->Get("Timestamp");
709 fEventType = fEventHeader->Get("Type");
715 fMaxTime = fParam->GetMaxTBin() + 6; // add 3 virtual time bins before and 3 after
716 const Int_t kNIS = fParam->GetNInnerSector();
717 const Int_t kNOS = fParam->GetNOuterSector();
718 const Int_t kNS = kNIS + kNOS;
719 fZWidth = fParam->GetZWidth();
720 Int_t zeroSup = fParam->GetZeroSup();
722 //alocate memory for sector - maximal case
724 Float_t** allBins = NULL;
725 Int_t** allSigBins = NULL;
726 Int_t* allNSigBins = NULL;
727 Int_t nRowsMax = roc->GetNRows(roc->GetNSector()-1);
728 Int_t nPadsMax = roc->GetNPads(roc->GetNSector()-1,nRowsMax-1);
729 allBins = new Float_t*[nRowsMax];
730 allSigBins = new Int_t*[nRowsMax];
731 allNSigBins = new Int_t[nRowsMax];
732 for (Int_t iRow = 0; iRow < nRowsMax; iRow++) {
734 Int_t maxBin = fMaxTime*(nPadsMax+6); // add 3 virtual pads before and 3 after
735 allBins[iRow] = new Float_t[maxBin];
736 memset(allBins[iRow],0,sizeof(Float_t)*maxBin);
737 allSigBins[iRow] = new Int_t[maxBin];
743 for(fSector = 0; fSector < kNS; fSector++) {
745 AliTPCCalROC * gainROC = gainTPC->GetCalROC(fSector); // pad gains per given sector
746 AliTPCCalROC * pedestalROC = pedestalTPC->GetCalROC(fSector); // pedestal per given sector
747 AliTPCCalROC * noiseROC = noiseTPC->GetCalROC(fSector); // noise per given sector
750 Int_t nDDLs = 0, indexDDL = 0;
751 if (fSector < kNIS) {
752 nRows = fParam->GetNRowLow();
753 fSign = (fSector < kNIS/2) ? 1 : -1;
755 indexDDL = fSector * 2;
758 nRows = fParam->GetNRowUp();
759 fSign = ((fSector-kNIS) < kNOS/2) ? 1 : -1;
761 indexDDL = (fSector-kNIS) * 4 + kNIS * 2;
764 for (Int_t iRow = 0; iRow < nRows; iRow++) {
767 maxPad = fParam->GetNPadsLow(iRow);
769 maxPad = fParam->GetNPadsUp(iRow);
771 Int_t maxBin = fMaxTime*(maxPad+6); // add 3 virtual pads before and 3 after
772 memset(allBins[iRow],0,sizeof(Float_t)*maxBin);
773 allNSigBins[iRow] = 0;
776 // Loas the raw data for corresponding DDLs
778 input.SetOldRCUFormat(fIsOldRCUFormat);
779 rawReader->Select("TPC",indexDDL,indexDDL+nDDLs-1);
781 // Begin loop over altro data
782 Bool_t calcPedestal = fRecoParam->GetCalcPedestal();
785 while (input.Next()) {
787 if (input.GetSector() != fSector)
788 AliFatal(Form("Sector index mismatch ! Expected (%d), but got (%d) !",fSector,input.GetSector()));
791 Int_t iRow = input.GetRow();
792 if (iRow < 0 || iRow >= nRows)
793 AliFatal(Form("Pad-row index (%d) outside the range (%d -> %d) !",
796 Int_t iPad = input.GetPad();
797 if (iPad < 0 || iPad >= nPadsMax)
798 AliFatal(Form("Pad index (%d) outside the range (%d -> %d) !",
799 iPad, 0, nPadsMax-1));
801 gain = gainROC->GetValue(iRow,iPad);
806 Int_t iTimeBin = input.GetTime();
807 if ( iTimeBin < 0 || iTimeBin >= fParam->GetMaxTBin())
808 AliFatal(Form("Timebin index (%d) outside the range (%d -> %d) !",
809 iTimeBin, 0, iTimeBin -1));
812 Float_t signal = input.GetSignal();
813 if (!calcPedestal && signal <= zeroSup) continue;
815 Int_t bin = iPad*fMaxTime+iTimeBin;
816 allBins[iRow][bin] = signal/gain;
817 allSigBins[iRow][allNSigBins[iRow]++] = bin;
819 allBins[iRow][iPad*fMaxTime+iTimeBin] = signal;
821 allBins[iRow][iPad*fMaxTime+0]=1.; // pad with signal
822 } // End of the loop over altro data
825 // Now loop over rows and perform pedestal subtraction
826 if (digCounter==0) continue;
827 // if (fPedSubtraction) {
829 for (Int_t iRow = 0; iRow < nRows; iRow++) {
832 maxPad = fParam->GetNPadsLow(iRow);
834 maxPad = fParam->GetNPadsUp(iRow);
836 for (Int_t iPad = 3; iPad < maxPad + 3; iPad++) {
837 if (allBins[iRow][iPad*fMaxTime+0] <1 ) continue; // no data
838 Float_t *p = &allBins[iRow][iPad*fMaxTime+3];
839 //Float_t pedestal = TMath::Median(fMaxTime, p);
840 Int_t id[3] = {fSector, iRow, iPad-3};
842 Double_t rmsCalib= noiseROC->GetValue(iRow,iPad-3);
843 Double_t pedestalCalib = pedestalROC->GetValue(iRow,iPad-3);
844 Double_t rmsEvent = rmsCalib;
845 Double_t pedestalEvent = pedestalCalib;
846 ProcesSignal(p, fMaxTime, id, rmsEvent, pedestalEvent);
847 if (rmsEvent<rmsCalib) rmsEvent = rmsCalib; // take worst scenario
848 if (TMath::Abs(pedestalEvent-pedestalCalib)<1.0) pedestalEvent = pedestalCalib;
851 for (Int_t iTimeBin = 0; iTimeBin < fMaxTime; iTimeBin++) {
852 Int_t bin = iPad*fMaxTime+iTimeBin;
853 allBins[iRow][bin] -= pedestalEvent;
854 if (iTimeBin < AliTPCReconstructor::GetRecoParam()->GetFirstBin())
855 allBins[iRow][bin] = 0;
856 if (iTimeBin > AliTPCReconstructor::GetRecoParam()->GetLastBin())
857 allBins[iRow][bin] = 0;
858 if (allBins[iRow][iPad*fMaxTime+iTimeBin] < zeroSup)
859 allBins[iRow][bin] = 0;
860 if (allBins[iRow][bin] < 3.0*rmsEvent) // 3 sigma cut on RMS
861 allBins[iRow][bin] = 0;
862 if (allBins[iRow][bin]) allSigBins[iRow][allNSigBins[iRow]++] = bin;
867 // Now loop over rows and find clusters
868 for (fRow = 0; fRow < nRows; fRow++) {
869 fRowCl = new AliTPCClustersRow;
870 fRowCl->SetClass("AliTPCclusterMI");
872 fRowCl->SetID(fParam->GetIndex(fSector, fRow));
873 fOutput->GetBranch("Segment")->SetAddress(&fRowCl);
875 fRx = fParam->GetPadRowRadii(fSector, fRow);
876 fPadLength = fParam->GetPadPitchLength(fSector, fRow);
877 fPadWidth = fParam->GetPadPitchWidth();
879 fMaxPad = fParam->GetNPadsLow(fRow);
881 fMaxPad = fParam->GetNPadsUp(fRow);
882 fMaxBin = fMaxTime*(fMaxPad+6); // add 3 virtual pads before and 3 after
884 fBins = allBins[fRow];
885 fSigBins = allSigBins[fRow];
886 fNSigBins = allNSigBins[fRow];
888 FindClusters(noiseROC);
892 nclusters += fNcluster;
893 } // End of loop to find clusters
894 } // End of loop over sectors
896 for (Int_t iRow = 0; iRow < nRowsMax; iRow++) {
897 delete [] allBins[iRow];
898 delete [] allSigBins[iRow];
901 delete [] allSigBins;
902 delete [] allNSigBins;
904 Info("Digits2Clusters", "File %s Event\t%d\tNumber of found clusters : %d\n", fOutput->GetName(),*(rawReader->GetEventId()), nclusters);
908 void AliTPCclustererMI::FindClusters(AliTPCCalROC * noiseROC)
912 // add virtual charge at the edge
914 Double_t kMaxDumpSize = 500000;
915 if (fRecoParam->GetCalcPedestal() && fOutput->GetZipBytes()< kMaxDumpSize) fBDumpSignal =kTRUE; //dump signal flag
919 Int_t crtime = Int_t((fParam->GetZLength(fSector)-fRecoParam->GetCtgRange()*fRx)/fZWidth-fParam->GetNTBinsL1()-5);
920 Float_t minMaxCutAbs = fRecoParam->GetMinMaxCutAbs();
921 Float_t minLeftRightCutAbs = fRecoParam->GetMinLeftRightCutAbs();
922 Float_t minUpDownCutAbs = fRecoParam->GetMinUpDownCutAbs();
923 Float_t minMaxCutSigma = fRecoParam->GetMinMaxCutSigma();
924 Float_t minLeftRightCutSigma = fRecoParam->GetMinLeftRightCutSigma();
925 Float_t minUpDownCutSigma = fRecoParam->GetMinUpDownCutSigma();
926 for (Int_t iSig = 0; iSig < fNSigBins; iSig++) {
927 Int_t i = fSigBins[iSig];
928 if (i%fMaxTime<=crtime) continue;
929 Float_t *b = &fBins[i];
931 if (b[0]<minMaxCutAbs) continue; //threshold for maxima
933 if (b[-1]+b[1]+b[-fMaxTime]+b[fMaxTime]<=0) continue; // cut on isolated clusters
934 // if (b[-1]+b[1]<=0) continue; // cut on isolated clusters
935 //if (b[-fMaxTime]+b[fMaxTime]<=0) continue; // cut on isolated clusters
937 if ((b[0]+b[-1]+b[1])<minUpDownCutAbs) continue; //threshold for up down (TRF)
938 if ((b[0]+b[-fMaxTime]+b[fMaxTime])<minLeftRightCutAbs) continue; //threshold for left right (PRF)
939 if (!IsMaximum(*b,fMaxTime,b)) continue;
941 Float_t noise = noiseROC->GetValue(fRow, i/fMaxTime);
943 if (b[0]<minMaxCutSigma*noise) continue; //threshold form maxima
944 if ((b[0]+b[-1]+b[1])<minUpDownCutSigma*noise) continue; //threshold for up town TRF
945 if ((b[0]+b[-fMaxTime]+b[fMaxTime])<minLeftRightCutSigma*noise) continue; //threshold for left right (PRF)
947 AliTPCclusterMI c(kFALSE); // default cosntruction without info
949 MakeCluster(i, fMaxTime, fBins, dummy,c);
956 Double_t AliTPCclustererMI::ProcesSignal(Float_t *signal, Int_t nchannels, Int_t id[3], Double_t &rmsEvent, Double_t &pedestalEvent){
958 // process signal on given pad - + streaming of additional information in special mode
965 // ESTIMATE pedestal and the noise
967 const Int_t kPedMax = 100;
968 Double_t kMaxDebugSize = 5000000.;
974 Float_t rmsCalib = rmsEvent; // backup initial value ( from calib)
975 Float_t pedestalCalib = pedestalEvent;// backup initial value ( from calib)
976 Int_t firstBin = AliTPCReconstructor::GetRecoParam()->GetFirstBin();
978 UShort_t histo[kPedMax];
979 memset(histo,0,kPedMax*sizeof(UShort_t));
980 for (Int_t i=0; i<fMaxTime; i++){
981 if (signal[i]<=0) continue;
982 if (signal[i]>max && i>firstBin) {
986 if (signal[i]>kPedMax-1) continue;
987 histo[int(signal[i]+0.5)]++;
991 for (Int_t i=1; i<kPedMax; i++){
992 if (count1<count0*0.5) median=i;
997 Float_t count10=histo[median] ,mean=histo[median]*median, rms=histo[median]*median*median ;
998 Float_t count06=histo[median] ,mean06=histo[median]*median, rms06=histo[median]*median*median ;
999 Float_t count09=histo[median] ,mean09=histo[median]*median, rms09=histo[median]*median*median ;
1001 for (Int_t idelta=1; idelta<10; idelta++){
1002 if (median-idelta<=0) continue;
1003 if (median+idelta>kPedMax) continue;
1004 if (count06<0.6*count1){
1005 count06+=histo[median-idelta];
1006 mean06 +=histo[median-idelta]*(median-idelta);
1007 rms06 +=histo[median-idelta]*(median-idelta)*(median-idelta);
1008 count06+=histo[median+idelta];
1009 mean06 +=histo[median+idelta]*(median+idelta);
1010 rms06 +=histo[median+idelta]*(median+idelta)*(median+idelta);
1012 if (count09<0.9*count1){
1013 count09+=histo[median-idelta];
1014 mean09 +=histo[median-idelta]*(median-idelta);
1015 rms09 +=histo[median-idelta]*(median-idelta)*(median-idelta);
1016 count09+=histo[median+idelta];
1017 mean09 +=histo[median+idelta]*(median+idelta);
1018 rms09 +=histo[median+idelta]*(median+idelta)*(median+idelta);
1020 if (count10<0.95*count1){
1021 count10+=histo[median-idelta];
1022 mean +=histo[median-idelta]*(median-idelta);
1023 rms +=histo[median-idelta]*(median-idelta)*(median-idelta);
1024 count10+=histo[median+idelta];
1025 mean +=histo[median+idelta]*(median+idelta);
1026 rms +=histo[median+idelta]*(median+idelta)*(median+idelta);
1032 rms = TMath::Sqrt(TMath::Abs(rms/count10-mean*mean));
1033 rms06 = TMath::Sqrt(TMath::Abs(rms06/count06-mean06*mean06));
1034 rms09 = TMath::Sqrt(TMath::Abs(rms09/count09-mean09*mean09));
1037 pedestalEvent = median;
1038 if (AliLog::GetDebugLevel("","AliTPCclustererMI")==0) return median;
1040 UInt_t uid[3] = {UInt_t(id[0]),UInt_t(id[1]),UInt_t(id[2])};
1042 // Dump mean signal info
1044 (*fDebugStreamer)<<"Signal"<<
1045 "TimeStamp="<<fTimeStamp<<
1046 "EventType="<<fEventType<<
1060 "RMSCalib="<<rmsCalib<<
1061 "PedCalib="<<pedestalCalib<<
1064 // fill pedestal histogram
1066 AliTPCROC * roc = AliTPCROC::Instance();
1067 if (!fAmplitudeHisto){
1068 fAmplitudeHisto = new TObjArray(72);
1071 if (uid[0]<roc->GetNSectors()
1072 && uid[1]< roc->GetNRows(uid[0]) &&
1073 uid[2] <roc->GetNPads(uid[0], uid[1])){
1074 TObjArray * sectorArray = (TObjArray*)fAmplitudeHisto->UncheckedAt(uid[0]);
1076 Int_t npads =roc->GetNChannels(uid[0]);
1077 sectorArray = new TObjArray(npads);
1078 fAmplitudeHisto->AddAt(sectorArray, uid[0]);
1080 Int_t position = uid[2]+roc->GetRowIndexes(uid[0])[uid[1]];
1081 TH1F * histo = (TH1F*)sectorArray->UncheckedAt(position);
1084 sprintf(hname,"Amp_%d_%d_%d",uid[0],uid[1],uid[2]);
1085 TFile * backup = gFile;
1086 fDebugStreamer->GetFile()->cd();
1087 histo = new TH1F(hname, hname, 100, 5,100);
1088 //histo->SetDirectory(0); // histogram not connected to directory -(File)
1089 sectorArray->AddAt(histo, position);
1090 if (backup) backup->cd();
1092 for (Int_t i=0; i<nchannels; i++){
1093 histo->Fill(signal[i]);
1099 Float_t kMin =fRecoParam->GetDumpAmplitudeMin(); // minimal signal to be dumped
1100 Float_t *dsignal = new Float_t[nchannels];
1101 Float_t *dtime = new Float_t[nchannels];
1102 for (Int_t i=0; i<nchannels; i++){
1104 dsignal[i] = signal[i];
1109 // if (max-median>30.*TMath::Max(1.,Double_t(rms06)) && (((*fDebugStreamer)<<"SignalDN").GetSize()<kMaxDebugSize)){
1112 // TGraph * graph =new TGraph(nchannels, dtime, dsignal);
1115 // // jumps left - right
1117 // Double_t deltaT0[2000];
1118 // Double_t deltaA0[2000];
1119 // Int_t lastJump0 = fRecoParam->GetFirstBin();
1121 // Double_t deltaT1[2000];
1122 // Double_t deltaA1[2000];
1123 // Int_t lastJump1 = fRecoParam->GetFirstBin();
1125 // Double_t deltaT2[2000];
1126 // Double_t deltaA2[2000];
1127 // Int_t lastJump2 = fRecoParam->GetFirstBin();
1129 // for (Int_t itime=fRecoParam->GetFirstBin()+1; itime<fRecoParam->GetLastBin()-1; itime++){
1130 // if (TMath::Abs(dsignal[itime]-dsignal[itime-1])>30.*TMath::Max(1.,Double_t(rms06)) &&
1131 // TMath::Abs(dsignal[itime]-dsignal[itime+1])>30.*TMath::Max(1.,Double_t(rms06)) &&
1132 // (dsignal[itime-1]-median<5.*rms06) &&
1133 // (dsignal[itime+1]-median<5.*rms06)
1135 // deltaA0[njumps0] = dsignal[itime]-dsignal[itime-1];
1136 // deltaT0[njumps0] = itime-lastJump0;
1137 // lastJump0 = itime;
1140 // if (TMath::Abs(dsignal[itime]-dsignal[itime-1])>30.*TMath::Max(1.,Double_t(rms06)) &&
1141 // (dsignal[itime-1]-median<5.*rms06)
1143 // deltaA1[njumps1] = dsignal[itime]-dsignal[itime-1];
1144 // deltaT1[njumps1] = itime-lastJump1;
1145 // lastJump1 = itime;
1148 // if (TMath::Abs(dsignal[itime]-dsignal[itime+1])>30.*TMath::Max(1.,Double_t(rms06)) &&
1149 // (dsignal[itime+1]-median<5.*rms06)
1151 // deltaA2[njumps2] = dsignal[itime]-dsignal[itime+1];
1152 // deltaT2[njumps2] = itime-lastJump2;
1153 // lastJump2 = itime;
1158 // if (njumps0>0 || njumps1>0 || njumps2>0){
1159 // TGraph *graphDN0 = new TGraph(njumps0, deltaT0, deltaA0);
1160 // TGraph *graphDN1 = new TGraph(njumps1, deltaT1, deltaA1);
1161 // TGraph *graphDN2 = new TGraph(njumps2, deltaT2, deltaA2);
1162 // (*fDebugStreamer)<<"SignalDN"<< //digital - noise pads - or random sample of pads
1163 // "TimeStamp="<<fTimeStamp<<
1164 // "EventType="<<fEventType<<
1165 // "Sector="<<uid[0]<<
1168 // "Graph="<<graph<<
1170 // "MaxPos="<<maxPos<<
1171 // "Graph.="<<graph<<
1172 // "P0GraphDN0.="<<graphDN0<<
1173 // "P1GraphDN1.="<<graphDN1<<
1174 // "P2GraphDN2.="<<graphDN2<<
1176 // "Median="<<median<<
1179 // "Mean06="<<mean06<<
1180 // "RMS06="<<rms06<<
1181 // "Mean09="<<mean09<<
1182 // "RMS09="<<rms09<<
1192 // NOISE STUDY Fourier transform
1195 Bool_t random = (gRandom->Rndm()<0.0003);
1196 if (((*fDebugStreamer)<<"SignalN").GetSize()<kMaxDebugSize)
1197 if (max-median>kMin || rms06>1.*fParam->GetZeroSup() || random){
1198 graph =new TGraph(nchannels, dtime, dsignal);
1199 if (rms06>1.*fParam->GetZeroSup() || random){
1200 //Double_t *input, Double_t threshold, Bool_t locMax, Double_t *freq, Double_t *re, Double_t *im, Double_t *mag, Double_t *phi);
1201 Float_t * input = &(dsignal[fRecoParam->GetFirstBin()]);
1202 Float_t freq[2000], re[2000], im[2000], mag[2000], phi[2000];
1203 Int_t npoints = TransformFFT(input, -1,kFALSE, freq, re, im, mag, phi);
1204 TGraph *graphMag0 = new TGraph(npoints, freq, mag);
1205 TGraph *graphPhi0 = new TGraph(npoints, freq, phi);
1206 npoints = TransformFFT(input, 0.5,kTRUE, freq, re, im, mag, phi);
1207 TGraph *graphMag1 = new TGraph(npoints, freq, mag);
1208 TGraph *graphPhi1 = new TGraph(npoints, freq, phi);
1210 (*fDebugStreamer)<<"SignalN"<< //noise pads - or random sample of pads
1211 "TimeStamp="<<fTimeStamp<<
1212 "EventType="<<fEventType<<
1228 "Mag0.="<<graphMag0<<
1229 "Mag1.="<<graphMag1<<
1230 "Phi0.="<<graphPhi0<<
1231 "Phi1.="<<graphPhi1<<
1239 // Big signals dumping
1242 if (max-median>kMin &&maxPos>AliTPCReconstructor::GetRecoParam()->GetFirstBin())
1243 (*fDebugStreamer)<<"SignalB"<< // pads with signal
1244 "TimeStamp="<<fTimeStamp<<
1245 "EventType="<<fEventType<<
1266 // Central Electrode signal analysis
1268 Float_t ceQmax =0, ceQsum=0, ceTime=0;
1269 Float_t cemean = mean06, cerms=rms06 ;
1271 Float_t ceThreshold=5.*cerms;
1272 Float_t ceSumThreshold=8.*cerms;
1273 const Int_t kCemin=5; // range for the analysis of the ce signal +- channels from the peak
1274 const Int_t kCemax=5;
1275 for (Int_t i=nchannels-2; i>nchannels/2; i--){
1276 if ( (dsignal[i]-mean06)>ceThreshold && dsignal[i]>=dsignal[i+1] && dsignal[i]>=dsignal[i-1] ){
1284 for (Int_t i=cemaxpos-20; i<cemaxpos+5; i++){
1285 if (i<0 || i>nchannels-1) continue;
1286 Double_t val=dsignal[i]- cemean;
1292 cemaxpos = cemaxpos2;
1294 for (Int_t i=cemaxpos-kCemin; i<cemaxpos+kCemax; i++){
1295 if (i>0 && i<nchannels&&dsignal[i]- cemean>0){
1296 Double_t val=dsignal[i]- cemean;
1297 ceTime+=val*dtime[i];
1299 if (val>ceQmax) ceQmax=val;
1302 if (ceQmax&&ceQsum>ceSumThreshold) {
1304 (*fDebugStreamer)<<"Signalce"<<
1305 "TimeStamp="<<fTimeStamp<<
1306 "EventType="<<fEventType<<
1318 // end of ce signal analysis
1322 // Gating grid signal analysis
1324 Double_t ggQmax =0, ggQsum=0, ggTime=0;
1325 Double_t ggmean = mean06, ggrms=rms06 ;
1327 Double_t ggThreshold=5.*ggrms;
1328 Double_t ggSumThreshold=8.*ggrms;
1330 for (Int_t i=1; i<nchannels/4; i++){
1331 if ( (dsignal[i]-mean06)>ggThreshold && dsignal[i]>=dsignal[i+1] && dsignal[i]>=dsignal[i-1] &&
1332 (dsignal[i]+dsignal[i+1]+dsignal[i-1]-3*mean06)>ggSumThreshold){
1334 if (dsignal[i-1]>dsignal[i+1]) ggmaxpos=i-1;
1339 for (Int_t i=ggmaxpos-1; i<ggmaxpos+3; i++){
1340 if (i>0 && i<nchannels && dsignal[i]-ggmean>0){
1341 Double_t val=dsignal[i]- ggmean;
1342 ggTime+=val*dtime[i];
1344 if (val>ggQmax) ggQmax=val;
1347 if (ggQmax&&ggQsum>ggSumThreshold) {
1349 (*fDebugStreamer)<<"Signalgg"<<
1350 "TimeStamp="<<fTimeStamp<<
1351 "EventType="<<fEventType<<
1363 // end of gg signal analysis
1368 if (rms06>fRecoParam->GetMaxNoise()) {
1369 pedestalEvent+=1024.;
1370 return 1024+median; // sign noisy channel in debug mode
1377 void AliTPCclustererMI::DumpHistos(){
1379 // Dump histogram information
1381 if (!fAmplitudeHisto) return;
1382 AliTPCROC * roc = AliTPCROC::Instance();
1383 for (UInt_t isector=0; isector<AliTPCROC::Instance()->GetNSectors(); isector++){
1384 TObjArray * array = (TObjArray*)fAmplitudeHisto->UncheckedAt(isector);
1385 if (!array) continue;
1386 for (UInt_t ipad = 0; ipad <(UInt_t)array->GetEntriesFast(); ipad++){
1387 TH1F * histo = (TH1F*) array->UncheckedAt(ipad);
1388 if (!histo) continue;
1389 if (histo->GetEntries()<100) continue;
1390 histo->Fit("gaus","q");
1391 Float_t mean = histo->GetMean();
1392 Float_t rms = histo->GetRMS();
1393 Float_t gmean = histo->GetFunction("gaus")->GetParameter(1);
1394 Float_t gsigma = histo->GetFunction("gaus")->GetParameter(2);
1395 Float_t gmeanErr = histo->GetFunction("gaus")->GetParError(1);
1396 Float_t gsigmaErr = histo->GetFunction("gaus")->GetParError(2);
1397 Float_t max = histo->GetFunction("gaus")->GetParameter(0);
1400 UInt_t row=0, pad =0;
1401 const UInt_t *indexes =roc->GetRowIndexes(isector);
1402 for (UInt_t irow=0; irow<roc->GetNRows(isector); irow++){
1403 if (indexes[irow]<=ipad){
1405 pad = ipad-indexes[irow];
1408 Int_t rpad = pad - (AliTPCROC::Instance()->GetNPads(isector,row))/2;
1410 (*fDebugStreamer)<<"Fit"<<
1411 "TimeStamp="<<fTimeStamp<<
1412 "EventType="<<fEventType<<
1413 "Sector="<<isector<<
1422 "GMeanErr="<<gmeanErr<<
1423 "GSigmaErr="<<gsigmaErr<<
1425 if (array->UncheckedAt(ipad)) fDebugStreamer->StoreObject(array->UncheckedAt(ipad));
1432 Int_t AliTPCclustererMI::TransformFFT(Float_t *input, Float_t threshold, Bool_t locMax, Float_t *freq, Float_t *re, Float_t *im, Float_t *mag, Float_t *phi)
1435 // calculate fourrie transform
1436 // return only frequncies with mag over threshold
1437 // if locMax is spectified only freque with local maxima over theshold is returned
1439 if (! fFFTr2c) return kFALSE;
1440 if (!freq) return kFALSE;
1443 Int_t nPoints = fRecoParam->GetLastBin()-fRecoParam->GetFirstBin();
1444 Double_t *in = new Double_t[nPoints];
1445 Double_t *rfft = new Double_t[nPoints];
1446 Double_t *ifft = new Double_t[nPoints];
1447 for (Int_t i=0; i<nPoints; i++){in[i]=input[i];}
1448 fFFTr2c->SetPoints(in);
1449 fFFTr2c->Transform();
1450 fFFTr2c->GetPointsComplex(rfft, ifft);
1451 for (Int_t i=3; i<nPoints/2-3; i++){
1452 Float_t lmag = TMath::Sqrt(rfft[i]*rfft[i]+ifft[i]*ifft[i])/nPoints;
1453 if (lmag<threshold) continue;
1455 if ( TMath::Sqrt(rfft[i-1]*rfft[i-1]+ifft[i-1]*ifft[i-1])/nPoints>lmag) continue;
1456 if ( TMath::Sqrt(rfft[i+1]*rfft[i+1]+ifft[i+1]*ifft[i+1])/nPoints>lmag) continue;
1457 if ( TMath::Sqrt(rfft[i-2]*rfft[i-2]+ifft[i-2]*ifft[i-2])/nPoints>lmag) continue;
1458 if ( TMath::Sqrt(rfft[i+2]*rfft[i+2]+ifft[i+2]*ifft[i+2])/nPoints>lmag) continue;
1459 if ( TMath::Sqrt(rfft[i-3]*rfft[i-3]+ifft[i-3]*ifft[i-3])/nPoints>lmag) continue;
1460 if ( TMath::Sqrt(rfft[i+3]*rfft[i+3]+ifft[i+3]*ifft[i+3])/nPoints>lmag) continue;
1463 freq[current] = Float_t(i)/Float_t(nPoints);
1465 re[current] = rfft[i];
1466 im[current] = ifft[i];
1468 phi[current]=TMath::ATan2(ifft[i],rfft[i]);