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