]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/AliTPCclustererMI.cxx
Record changes.
[u/mrichter/AliRoot.git] / TPC / AliTPCclustererMI.cxx
CommitLineData
1c53abe2 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
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 **************************************************************************/
15
88cb7938 16/* $Id$ */
1c53abe2 17
18//-------------------------------------------------------
19// Implementation of the TPC clusterer
20//
21// Origin: Marian Ivanov
22//-------------------------------------------------------
23
d7a11555 24#include "AliTPCReconstructor.h"
1c53abe2 25#include "AliTPCclustererMI.h"
26#include "AliTPCclusterMI.h"
27#include <TObjArray.h>
28#include <TFile.h>
65c045f0 29#include "TGraph.h"
30#include "TF1.h"
31#include "TRandom.h"
32#include "AliMathBase.h"
33
1c53abe2 34#include "AliTPCClustersArray.h"
35#include "AliTPCClustersRow.h"
36#include "AliDigits.h"
37#include "AliSimDigits.h"
38#include "AliTPCParam.h"
f8aae377 39#include "AliRawReader.h"
40#include "AliTPCRawStream.h"
41#include "AliRunLoader.h"
42#include "AliLoader.h"
cc5e9db0 43#include "Riostream.h"
1c53abe2 44#include <TTree.h>
13116aec 45#include "AliTPCcalibDB.h"
46#include "AliTPCCalPad.h"
47#include "AliTPCCalROC.h"
65c045f0 48#include "TTreeStream.h"
49#include "AliLog.h"
13116aec 50
51
1c53abe2 52ClassImp(AliTPCclustererMI)
53
54
55
f8aae377 56AliTPCclustererMI::AliTPCclustererMI(const AliTPCParam* par)
1c53abe2 57{
22c352f8 58 fPedSubtraction = kFALSE;
59 fIsOldRCUFormat = kFALSE;
1c53abe2 60 fInput =0;
61 fOutput=0;
f8aae377 62 fParam = par;
65c045f0 63 fDebugStreamer = new TTreeSRedirector("TPCsignal.root");
64 fAmplitudeHisto = 0;
65}
66
67AliTPCclustererMI::~AliTPCclustererMI(){
68 DumpHistos();
69 if (fAmplitudeHisto) delete fAmplitudeHisto;
70 if (fDebugStreamer) delete fDebugStreamer;
1c53abe2 71}
22c352f8 72
1c53abe2 73void AliTPCclustererMI::SetInput(TTree * tree)
74{
75 //
76 // set input tree with digits
77 //
78 fInput = tree;
79 if (!fInput->GetBranch("Segment")){
80 cerr<<"AliTPC::Digits2Clusters(): no porper input tree !\n";
81 fInput=0;
82 return;
83 }
84}
85
86void AliTPCclustererMI::SetOutput(TTree * tree)
87{
88 //
89 //
90 fOutput= tree;
91 AliTPCClustersRow clrow;
92 AliTPCClustersRow *pclrow=&clrow;
93 clrow.SetClass("AliTPCclusterMI");
94 clrow.SetArray(1); // to make Clones array
95 fOutput->Branch("Segment","AliTPCClustersRow",&pclrow,32000,200);
96}
97
98
99Float_t AliTPCclustererMI::GetSigmaY2(Int_t iz){
100 // sigma y2 = in digits - we don't know the angle
753797ce 101 Float_t z = iz*fParam->GetZWidth()+fParam->GetNTBinsL1()*fParam->GetZWidth();
1c53abe2 102 Float_t sd2 = (z*fParam->GetDiffL()*fParam->GetDiffL())/
103 (fPadWidth*fPadWidth);
104 Float_t sres = 0.25;
105 Float_t res = sd2+sres;
106 return res;
107}
108
109
110Float_t AliTPCclustererMI::GetSigmaZ2(Int_t iz){
111 //sigma z2 = in digits - angle estimated supposing vertex constraint
753797ce 112 Float_t z = iz*fZWidth+fParam->GetNTBinsL1()*fParam->GetZWidth();
1c53abe2 113 Float_t sd2 = (z*fParam->GetDiffL()*fParam->GetDiffL())/(fZWidth*fZWidth);
114 Float_t angular = fPadLength*(fParam->GetZLength()-z)/(fRx*fZWidth);
115 angular*=angular;
116 angular/=12.;
117 Float_t sres = fParam->GetZSigma()/fZWidth;
118 sres *=sres;
119 Float_t res = angular +sd2+sres;
120 return res;
121}
122
13116aec 123void AliTPCclustererMI::MakeCluster(Int_t k,Int_t max,Float_t *bins, UInt_t /*m*/,
1c53abe2 124AliTPCclusterMI &c)
125{
126 Int_t i0=k/max; //central pad
127 Int_t j0=k%max; //central time bin
128
129 // set pointers to data
130 //Int_t dummy[5] ={0,0,0,0,0};
13116aec 131 Float_t * matrix[5]; //5x5 matrix with digits - indexing i = 0 ..4 j = -2..2
132 Float_t * resmatrix[5];
1c53abe2 133 for (Int_t di=-2;di<=2;di++){
134 matrix[di+2] = &bins[k+di*max];
135 resmatrix[di+2] = &fResBins[k+di*max];
136 }
137 //build matrix with virtual charge
138 Float_t sigmay2= GetSigmaY2(j0);
139 Float_t sigmaz2= GetSigmaZ2(j0);
140
141 Float_t vmatrix[5][5];
142 vmatrix[2][2] = matrix[2][0];
143 c.SetType(0);
6c024a0e 144 c.SetMax((UShort_t)(vmatrix[2][2])); // write maximal amplitude
1c53abe2 145 for (Int_t di =-1;di <=1;di++)
146 for (Int_t dj =-1;dj <=1;dj++){
147 Float_t amp = matrix[di+2][dj];
148 if ( (amp<2) && (fLoop<2)){
149 // if under threshold - calculate virtual charge
150 Float_t ratio = TMath::Exp(-1.2*TMath::Abs(di)/sigmay2)*TMath::Exp(-1.2*TMath::Abs(dj)/sigmaz2);
151 amp = ((matrix[2][0]-2)*(matrix[2][0]-2)/(matrix[-di+2][-dj]+2))*ratio;
152 if (amp>2) amp = 2;
153 vmatrix[2+di][2+dj]=amp;
154 vmatrix[2+2*di][2+2*dj]=0;
155 if ( (di*dj)!=0){
156 //DIAGONAL ELEMENTS
157 vmatrix[2+2*di][2+dj] =0;
158 vmatrix[2+di][2+2*dj] =0;
159 }
160 continue;
161 }
162 if (amp<4){
163 //if small amplitude - below 2 x threshold - don't consider other one
164 vmatrix[2+di][2+dj]=amp;
165 vmatrix[2+2*di][2+2*dj]=0; // don't take to the account next bin
166 if ( (di*dj)!=0){
167 //DIAGONAL ELEMENTS
168 vmatrix[2+2*di][2+dj] =0;
169 vmatrix[2+di][2+2*dj] =0;
170 }
171 continue;
172 }
173 //if bigger then take everything
174 vmatrix[2+di][2+dj]=amp;
175 vmatrix[2+2*di][2+2*dj]= matrix[2*di+2][2*dj] ;
176 if ( (di*dj)!=0){
177 //DIAGONAL ELEMENTS
178 vmatrix[2+2*di][2+dj] = matrix[2*di+2][dj];
179 vmatrix[2+di][2+2*dj] = matrix[2+di][dj*2];
180 }
181 }
182
183
184
185 Float_t sumw=0;
186 Float_t sumiw=0;
187 Float_t sumi2w=0;
188 Float_t sumjw=0;
189 Float_t sumj2w=0;
190 //
191 for (Int_t i=-2;i<=2;i++)
192 for (Int_t j=-2;j<=2;j++){
193 Float_t amp = vmatrix[i+2][j+2];
194
195 sumw += amp;
196 sumiw += i*amp;
197 sumi2w += i*i*amp;
198 sumjw += j*amp;
199 sumj2w += j*j*amp;
200 }
201 //
202 Float_t meani = sumiw/sumw;
203 Float_t mi2 = sumi2w/sumw-meani*meani;
204 Float_t meanj = sumjw/sumw;
205 Float_t mj2 = sumj2w/sumw-meanj*meanj;
206 //
207 Float_t ry = mi2/sigmay2;
208 Float_t rz = mj2/sigmaz2;
209
210 //
211 if ( ( (ry<0.6) || (rz<0.6) ) && fLoop==2) return;
212 if ( (ry <1.2) && (rz<1.2) ) {
213 //if cluster looks like expected
214 //+1.2 deviation from expected sigma accepted
215 // c.fMax = FitMax(vmatrix,meani,meanj,TMath::Sqrt(sigmay2),TMath::Sqrt(sigmaz2));
216
217 meani +=i0;
218 meanj +=j0;
219 //set cluster parameters
220 c.SetQ(sumw);
221 c.SetY(meani*fPadWidth);
222 c.SetZ(meanj*fZWidth);
223 c.SetSigmaY2(mi2);
224 c.SetSigmaZ2(mj2);
225 AddCluster(c);
226 //remove cluster data from data
227 for (Int_t di=-2;di<=2;di++)
228 for (Int_t dj=-2;dj<=2;dj++){
13116aec 229 resmatrix[di+2][dj] -= vmatrix[di+2][dj+2];
1c53abe2 230 if (resmatrix[di+2][dj]<0) resmatrix[di+2][dj]=0;
231 }
232 resmatrix[2][0] =0;
233 return;
234 }
235 //
236 //unfolding when neccessary
237 //
238
13116aec 239 Float_t * matrix2[7]; //7x7 matrix with digits - indexing i = 0 ..6 j = -3..3
240 Float_t dummy[7]={0,0,0,0,0,0};
1c53abe2 241 for (Int_t di=-3;di<=3;di++){
242 matrix2[di+3] = &bins[k+di*max];
243 if ((k+di*max)<3) matrix2[di+3] = &dummy[3];
244 if ((k+di*max)>fMaxBin-3) matrix2[di+3] = &dummy[3];
245 }
246 Float_t vmatrix2[5][5];
247 Float_t sumu;
248 Float_t overlap;
249 UnfoldCluster(matrix2,vmatrix2,meani,meanj,sumu,overlap);
250 //
251 // c.fMax = FitMax(vmatrix2,meani,meanj,TMath::Sqrt(sigmay2),TMath::Sqrt(sigmaz2));
252 meani +=i0;
253 meanj +=j0;
254 //set cluster parameters
255 c.SetQ(sumu);
256 c.SetY(meani*fPadWidth);
257 c.SetZ(meanj*fZWidth);
258 c.SetSigmaY2(mi2);
259 c.SetSigmaZ2(mj2);
260 c.SetType(Char_t(overlap)+1);
261 AddCluster(c);
262
263 //unfolding 2
264 meani-=i0;
265 meanj-=j0;
266 if (gDebug>4)
267 printf("%f\t%f\n", vmatrix2[2][2], vmatrix[2][2]);
268}
269
270
271
13116aec 272void AliTPCclustererMI::UnfoldCluster(Float_t * matrix2[7], Float_t recmatrix[5][5], Float_t & meani, Float_t & meanj,
1c53abe2 273 Float_t & sumu, Float_t & overlap )
274{
275 //
276 //unfold cluster from input matrix
277 //data corresponding to cluster writen in recmatrix
278 //output meani and meanj
279
280 //take separatelly y and z
281
282 Float_t sum3i[7] = {0,0,0,0,0,0,0};
283 Float_t sum3j[7] = {0,0,0,0,0,0,0};
284
285 for (Int_t k =0;k<7;k++)
286 for (Int_t l = -1; l<=1;l++){
287 sum3i[k]+=matrix2[k][l];
288 sum3j[k]+=matrix2[l+3][k-3];
289 }
290 Float_t mratio[3][3]={{1,1,1},{1,1,1},{1,1,1}};
291 //
292 //unfold y
293 Float_t sum3wi = 0; //charge minus overlap
294 Float_t sum3wio = 0; //full charge
295 Float_t sum3iw = 0; //sum for mean value
296 for (Int_t dk=-1;dk<=1;dk++){
297 sum3wio+=sum3i[dk+3];
298 if (dk==0){
299 sum3wi+=sum3i[dk+3];
300 }
301 else{
302 Float_t ratio =1;
303 if ( ( ((sum3i[dk+3]+3)/(sum3i[3]-3))+1 < (sum3i[2*dk+3]-3)/(sum3i[dk+3]+3))||
304 sum3i[dk+3]<=sum3i[2*dk+3] && sum3i[dk+3]>2 ){
305 Float_t xm2 = sum3i[-dk+3];
306 Float_t xm1 = sum3i[+3];
307 Float_t x1 = sum3i[2*dk+3];
308 Float_t x2 = sum3i[3*dk+3];
cc5e9db0 309 Float_t w11 = TMath::Max((Float_t)(4.*xm1-xm2),(Float_t)0.000001);
310 Float_t w12 = TMath::Max((Float_t)(4 *x1 -x2),(Float_t)0.);
1c53abe2 311 ratio = w11/(w11+w12);
312 for (Int_t dl=-1;dl<=1;dl++)
313 mratio[dk+1][dl+1] *= ratio;
314 }
315 Float_t amp = sum3i[dk+3]*ratio;
316 sum3wi+=amp;
317 sum3iw+= dk*amp;
318 }
319 }
320 meani = sum3iw/sum3wi;
321 Float_t overlapi = (sum3wio-sum3wi)/sum3wio;
322
323
324
325 //unfold z
326 Float_t sum3wj = 0; //charge minus overlap
327 Float_t sum3wjo = 0; //full charge
328 Float_t sum3jw = 0; //sum for mean value
329 for (Int_t dk=-1;dk<=1;dk++){
330 sum3wjo+=sum3j[dk+3];
331 if (dk==0){
332 sum3wj+=sum3j[dk+3];
333 }
334 else{
335 Float_t ratio =1;
336 if ( ( ((sum3j[dk+3]+3)/(sum3j[3]-3))+1 < (sum3j[2*dk+3]-3)/(sum3j[dk+3]+3)) ||
337 (sum3j[dk+3]<=sum3j[2*dk+3] && sum3j[dk+3]>2)){
338 Float_t xm2 = sum3j[-dk+3];
339 Float_t xm1 = sum3j[+3];
340 Float_t x1 = sum3j[2*dk+3];
341 Float_t x2 = sum3j[3*dk+3];
cc5e9db0 342 Float_t w11 = TMath::Max((Float_t)(4.*xm1-xm2),(Float_t)0.000001);
343 Float_t w12 = TMath::Max((Float_t)(4 *x1 -x2),(Float_t)0.);
1c53abe2 344 ratio = w11/(w11+w12);
345 for (Int_t dl=-1;dl<=1;dl++)
346 mratio[dl+1][dk+1] *= ratio;
347 }
348 Float_t amp = sum3j[dk+3]*ratio;
349 sum3wj+=amp;
350 sum3jw+= dk*amp;
351 }
352 }
353 meanj = sum3jw/sum3wj;
354 Float_t overlapj = (sum3wjo-sum3wj)/sum3wjo;
355 overlap = Int_t(100*TMath::Max(overlapi,overlapj)+3);
356 sumu = (sum3wj+sum3wi)/2.;
357
358 if (overlap ==3) {
359 //if not overlap detected remove everything
360 for (Int_t di =-2; di<=2;di++)
361 for (Int_t dj =-2; dj<=2;dj++){
362 recmatrix[di+2][dj+2] = matrix2[3+di][dj];
363 }
364 }
365 else{
366 for (Int_t di =-1; di<=1;di++)
367 for (Int_t dj =-1; dj<=1;dj++){
368 Float_t ratio =1;
369 if (mratio[di+1][dj+1]==1){
370 recmatrix[di+2][dj+2] = matrix2[3+di][dj];
371 if (TMath::Abs(di)+TMath::Abs(dj)>1){
372 recmatrix[2*di+2][dj+2] = matrix2[3+2*di][dj];
373 recmatrix[di+2][2*dj+2] = matrix2[3+di][2*dj];
374 }
375 recmatrix[2*di+2][2*dj+2] = matrix2[3+2*di][2*dj];
376 }
377 else
378 {
379 //if we have overlap in direction
380 recmatrix[di+2][dj+2] = mratio[di+1][dj+1]* matrix2[3+di][dj];
381 if (TMath::Abs(di)+TMath::Abs(dj)>1){
cc5e9db0 382 ratio = TMath::Min((Float_t)(recmatrix[di+2][dj+2]/(matrix2[3+0*di][1*dj]+1)),(Float_t)1.);
1c53abe2 383 recmatrix[2*di+2][dj+2] = ratio*recmatrix[di+2][dj+2];
384 //
cc5e9db0 385 ratio = TMath::Min((Float_t)(recmatrix[di+2][dj+2]/(matrix2[3+1*di][0*dj]+1)),(Float_t)1.);
1c53abe2 386 recmatrix[di+2][2*dj+2] = ratio*recmatrix[di+2][dj+2];
387 }
388 else{
389 ratio = recmatrix[di+2][dj+2]/matrix2[3][0];
390 recmatrix[2*di+2][2*dj+2] = ratio*recmatrix[di+2][dj+2];
391 }
392 }
393 }
394 }
395 if (gDebug>4)
396 printf("%f\n", recmatrix[2][2]);
397
398}
399
400Float_t AliTPCclustererMI::FitMax(Float_t vmatrix[5][5], Float_t y, Float_t z, Float_t sigmay, Float_t sigmaz)
401{
402 //
403 // estimate max
404 Float_t sumteor= 0;
405 Float_t sumamp = 0;
406
407 for (Int_t di = -1;di<=1;di++)
408 for (Int_t dj = -1;dj<=1;dj++){
409 if (vmatrix[2+di][2+dj]>2){
410 Float_t teor = TMath::Gaus(di,y,sigmay*1.2)*TMath::Gaus(dj,z,sigmaz*1.2);
411 sumteor += teor*vmatrix[2+di][2+dj];
412 sumamp += vmatrix[2+di][2+dj]*vmatrix[2+di][2+dj];
413 }
414 }
415 Float_t max = sumamp/sumteor;
416 return max;
417}
418
419void AliTPCclustererMI::AddCluster(AliTPCclusterMI &c){
420 //
421 // transform cluster to the global coordinata
422 // add the cluster to the array
423 //
424 Float_t meani = c.GetY()/fPadWidth;
425 Float_t meanj = c.GetZ()/fZWidth;
426
427 Int_t ki = TMath::Nint(meani-3);
428 if (ki<0) ki=0;
429 if (ki>=fMaxPad) ki = fMaxPad-1;
430 Int_t kj = TMath::Nint(meanj-3);
431 if (kj<0) kj=0;
432 if (kj>=fMaxTime-3) kj=fMaxTime-4;
433 // ki and kj shifted to "real" coordinata
f8aae377 434 if (fRowDig) {
435 c.SetLabel(fRowDig->GetTrackIDFast(kj,ki,0)-2,0);
436 c.SetLabel(fRowDig->GetTrackIDFast(kj,ki,1)-2,1);
437 c.SetLabel(fRowDig->GetTrackIDFast(kj,ki,2)-2,2);
438 }
1c53abe2 439
440
441 Float_t s2 = c.GetSigmaY2();
442 Float_t w=fParam->GetPadPitchWidth(fSector);
443
444 c.SetSigmaY2(s2*w*w);
445 s2 = c.GetSigmaZ2();
446 w=fZWidth;
447 c.SetSigmaZ2(s2*w*w);
448 c.SetY((meani - 2.5 - 0.5*fMaxPad)*fParam->GetPadPitchWidth(fSector));
449 c.SetZ(fZWidth*(meanj-3));
753797ce 450 c.SetZ(c.GetZ() - 3.*fParam->GetZSigma() + fParam->GetNTBinsL1()*fParam->GetZWidth()); // PASA delay + L1 delay
1c53abe2 451 c.SetZ(fSign*(fParam->GetZLength() - c.GetZ()));
508541c7 452 c.SetX(fRx);
453 c.SetDetector(fSector);
454 c.SetRow(fRow);
455
1c53abe2 456 if (ki<=1 || ki>=fMaxPad-1 || kj==1 || kj==fMaxTime-2) {
457 //c.SetSigmaY2(c.GetSigmaY2()*25.);
458 //c.SetSigmaZ2(c.GetSigmaZ2()*4.);
459 c.SetType(-(c.GetType()+3)); //edge clusters
460 }
461 if (fLoop==2) c.SetType(100);
462
463 TClonesArray * arr = fRowCl->GetArray();
45bcf167 464 // AliTPCclusterMI * cl =
465 new ((*arr)[fNcluster]) AliTPCclusterMI(c);
1c53abe2 466
467 fNcluster++;
468}
469
470
471//_____________________________________________________________________________
f8aae377 472void AliTPCclustererMI::Digits2Clusters()
1c53abe2 473{
474 //-----------------------------------------------------------------
475 // This is a simple cluster finder.
476 //-----------------------------------------------------------------
1c53abe2 477
f8aae377 478 if (!fInput) {
479 Error("Digits2Clusters", "input tree not initialised");
1c53abe2 480 return;
481 }
482
f8aae377 483 if (!fOutput) {
484 Error("Digits2Clusters", "output tree not initialised");
485 return;
1c53abe2 486 }
487
13116aec 488 AliTPCCalPad * gainTPC = AliTPCcalibDB::Instance()->GetPadGainFactor();
489
1c53abe2 490 AliSimDigits digarr, *dummy=&digarr;
491 fRowDig = dummy;
492 fInput->GetBranch("Segment")->SetAddress(&dummy);
493 Stat_t nentries = fInput->GetEntries();
494
f8aae377 495 fMaxTime=fParam->GetMaxTBin()+6; // add 3 virtual time bins before and 3 after
1c53abe2 496
1c53abe2 497 Int_t nclusters = 0;
13116aec 498
1c53abe2 499 for (Int_t n=0; n<nentries; n++) {
500 fInput->GetEvent(n);
508541c7 501 if (!fParam->AdjustSectorRow(digarr.GetID(),fSector,fRow)) {
1c53abe2 502 cerr<<"AliTPC warning: invalid segment ID ! "<<digarr.GetID()<<endl;
503 continue;
504 }
508541c7 505 Int_t row = fRow;
13116aec 506 AliTPCCalROC * gainROC = gainTPC->GetCalROC(fSector); // pad gains per given sector
507
508 //
1c53abe2 509 AliTPCClustersRow *clrow= new AliTPCClustersRow();
510 fRowCl = clrow;
511 clrow->SetClass("AliTPCclusterMI");
512 clrow->SetArray(1);
513
514 clrow->SetID(digarr.GetID());
515 fOutput->GetBranch("Segment")->SetAddress(&clrow);
f8aae377 516 fRx=fParam->GetPadRowRadii(fSector,row);
1c53abe2 517
518
f8aae377 519 const Int_t kNIS=fParam->GetNInnerSector(), kNOS=fParam->GetNOuterSector();
1c53abe2 520 fZWidth = fParam->GetZWidth();
521 if (fSector < kNIS) {
f8aae377 522 fMaxPad = fParam->GetNPadsLow(row);
1c53abe2 523 fSign = (fSector < kNIS/2) ? 1 : -1;
f8aae377 524 fPadLength = fParam->GetPadPitchLength(fSector,row);
525 fPadWidth = fParam->GetPadPitchWidth();
1c53abe2 526 } else {
f8aae377 527 fMaxPad = fParam->GetNPadsUp(row);
1c53abe2 528 fSign = ((fSector-kNIS) < kNOS/2) ? 1 : -1;
f8aae377 529 fPadLength = fParam->GetPadPitchLength(fSector,row);
530 fPadWidth = fParam->GetPadPitchWidth();
1c53abe2 531 }
532
533
534 fMaxBin=fMaxTime*(fMaxPad+6); // add 3 virtual pads before and 3 after
13116aec 535 fBins =new Float_t[fMaxBin];
536 fResBins =new Float_t[fMaxBin]; //fBins with residuals after 1 finder loop
537 memset(fBins,0,sizeof(Float_t)*fMaxBin);
538 memset(fResBins,0,sizeof(Float_t)*fMaxBin);
1c53abe2 539
540 if (digarr.First()) //MI change
541 do {
13116aec 542 Float_t dig=digarr.CurrentDigit();
f8aae377 543 if (dig<=fParam->GetZeroSup()) continue;
1c53abe2 544 Int_t j=digarr.CurrentRow()+3, i=digarr.CurrentColumn()+3;
9d4f75a9 545 Float_t gain = gainROC->GetValue(row,digarr.CurrentColumn());
13116aec 546 fBins[i*fMaxTime+j]=dig/gain;
1c53abe2 547 } while (digarr.Next());
548 digarr.ExpandTrackBuffer();
549
f8aae377 550 FindClusters();
8569a2b0 551
552 fOutput->Fill();
88cb7938 553 delete clrow;
554 nclusters+=fNcluster;
8569a2b0 555 delete[] fBins;
556 delete[] fResBins;
88cb7938 557 }
f8aae377 558
19dd5b2f 559 Info("Digits2Clusters", "Number of found clusters : %d", nclusters);
f8aae377 560}
561
562void AliTPCclustererMI::Digits2Clusters(AliRawReader* rawReader)
563{
564//-----------------------------------------------------------------
22c352f8 565// This is a cluster finder for the TPC raw data.
566// The method assumes NO ordering of the altro channels.
567// The pedestal subtraction can be switched on and off
568// using an option of the TPC reconstructor
f8aae377 569//-----------------------------------------------------------------
570
571 if (!fOutput) {
572 Error("Digits2Clusters", "output tree not initialised");
573 return;
574 }
575
f8aae377 576 fRowDig = NULL;
577
22c352f8 578 AliTPCCalPad * gainTPC = AliTPCcalibDB::Instance()->GetPadGainFactor();
579
f8aae377 580 Int_t nclusters = 0;
88cb7938 581
f8aae377 582 fMaxTime = fParam->GetMaxTBin() + 6; // add 3 virtual time bins before and 3 after
583 const Int_t kNIS = fParam->GetNInnerSector();
584 const Int_t kNOS = fParam->GetNOuterSector();
585 const Int_t kNS = kNIS + kNOS;
586 fZWidth = fParam->GetZWidth();
587 Int_t zeroSup = fParam->GetZeroSup();
588
22c352f8 589 Float_t** allBins = NULL;
590 Float_t** allBinsRes = NULL;
f8aae377 591
22c352f8 592 // Loop over sectors
593 for(fSector = 0; fSector < kNS; fSector++) {
f8aae377 594
22c352f8 595 AliTPCCalROC * gainROC = gainTPC->GetCalROC(fSector); // pad gains per given sector
596
597 Int_t nRows = 0;
598 Int_t nDDLs = 0, indexDDL = 0;
599 if (fSector < kNIS) {
600 nRows = fParam->GetNRowLow();
601 fSign = (fSector < kNIS/2) ? 1 : -1;
602 nDDLs = 2;
603 indexDDL = fSector * 2;
604 }
605 else {
606 nRows = fParam->GetNRowUp();
607 fSign = ((fSector-kNIS) < kNOS/2) ? 1 : -1;
608 nDDLs = 4;
609 indexDDL = (fSector-kNIS) * 4 + kNIS * 2;
610 }
611
612 allBins = new Float_t*[nRows];
613 allBinsRes = new Float_t*[nRows];
614
615 for (Int_t iRow = 0; iRow < nRows; iRow++) {
616 Int_t maxPad;
617 if (fSector < kNIS)
618 maxPad = fParam->GetNPadsLow(iRow);
619 else
620 maxPad = fParam->GetNPadsUp(iRow);
621
622 Int_t maxBin = fMaxTime*(maxPad+6); // add 3 virtual pads before and 3 after
623 allBins[iRow] = new Float_t[maxBin];
624 allBinsRes[iRow] = new Float_t[maxBin];
625 memset(allBins[iRow],0,sizeof(Float_t)*maxBin);
626 memset(allBinsRes[iRow],0,sizeof(Float_t)*maxBin);
627 }
628
629 // Loas the raw data for corresponding DDLs
630 rawReader->Reset();
631 AliTPCRawStream input(rawReader);
632 input.SetOldRCUFormat(fIsOldRCUFormat);
362c9d61 633 rawReader->Select("TPC",indexDDL,indexDDL+nDDLs-1);
f8aae377 634
22c352f8 635 // Begin loop over altro data
636 while (input.Next()) {
637
638 if (input.GetSector() != fSector)
639 AliFatal(Form("Sector index mismatch ! Expected (%d), but got (%d) !",fSector,input.GetSector()));
640
641 if (input.GetTime() < 40) continue;
642
643 Int_t iRow = input.GetRow();
644 if (iRow < 0 || iRow >= nRows)
645 AliFatal(Form("Pad-row index (%d) outside the range (%d -> %d) !",
646 iRow, 0, nRows -1));
647
648 Int_t iPad = input.GetPad() + 3;
649
650 Int_t maxPad;
651 if (fSector < kNIS)
652 maxPad = fParam->GetNPadsLow(iRow);
653 else
654 maxPad = fParam->GetNPadsUp(iRow);
655
656 if (input.GetPad() < 0 || input.GetPad() >= maxPad)
657 AliFatal(Form("Pad index (%d) outside the range (%d -> %d) !",
658 input.GetPad(), 0, maxPad -1));
659
660 Int_t iTimeBin = input.GetTime() + 3;
661 if ( input.GetTime() < 0 || input.GetTime() >= fParam->GetMaxTBin())
662 AliFatal(Form("Timebin index (%d) outside the range (%d -> %d) !",
663 input.GetTime(), 0, fParam->GetMaxTBin() -1));
664
665 Int_t maxBin = fMaxTime*(maxPad+6); // add 3 virtual pads before and 3 after
666
667 if (((iPad*fMaxTime+iTimeBin) >= maxBin) ||
668 ((iPad*fMaxTime+iTimeBin) < 0))
669 AliFatal(Form("Index outside the allowed range"
670 " Sector=%d Row=%d Pad=%d Timebin=%d"
671 " (Max.index=%d)",fSector,iRow,iPad,iTimeBin,maxBin));
672
673 Float_t signal = input.GetSignal();
674 if (!fPedSubtraction && signal <= zeroSup) continue;
675
676 Float_t gain = gainROC->GetValue(iRow,input.GetPad());
677 allBins[iRow][iPad*fMaxTime+iTimeBin] = signal/gain;
65c045f0 678 allBins[iRow][iPad*fMaxTime+0] = 1; // pad with signal
22c352f8 679 } // End of the loop over altro data
680
681 // Now loop over rows and perform pedestal subtraction
682 if (fPedSubtraction) {
683 for (Int_t iRow = 0; iRow < nRows; iRow++) {
684 Int_t maxPad;
685 if (fSector < kNIS)
686 maxPad = fParam->GetNPadsLow(iRow);
687 else
688 maxPad = fParam->GetNPadsUp(iRow);
689
690 for (Int_t iPad = 0; iPad < maxPad + 6; iPad++) {
65c045f0 691 if (allBins[iRow][iPad*fMaxTime+0] !=1) continue; // no data
22c352f8 692 Float_t *p = &allBins[iRow][iPad*fMaxTime+3];
65c045f0 693 //Float_t pedestal = TMath::Median(fMaxTime, p);
694 Int_t id[3] = {fSector, iRow, iPad-3};
695 Float_t pedestal = ProcesSignal(p, fMaxTime, id);
22c352f8 696 for (Int_t iTimeBin = 0; iTimeBin < fMaxTime; iTimeBin++) {
697 allBins[iRow][iPad*fMaxTime+iTimeBin] -= pedestal;
698 if (allBins[iRow][iPad*fMaxTime+iTimeBin] < zeroSup)
699 allBins[iRow][iPad*fMaxTime+iTimeBin] = 0;
700 }
701 }
f8aae377 702 }
22c352f8 703 }
704
705 // Now loop over rows and find clusters
706 for (fRow = 0; fRow < nRows; fRow++) {
707 fRowCl = new AliTPCClustersRow;
708 fRowCl->SetClass("AliTPCclusterMI");
709 fRowCl->SetArray(1);
710 fRowCl->SetID(fParam->GetIndex(fSector, fRow));
711 fOutput->GetBranch("Segment")->SetAddress(&fRowCl);
712
713 fRx = fParam->GetPadRowRadii(fSector, fRow);
714 fPadLength = fParam->GetPadPitchLength(fSector, fRow);
f8aae377 715 fPadWidth = fParam->GetPadPitchWidth();
22c352f8 716 if (fSector < kNIS)
717 fMaxPad = fParam->GetNPadsLow(fRow);
718 else
719 fMaxPad = fParam->GetNPadsUp(fRow);
f8aae377 720 fMaxBin = fMaxTime*(fMaxPad+6); // add 3 virtual pads before and 3 after
22c352f8 721
722 fBins = allBins[fRow];
723 fResBins = allBinsRes[fRow];
724
725 FindClusters();
726
727 fOutput->Fill();
728 delete fRowCl;
729 nclusters += fNcluster;
730 } // End of loop to find clusters
731
732
733 for (Int_t iRow = 0; iRow < nRows; iRow++) {
734 delete [] allBins[iRow];
735 delete [] allBinsRes[iRow];
f8aae377 736 }
737
22c352f8 738 delete [] allBins;
739 delete [] allBinsRes;
f8aae377 740
22c352f8 741 } // End of loop over sectors
f8aae377 742
f8aae377 743 Info("Digits2Clusters", "Number of found clusters : %d\n", nclusters);
22c352f8 744
f8aae377 745}
746
747void AliTPCclustererMI::FindClusters()
748{
749 //add virtual charge at the edge
750 for (Int_t i=0; i<fMaxTime; i++){
751 Float_t amp1 = fBins[i+3*fMaxTime];
752 Float_t amp0 =0;
753 if (amp1>0){
754 Float_t amp2 = fBins[i+4*fMaxTime];
755 if (amp2==0) amp2=0.5;
756 Float_t sigma2 = GetSigmaY2(i);
757 amp0 = (amp1*amp1/amp2)*TMath::Exp(-1./sigma2);
758 if (gDebug>4) printf("\n%f\n",amp0);
759 }
13116aec 760 fBins[i+2*fMaxTime] = amp0;
f8aae377 761 amp0 = 0;
762 amp1 = fBins[(fMaxPad+2)*fMaxTime+i];
763 if (amp1>0){
764 Float_t amp2 = fBins[i+(fMaxPad+1)*fMaxTime];
765 if (amp2==0) amp2=0.5;
766 Float_t sigma2 = GetSigmaY2(i);
767 amp0 = (amp1*amp1/amp2)*TMath::Exp(-1./sigma2);
768 if (gDebug>4) printf("\n%f\n",amp0);
769 }
13116aec 770 fBins[(fMaxPad+3)*fMaxTime+i] = amp0;
f8aae377 771 }
772
773// memcpy(fResBins,fBins, fMaxBin*2);
774 memcpy(fResBins,fBins, fMaxBin);
775 //
776 fNcluster=0;
777 //first loop - for "gold cluster"
778 fLoop=1;
13116aec 779 Float_t *b=&fBins[-1]+2*fMaxTime;
753797ce 780 Int_t crtime = Int_t((fParam->GetZLength()-AliTPCReconstructor::GetCtgRange()*fRx)/fZWidth-fParam->GetNTBinsL1()-5);
f8aae377 781
782 for (Int_t i=2*fMaxTime; i<fMaxBin-2*fMaxTime; i++) {
783 b++;
784 if (*b<8) continue; //threshold form maxima
785 if (i%fMaxTime<crtime) {
786 Int_t delta = -(i%fMaxTime)+crtime;
787 b+=delta;
788 i+=delta;
789 continue;
790 }
791
792 if (!IsMaximum(*b,fMaxTime,b)) continue;
793 AliTPCclusterMI c;
794 Int_t dummy=0;
795 MakeCluster(i, fMaxTime, fBins, dummy,c);
796 //}
797 }
798 //memcpy(fBins,fResBins, fMaxBin*2);
799 //second loop - for rest cluster
800 /*
801 fLoop=2;
802 b=&fResBins[-1]+2*fMaxTime;
803 for (Int_t i=2*fMaxTime; i<fMaxBin-2*fMaxTime; i++) {
804 b++;
805 if (*b<25) continue; // bigger threshold for maxima
806 if (!IsMaximum(*b,fMaxTime,b)) continue;
807 AliTPCclusterMI c;
808 Int_t dummy;
809 MakeCluster(i, fMaxTime, fResBins, dummy,c);
810 //}
811 }
812 */
8569a2b0 813}
65c045f0 814
815
816Double_t AliTPCclustererMI::ProcesSignal(Float_t *signal, Int_t nchannels, Int_t id[3]){
817 //
818 // process signal on given pad - + streaming of additional information in special mode
819 //
820 // id[0] - sector
821 // id[1] - row
822 // id[2] - pad
823 Int_t offset =100;
824 Float_t kMin =50;
825 Float_t kMaxNoise = 3;
826 Double_t median = TMath::Median(nchannels-offset, &(signal[offset]));
827 if (AliLog::GetDebugLevel("","AliTPCclustererMI")==0) return median;
828 //
829
830 Double_t mean = TMath::Mean(nchannels-offset, &(signal[offset]));
831 Double_t rms = TMath::RMS(nchannels-offset, &(signal[offset]));
832 Double_t *dsignal = new Double_t[nchannels];
833 Double_t *dtime = new Double_t[nchannels];
834 Float_t max = 0;
835 Float_t maxPos = 0;
836 for (Int_t i=0; i<nchannels; i++){
837 dtime[i] = i;
838 dsignal[i] = signal[i];
839 if (signal[i]>max && i <fMaxTime-100) { // temporary remove spike signals at the end
840 max = signal[i];
841 maxPos = i;
842 }
843 }
844
845 Double_t mean06=0,mean07=0, mean08=0, mean09=0;
846 Double_t rms06=0, rms07=0, rms08, rms09=0;
847 AliMathBase::EvaluateUni(nchannels-offset, &(dsignal[offset]), mean06, rms06, int(0.6*(nchannels-offset)));
848 AliMathBase::EvaluateUni(nchannels-offset, &(dsignal[offset]), mean07, rms07, int(0.7*(nchannels-offset)));
849 AliMathBase::EvaluateUni(nchannels-offset, &(dsignal[offset]), mean08, rms08, int(0.8*(nchannels-offset)));
850 AliMathBase::EvaluateUni(nchannels-offset, &(dsignal[offset]), mean09, rms09, int(0.9*(nchannels-offset)));
851
852 //
853 UInt_t uid[3] = {UInt_t(id[0]),UInt_t(id[1]),UInt_t(id[2])};
854 if (uid[0]< AliTPCROC::Instance()->GetNSectors()
855 && uid[1]< AliTPCROC::Instance()->GetNRows(uid[0]) &&
856 uid[2] < AliTPCROC::Instance()->GetNPads(uid[0], uid[1])){
857 if (!fAmplitudeHisto){
858 fAmplitudeHisto = new TObjArray(72);
859 }
860 TObjArray * sectorArray = (TObjArray*)fAmplitudeHisto->UncheckedAt(uid[0]);
861 if (!sectorArray){
862 Int_t npads = AliTPCROC::Instance()->GetNChannels(uid[0]);
863 sectorArray = new TObjArray(npads);
864 fAmplitudeHisto->AddAt(sectorArray, uid[0]);
865 }
866 Int_t position = uid[2]+ AliTPCROC::Instance()->GetRowIndexes(uid[0])[uid[1]];
867 TH1F * histo = (TH1F*)sectorArray->UncheckedAt(position);
868 if (!histo){
869 char hname[100];
870 sprintf(hname,"Amp_%d_%d_%d",uid[0],uid[1],uid[2]);
871 TFile * backup = gFile;
872 fDebugStreamer->GetFile()->cd();
873 histo = new TH1F(hname, hname, 100, 1,100);
874 //histo->SetDirectory(0); // histogram not connected to directory -(File)
875 sectorArray->AddAt(histo, position);
876 if (backup) backup->cd();
877 }
878 for (Int_t i=0; i<nchannels; i++){
879 if (signal[i]>0) histo->Fill(signal[i]);
880 }
881 }
882 //
883 TGraph * graph;
884 Bool_t random = (gRandom->Rndm()<0.0001);
885 if (max-median>kMin || rms06>2.*fParam->GetZeroSup() || random){
886 graph =new TGraph(nchannels, dtime, dsignal);
887 if (rms06>2.*fParam->GetZeroSup() || random)
888 (*fDebugStreamer)<<"SignalN"<< //noise pads - or random sample of pads
889 "Sector="<<uid[0]<<
890 "Row="<<uid[1]<<
891 "Pad="<<uid[2]<<
892 "Graph.="<<graph<<
893 "Max="<<max<<
894 "MaxPos="<<maxPos<<
895 //
896 "Median="<<median<<
897 "Mean="<<mean<<
898 "RMS="<<rms<<
899 "Mean06="<<mean06<<
900 "RMS06="<<rms06<<
901 "Mean07="<<mean07<<
902 "RMS07="<<rms07<<
903 "Mean08="<<mean08<<
904 "RMS08="<<rms08<<
905 "Mean09="<<mean09<<
906 "RMS09="<<rms09<<
907 "\n";
908 if (max-median>kMin)
909 (*fDebugStreamer)<<"SignalB"<< // pads with signal
910 "Sector="<<uid[0]<<
911 "Row="<<uid[1]<<
912 "Pad="<<uid[2]<<
913 "Graph.="<<graph<<
914 "Max="<<max<<
915 "MaxPos="<<maxPos<<
916 //
917 "Median="<<median<<
918 "Mean="<<mean<<
919 "RMS="<<rms<<
920 "Mean06="<<mean06<<
921 "RMS06="<<rms06<<
922 "Mean07="<<mean07<<
923 "RMS07="<<rms07<<
924 "Mean08="<<mean08<<
925 "RMS08="<<rms08<<
926 "Mean09="<<mean09<<
927 "RMS09="<<rms09<<
928 "\n";
929 delete graph;
930 }
931
932 (*fDebugStreamer)<<"Signal"<<
933 "Sector="<<uid[0]<<
934 "Row="<<uid[1]<<
935 "Pad="<<uid[2]<<
936 "Max="<<max<<
937 "MaxPos="<<maxPos<<
938 //
939 "Median="<<median<<
940 "Mean="<<mean<<
941 "RMS="<<rms<<
942 "Mean06="<<mean06<<
943 "RMS06="<<rms06<<
944 "Mean07="<<mean07<<
945 "RMS07="<<rms07<<
946 "Mean08="<<mean08<<
947 "RMS08="<<rms08<<
948 "Mean09="<<mean09<<
949 "RMS09="<<rms09<<
950 "\n";
951 delete [] dsignal;
952 delete [] dtime;
953 if (rms06>kMaxNoise) return 1024+median; // sign noisy channel in debug mode
954 return median;
955}
956
957
958
959void AliTPCclustererMI::DumpHistos(){
960 if (!fAmplitudeHisto) return;
961 for (UInt_t isector=0; isector<AliTPCROC::Instance()->GetNSectors(); isector++){
962 TObjArray * array = (TObjArray*)fAmplitudeHisto->UncheckedAt(isector);
963 if (!array) continue;
964 for (UInt_t ipad = 0; ipad <(UInt_t)array->GetEntriesFast(); ipad++){
965 TH1F * histo = (TH1F*) array->UncheckedAt(ipad);
966 if (!histo) continue;
967 if (histo->GetEntries()<100) continue;
968 histo->Fit("gaus","q");
969 Float_t mean = histo->GetMean();
970 Float_t rms = histo->GetRMS();
971 Float_t gmean = histo->GetFunction("gaus")->GetParameter(1);
972 Float_t gsigma = histo->GetFunction("gaus")->GetParameter(2);
973 Float_t max = histo->GetFunction("gaus")->GetParameter(0);
974
975 // get pad number
976 UInt_t row=0, pad =0;
977 const UInt_t *indexes = AliTPCROC::Instance()->GetRowIndexes(isector);
978 for (UInt_t irow=0; irow< AliTPCROC::Instance()->GetNRows(isector); irow++){
979 if (indexes[irow]<ipad){
980 row = irow;
981 pad = ipad-indexes[irow];
982 }
983 }
984 Int_t rpad = pad - (AliTPCROC::Instance()->GetNPads(isector,row))/2;
985 //
986 (*fDebugStreamer)<<"Fit"<<
987 "Sector="<<isector<<
988 "Row="<<row<<
989 "Pad="<<pad<<
990 "RPad="<<rpad<<
991 "Max="<<max<<
992 "Mean="<<mean<<
993 "RMS="<<rms<<
994 "GMean="<<gmean<<
995 "GSigma="<<gsigma<<
996 "\n";
997 if (array->UncheckedAt(ipad)) fDebugStreamer->StoreObject(array->UncheckedAt(ipad));
998 }
999 }
1000}