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 **************************************************************************/
23 //-------------------------------------------------------
24 // Implementation of the TPC clusterer
26 // Origin: Marian Ivanov
27 //-------------------------------------------------------
29 #include "AliTPCclustererMI.h"
30 #include "AliTPCclusterMI.h"
31 #include <TObjArray.h>
33 #include "AliTPCClustersArray.h"
34 #include "AliTPCClustersRow.h"
35 #include "AliDigits.h"
36 #include "AliSimDigits.h"
37 #include "AliTPCParam.h"
38 #include "Riostream.h"
41 ClassImp(AliTPCclustererMI)
45 AliTPCclustererMI::AliTPCclustererMI()
50 void AliTPCclustererMI::SetInput(TTree * tree)
53 // set input tree with digits
56 if (!fInput->GetBranch("Segment")){
57 cerr<<"AliTPC::Digits2Clusters(): no porper input tree !\n";
63 void AliTPCclustererMI::SetOutput(TTree * tree)
68 AliTPCClustersRow clrow;
69 AliTPCClustersRow *pclrow=&clrow;
70 clrow.SetClass("AliTPCclusterMI");
71 clrow.SetArray(1); // to make Clones array
72 fOutput->Branch("Segment","AliTPCClustersRow",&pclrow,32000,200);
76 Float_t AliTPCclustererMI::GetSigmaY2(Int_t iz){
77 // sigma y2 = in digits - we don't know the angle
78 Float_t z = iz*fParam->GetZWidth();
79 Float_t sd2 = (z*fParam->GetDiffL()*fParam->GetDiffL())/
80 (fPadWidth*fPadWidth);
82 Float_t res = sd2+sres;
87 Float_t AliTPCclustererMI::GetSigmaZ2(Int_t iz){
88 //sigma z2 = in digits - angle estimated supposing vertex constraint
89 Float_t z = iz*fZWidth;
90 Float_t sd2 = (z*fParam->GetDiffL()*fParam->GetDiffL())/(fZWidth*fZWidth);
91 Float_t angular = fPadLength*(fParam->GetZLength()-z)/(fRx*fZWidth);
94 Float_t sres = fParam->GetZSigma()/fZWidth;
96 Float_t res = angular +sd2+sres;
100 void AliTPCclustererMI::MakeCluster(Int_t k,Int_t max,Int_t *bins, UInt_t m,
103 Int_t i0=k/max; //central pad
104 Int_t j0=k%max; //central time bin
106 // set pointers to data
107 //Int_t dummy[5] ={0,0,0,0,0};
108 Int_t * matrix[5]; //5x5 matrix with digits - indexing i = 0 ..4 j = -2..2
109 Int_t * resmatrix[5];
110 for (Int_t di=-2;di<=2;di++){
111 matrix[di+2] = &bins[k+di*max];
112 resmatrix[di+2] = &fResBins[k+di*max];
114 //build matrix with virtual charge
115 Float_t sigmay2= GetSigmaY2(j0);
116 Float_t sigmaz2= GetSigmaZ2(j0);
118 Float_t vmatrix[5][5];
119 vmatrix[2][2] = matrix[2][0];
121 c.SetMax(Short_t(vmatrix[2][2])); // write maximal amplitude
122 for (Int_t di =-1;di <=1;di++)
123 for (Int_t dj =-1;dj <=1;dj++){
124 Float_t amp = matrix[di+2][dj];
125 if ( (amp<2) && (fLoop<2)){
126 // if under threshold - calculate virtual charge
127 Float_t ratio = TMath::Exp(-1.2*TMath::Abs(di)/sigmay2)*TMath::Exp(-1.2*TMath::Abs(dj)/sigmaz2);
128 amp = ((matrix[2][0]-2)*(matrix[2][0]-2)/(matrix[-di+2][-dj]+2))*ratio;
130 vmatrix[2+di][2+dj]=amp;
131 vmatrix[2+2*di][2+2*dj]=0;
134 vmatrix[2+2*di][2+dj] =0;
135 vmatrix[2+di][2+2*dj] =0;
140 //if small amplitude - below 2 x threshold - don't consider other one
141 vmatrix[2+di][2+dj]=amp;
142 vmatrix[2+2*di][2+2*dj]=0; // don't take to the account next bin
145 vmatrix[2+2*di][2+dj] =0;
146 vmatrix[2+di][2+2*dj] =0;
150 //if bigger then take everything
151 vmatrix[2+di][2+dj]=amp;
152 vmatrix[2+2*di][2+2*dj]= matrix[2*di+2][2*dj] ;
155 vmatrix[2+2*di][2+dj] = matrix[2*di+2][dj];
156 vmatrix[2+di][2+2*dj] = matrix[2+di][dj*2];
168 for (Int_t i=-2;i<=2;i++)
169 for (Int_t j=-2;j<=2;j++){
170 Float_t amp = vmatrix[i+2][j+2];
179 Float_t meani = sumiw/sumw;
180 Float_t mi2 = sumi2w/sumw-meani*meani;
181 Float_t meanj = sumjw/sumw;
182 Float_t mj2 = sumj2w/sumw-meanj*meanj;
184 Float_t ry = mi2/sigmay2;
185 Float_t rz = mj2/sigmaz2;
188 if ( ( (ry<0.6) || (rz<0.6) ) && fLoop==2) return;
189 if ( (ry <1.2) && (rz<1.2) ) {
190 //if cluster looks like expected
191 //+1.2 deviation from expected sigma accepted
192 // c.fMax = FitMax(vmatrix,meani,meanj,TMath::Sqrt(sigmay2),TMath::Sqrt(sigmaz2));
196 //set cluster parameters
198 c.SetY(meani*fPadWidth);
199 c.SetZ(meanj*fZWidth);
203 //remove cluster data from data
204 for (Int_t di=-2;di<=2;di++)
205 for (Int_t dj=-2;dj<=2;dj++){
206 resmatrix[di+2][dj] -= Int_t(vmatrix[di+2][dj+2]);
207 if (resmatrix[di+2][dj]<0) resmatrix[di+2][dj]=0;
213 //unfolding when neccessary
216 Int_t * matrix2[7]; //7x7 matrix with digits - indexing i = 0 ..6 j = -3..3
217 Int_t dummy[7]={0,0,0,0,0,0};
218 for (Int_t di=-3;di<=3;di++){
219 matrix2[di+3] = &bins[k+di*max];
220 if ((k+di*max)<3) matrix2[di+3] = &dummy[3];
221 if ((k+di*max)>fMaxBin-3) matrix2[di+3] = &dummy[3];
223 Float_t vmatrix2[5][5];
226 UnfoldCluster(matrix2,vmatrix2,meani,meanj,sumu,overlap);
228 // c.fMax = FitMax(vmatrix2,meani,meanj,TMath::Sqrt(sigmay2),TMath::Sqrt(sigmaz2));
231 //set cluster parameters
233 c.SetY(meani*fPadWidth);
234 c.SetZ(meanj*fZWidth);
237 c.SetType(Char_t(overlap)+1);
244 printf("%f\t%f\n", vmatrix2[2][2], vmatrix[2][2]);
249 void AliTPCclustererMI::UnfoldCluster(Int_t * matrix2[7], Float_t recmatrix[5][5], Float_t & meani, Float_t & meanj,
250 Float_t & sumu, Float_t & overlap )
253 //unfold cluster from input matrix
254 //data corresponding to cluster writen in recmatrix
255 //output meani and meanj
257 //take separatelly y and z
259 Float_t sum3i[7] = {0,0,0,0,0,0,0};
260 Float_t sum3j[7] = {0,0,0,0,0,0,0};
262 for (Int_t k =0;k<7;k++)
263 for (Int_t l = -1; l<=1;l++){
264 sum3i[k]+=matrix2[k][l];
265 sum3j[k]+=matrix2[l+3][k-3];
267 Float_t mratio[3][3]={{1,1,1},{1,1,1},{1,1,1}};
270 Float_t sum3wi = 0; //charge minus overlap
271 Float_t sum3wio = 0; //full charge
272 Float_t sum3iw = 0; //sum for mean value
273 for (Int_t dk=-1;dk<=1;dk++){
274 sum3wio+=sum3i[dk+3];
280 if ( ( ((sum3i[dk+3]+3)/(sum3i[3]-3))+1 < (sum3i[2*dk+3]-3)/(sum3i[dk+3]+3))||
281 sum3i[dk+3]<=sum3i[2*dk+3] && sum3i[dk+3]>2 ){
282 Float_t xm2 = sum3i[-dk+3];
283 Float_t xm1 = sum3i[+3];
284 Float_t x1 = sum3i[2*dk+3];
285 Float_t x2 = sum3i[3*dk+3];
286 Float_t w11 = TMath::Max((Float_t)(4.*xm1-xm2),(Float_t)0.000001);
287 Float_t w12 = TMath::Max((Float_t)(4 *x1 -x2),(Float_t)0.);
288 ratio = w11/(w11+w12);
289 for (Int_t dl=-1;dl<=1;dl++)
290 mratio[dk+1][dl+1] *= ratio;
292 Float_t amp = sum3i[dk+3]*ratio;
297 meani = sum3iw/sum3wi;
298 Float_t overlapi = (sum3wio-sum3wi)/sum3wio;
303 Float_t sum3wj = 0; //charge minus overlap
304 Float_t sum3wjo = 0; //full charge
305 Float_t sum3jw = 0; //sum for mean value
306 for (Int_t dk=-1;dk<=1;dk++){
307 sum3wjo+=sum3j[dk+3];
313 if ( ( ((sum3j[dk+3]+3)/(sum3j[3]-3))+1 < (sum3j[2*dk+3]-3)/(sum3j[dk+3]+3)) ||
314 (sum3j[dk+3]<=sum3j[2*dk+3] && sum3j[dk+3]>2)){
315 Float_t xm2 = sum3j[-dk+3];
316 Float_t xm1 = sum3j[+3];
317 Float_t x1 = sum3j[2*dk+3];
318 Float_t x2 = sum3j[3*dk+3];
319 Float_t w11 = TMath::Max((Float_t)(4.*xm1-xm2),(Float_t)0.000001);
320 Float_t w12 = TMath::Max((Float_t)(4 *x1 -x2),(Float_t)0.);
321 ratio = w11/(w11+w12);
322 for (Int_t dl=-1;dl<=1;dl++)
323 mratio[dl+1][dk+1] *= ratio;
325 Float_t amp = sum3j[dk+3]*ratio;
330 meanj = sum3jw/sum3wj;
331 Float_t overlapj = (sum3wjo-sum3wj)/sum3wjo;
332 overlap = Int_t(100*TMath::Max(overlapi,overlapj)+3);
333 sumu = (sum3wj+sum3wi)/2.;
336 //if not overlap detected remove everything
337 for (Int_t di =-2; di<=2;di++)
338 for (Int_t dj =-2; dj<=2;dj++){
339 recmatrix[di+2][dj+2] = matrix2[3+di][dj];
343 for (Int_t di =-1; di<=1;di++)
344 for (Int_t dj =-1; dj<=1;dj++){
346 if (mratio[di+1][dj+1]==1){
347 recmatrix[di+2][dj+2] = matrix2[3+di][dj];
348 if (TMath::Abs(di)+TMath::Abs(dj)>1){
349 recmatrix[2*di+2][dj+2] = matrix2[3+2*di][dj];
350 recmatrix[di+2][2*dj+2] = matrix2[3+di][2*dj];
352 recmatrix[2*di+2][2*dj+2] = matrix2[3+2*di][2*dj];
356 //if we have overlap in direction
357 recmatrix[di+2][dj+2] = mratio[di+1][dj+1]* matrix2[3+di][dj];
358 if (TMath::Abs(di)+TMath::Abs(dj)>1){
359 ratio = TMath::Min((Float_t)(recmatrix[di+2][dj+2]/(matrix2[3+0*di][1*dj]+1)),(Float_t)1.);
360 recmatrix[2*di+2][dj+2] = ratio*recmatrix[di+2][dj+2];
362 ratio = TMath::Min((Float_t)(recmatrix[di+2][dj+2]/(matrix2[3+1*di][0*dj]+1)),(Float_t)1.);
363 recmatrix[di+2][2*dj+2] = ratio*recmatrix[di+2][dj+2];
366 ratio = recmatrix[di+2][dj+2]/matrix2[3][0];
367 recmatrix[2*di+2][2*dj+2] = ratio*recmatrix[di+2][dj+2];
373 printf("%f\n", recmatrix[2][2]);
377 Float_t AliTPCclustererMI::FitMax(Float_t vmatrix[5][5], Float_t y, Float_t z, Float_t sigmay, Float_t sigmaz)
384 for (Int_t di = -1;di<=1;di++)
385 for (Int_t dj = -1;dj<=1;dj++){
386 if (vmatrix[2+di][2+dj]>2){
387 Float_t teor = TMath::Gaus(di,y,sigmay*1.2)*TMath::Gaus(dj,z,sigmaz*1.2);
388 sumteor += teor*vmatrix[2+di][2+dj];
389 sumamp += vmatrix[2+di][2+dj]*vmatrix[2+di][2+dj];
392 Float_t max = sumamp/sumteor;
396 void AliTPCclustererMI::AddCluster(AliTPCclusterMI &c){
398 // transform cluster to the global coordinata
399 // add the cluster to the array
401 Float_t meani = c.GetY()/fPadWidth;
402 Float_t meanj = c.GetZ()/fZWidth;
404 Int_t ki = TMath::Nint(meani-3);
406 if (ki>=fMaxPad) ki = fMaxPad-1;
407 Int_t kj = TMath::Nint(meanj-3);
409 if (kj>=fMaxTime-3) kj=fMaxTime-4;
410 // ki and kj shifted to "real" coordinata
411 c.SetLabel(fRowDig->GetTrackIDFast(kj,ki,0)-2,0);
412 c.SetLabel(fRowDig->GetTrackIDFast(kj,ki,1)-2,1);
413 c.SetLabel(fRowDig->GetTrackIDFast(kj,ki,2)-2,2);
416 Float_t s2 = c.GetSigmaY2();
417 Float_t w=fParam->GetPadPitchWidth(fSector);
419 c.SetSigmaY2(s2*w*w);
422 c.SetSigmaZ2(s2*w*w);
423 c.SetY((meani - 2.5 - 0.5*fMaxPad)*fParam->GetPadPitchWidth(fSector));
424 c.SetZ(fZWidth*(meanj-3));
425 c.SetZ(c.GetZ() - 3.*fParam->GetZSigma()); // PASA delay
426 c.SetZ(fSign*(fParam->GetZLength() - c.GetZ()));
428 if (ki<=1 || ki>=fMaxPad-1 || kj==1 || kj==fMaxTime-2) {
429 //c.SetSigmaY2(c.GetSigmaY2()*25.);
430 //c.SetSigmaZ2(c.GetSigmaZ2()*4.);
431 c.SetType(-(c.GetType()+3)); //edge clusters
433 if (fLoop==2) c.SetType(100);
435 TClonesArray * arr = fRowCl->GetArray();
436 // AliTPCclusterMI * cl =
437 new ((*arr)[fNcluster]) AliTPCclusterMI(c);
443 //_____________________________________________________________________________
444 void AliTPCclustererMI::Digits2Clusters(const AliTPCParam *par, Int_t eventn)
446 //-----------------------------------------------------------------
447 // This is a simple cluster finder.
448 //-----------------------------------------------------------------
449 TDirectory *savedir=gDirectory;
452 cerr<<"AliTPC::Digits2Clusters(): input tree not initialised !\n";
457 cerr<<"AliTPC::Digits2Clusters(): output tree not initialised !\n";
461 AliSimDigits digarr, *dummy=&digarr;
463 fInput->GetBranch("Segment")->SetAddress(&dummy);
464 Stat_t nentries = fInput->GetEntries();
466 fMaxTime=par->GetMaxTBin()+6; // add 3 virtual time bins before and 3 after
469 ((AliTPCParam*)par)->Write(par->GetTitle());
473 for (Int_t n=0; n<nentries; n++) {
476 if (!par->AdjustSectorRow(digarr.GetID(),fSector,row)) {
477 cerr<<"AliTPC warning: invalid segment ID ! "<<digarr.GetID()<<endl;
481 AliTPCClustersRow *clrow= new AliTPCClustersRow();
483 clrow->SetClass("AliTPCclusterMI");
486 clrow->SetID(digarr.GetID());
487 fOutput->GetBranch("Segment")->SetAddress(&clrow);
488 fRx=par->GetPadRowRadii(fSector,row);
491 const Int_t kNIS=par->GetNInnerSector(), kNOS=par->GetNOuterSector();
492 fZWidth = fParam->GetZWidth();
493 if (fSector < kNIS) {
494 fMaxPad = par->GetNPadsLow(row);
495 fSign = (fSector < kNIS/2) ? 1 : -1;
496 fPadLength = par->GetPadPitchLength(fSector,row);
497 fPadWidth = par->GetPadPitchWidth();
499 fMaxPad = par->GetNPadsUp(row);
500 fSign = ((fSector-kNIS) < kNOS/2) ? 1 : -1;
501 fPadLength = par->GetPadPitchLength(fSector,row);
502 fPadWidth = par->GetPadPitchWidth();
506 fMaxBin=fMaxTime*(fMaxPad+6); // add 3 virtual pads before and 3 after
507 fBins =new Int_t[fMaxBin];
508 fResBins =new Int_t[fMaxBin]; //fBins with residuals after 1 finder loop
509 memset(fBins,0,sizeof(Int_t)*fMaxBin);
511 if (digarr.First()) //MI change
513 Short_t dig=digarr.CurrentDigit();
514 if (dig<=par->GetZeroSup()) continue;
515 Int_t j=digarr.CurrentRow()+3, i=digarr.CurrentColumn()+3;
516 fBins[i*fMaxTime+j]=dig;
517 } while (digarr.Next());
518 digarr.ExpandTrackBuffer();
520 //add virtual charge at the edge
521 for (Int_t i=0; i<fMaxTime; i++){
522 Float_t amp1 = fBins[i+3*fMaxTime];
525 Float_t amp2 = fBins[i+4*fMaxTime];
526 if (amp2==0) amp2=0.5;
527 Float_t sigma2 = GetSigmaY2(i);
528 amp0 = (amp1*amp1/amp2)*TMath::Exp(-1./sigma2);
529 if (gDebug>4) printf("\n%f\n",amp0);
531 fBins[i+2*fMaxTime] = Int_t(amp0);
533 amp1 = fBins[(fMaxPad+2)*fMaxTime+i];
535 Float_t amp2 = fBins[i+(fMaxPad+1)*fMaxTime];
536 if (amp2==0) amp2=0.5;
537 Float_t sigma2 = GetSigmaY2(i);
538 amp0 = (amp1*amp1/amp2)*TMath::Exp(-1./sigma2);
539 if (gDebug>4) printf("\n%f\n",amp0);
541 fBins[(fMaxPad+3)*fMaxTime+i] = Int_t(amp0);
544 memcpy(fResBins,fBins, fMaxBin*2);
547 //first loop - for "gold cluster"
549 Int_t *b=&fBins[-1]+2*fMaxTime;
550 Int_t crtime = Int_t((fParam->GetZLength()-1.05*fRx)/fZWidth-5);
552 for (Int_t i=2*fMaxTime; i<fMaxBin-2*fMaxTime; i++) {
554 if (*b<8) continue; //threshold form maxima
555 if (i%fMaxTime<crtime) {
556 Int_t delta = -(i%fMaxTime)+crtime;
562 if (!IsMaximum(*b,fMaxTime,b)) continue;
565 MakeCluster(i, fMaxTime, fBins, dummy,c);
568 //memcpy(fBins,fResBins, fMaxBin*2);
569 //second loop - for rest cluster
572 b=&fResBins[-1]+2*fMaxTime;
573 for (Int_t i=2*fMaxTime; i<fMaxBin-2*fMaxTime; i++) {
575 if (*b<25) continue; // bigger threshold for maxima
576 if (!IsMaximum(*b,fMaxTime,b)) continue;
579 MakeCluster(i, fMaxTime, fResBins, dummy,c);
586 nclusters+=fNcluster;
590 cerr<<"Number of found clusters : "<<nclusters<<" \n";