]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSClusterFinderSDD.cxx
Change needed on alphacxx6
[u/mrichter/AliRoot.git] / ITS / AliITSClusterFinderSDD.cxx
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 /*
16   $Id$
17  */
18 /////////////////////////////////////////////////////////////////////////// 
19 //  Cluster finder                                                       //
20 //  for Silicon                                                          //
21 //  Drift Detector                                                       //
22 ////////////////////////////////////////////////////////////////////////// 
23 #include <Riostream.h>
24
25 #include <TMath.h>
26 #include <math.h>
27
28 #include "AliITS.h"
29 #include "AliITSClusterFinderSDD.h"
30 #include "AliITSMapA1.h"
31 #include "AliITSRawClusterSDD.h"
32 #include "AliITSRecPoint.h"
33 #include "AliITSdigitSDD.h"
34 #include "AliITSresponseSDD.h"
35 #include "AliITSsegmentationSDD.h"
36 #include "AliLog.h"
37 #include "AliRun.h"
38
39 ClassImp(AliITSClusterFinderSDD)
40
41 //______________________________________________________________________
42 AliITSClusterFinderSDD::AliITSClusterFinderSDD():
43 AliITSClusterFinder(),
44 fNclusters(0),
45 fDAnode(0.0),
46 fDTime(0.0),
47 fTimeCorr(0.0),
48 fCutAmplitude(0),
49 fMinPeak(0),
50 fMinCharge(0),
51 fMinNCells(0),
52 fMaxNCells(0){
53     // default constructor
54 }
55 //______________________________________________________________________
56 AliITSClusterFinderSDD::AliITSClusterFinderSDD(AliITSsegmentation *seg,
57                                                AliITSresponse *response,
58                                                TClonesArray *digits,
59                                                TClonesArray *recp):
60 AliITSClusterFinder(seg,response),
61 fNclusters(0),
62 fDAnode(0.0),
63 fDTime(0.0),
64 fTimeCorr(0.0),
65 fCutAmplitude(0),
66 fMinPeak(0),
67 fMinCharge(0),
68 fMinNCells(0),
69 fMaxNCells(0){
70     // standard constructor
71
72     SetDigits(digits);
73     SetClusters(recp);
74     SetCutAmplitude();
75     SetDAnode();
76     SetDTime();
77     SetMinPeak((Int_t)(((AliITSresponseSDD*)GetResp())->
78                        GetNoiseAfterElectronics()*5));
79     //    SetMinPeak();
80     SetMinNCells();
81     SetMaxNCells();
82     SetTimeCorr();
83     SetMinCharge();
84     SetMap(new AliITSMapA1(GetSeg(),Digits(),fCutAmplitude));
85 }
86 //______________________________________________________________________
87 void AliITSClusterFinderSDD::SetCutAmplitude(Double_t nsigma){
88     // set the signal threshold for cluster finder
89     Double_t baseline,noise,noiseAfterEl;
90
91     GetResp()->GetNoiseParam(noise,baseline);
92     noiseAfterEl = ((AliITSresponseSDD*)GetResp())->GetNoiseAfterElectronics();
93     fCutAmplitude = (Int_t)((baseline + nsigma*noiseAfterEl));
94 }
95 //______________________________________________________________________
96 void AliITSClusterFinderSDD::Find1DClusters(){
97     // find 1D clusters
98   
99     // retrieve the parameters 
100     Int_t fNofMaps       = GetSeg()->Npz();
101     Int_t fMaxNofSamples = GetSeg()->Npx();
102     Int_t fNofAnodes     = fNofMaps/2;
103     Int_t dummy          = 0;
104     Double_t fTimeStep    = GetSeg()->Dpx(dummy);
105     Double_t fSddLength   = GetSeg()->Dx();
106     Double_t fDriftSpeed  = GetResp()->DriftSpeed();  
107     Double_t anodePitch   = GetSeg()->Dpz(dummy);
108
109     // map the signal
110     Map()->ClearMap();
111     Map()->SetThreshold(fCutAmplitude);
112     Map()->FillMap();
113   
114     Double_t noise;
115     Double_t baseline;
116     GetResp()->GetNoiseParam(noise,baseline);
117   
118     Int_t nofFoundClusters = 0;
119     Int_t i;
120     Double_t **dfadc = new Double_t*[fNofAnodes];
121     for(i=0;i<fNofAnodes;i++) dfadc[i] = new Double_t[fMaxNofSamples];
122     Double_t fadc  = 0.;
123     Double_t fadc1 = 0.;
124     Double_t fadc2 = 0.;
125     Int_t j,k,idx,l,m;
126     for(j=0;j<2;j++) {
127         for(k=0;k<fNofAnodes;k++) {
128             idx = j*fNofAnodes+k;
129             // signal (fadc) & derivative (dfadc)
130             dfadc[k][255]=0.;
131             for(l=0; l<fMaxNofSamples; l++) {
132                 fadc2=(Double_t)Map()->GetSignal(idx,l);
133                 if(l>0) fadc1=(Double_t)Map()->GetSignal(idx,l-1);
134                 if(l>0) dfadc[k][l-1] = fadc2-fadc1;
135             } // samples
136         } // anodes
137
138         for(k=0;k<fNofAnodes;k++) {
139             AliDebug(5,Form("Anode: %d, Wing: %d",k+1,j+1));
140             idx = j*fNofAnodes+k;
141             Int_t imax  = 0;
142             Int_t imaxd = 0;
143             Int_t it    = 0;
144             while(it <= fMaxNofSamples-3) {
145                 imax  = it;
146                 imaxd = it;
147                 // maximum of signal          
148                 Double_t fadcmax  = 0.;
149                 Double_t dfadcmax = 0.;
150                 Int_t lthrmina   = 1;
151                 Int_t lthrmint   = 3;
152                 Int_t lthra      = 1;
153                 Int_t lthrt      = 0;
154                 for(m=0;m<20;m++) {
155                     Int_t id = it+m;
156                     if(id>=fMaxNofSamples) break;
157                     fadc=(float)Map()->GetSignal(idx,id);
158                     if(fadc > fadcmax) { fadcmax = fadc; imax = id;}
159                     if(fadc > (float)fCutAmplitude)lthrt++; 
160                     if(dfadc[k][id] > dfadcmax) {
161                         dfadcmax = dfadc[k][id];
162                         imaxd    = id;
163                     } // end if
164                 } // end for m
165                 it = imaxd;
166                 if(Map()->TestHit(idx,imax) == kEmpty) {it++; continue;}
167                 // cluster charge
168                 Int_t tstart = it-2;
169                 if(tstart < 0) tstart = 0;
170                 Bool_t ilcl = 0;
171                 if(lthrt >= lthrmint && lthra >= lthrmina) ilcl = 1;
172                 if(ilcl) {
173                     nofFoundClusters++;
174                     Int_t tstop      = tstart;
175                     Double_t dfadcmin = 10000.;
176                     Int_t ij;
177                     for(ij=0; ij<20; ij++) {
178                         if(tstart+ij > 255) { tstop = 255; break; }
179                         fadc=(float)Map()->GetSignal(idx,tstart+ij);
180                         if((dfadc[k][tstart+ij] < dfadcmin) && 
181                            (fadc > fCutAmplitude)) {
182                             tstop = tstart+ij+5;
183                             if(tstop > 255) tstop = 255;
184                             dfadcmin = dfadc[k][it+ij];
185                         } // end if
186                     } // end for ij
187
188                     Double_t clusterCharge = 0.;
189                     Double_t clusterAnode  = k+0.5;
190                     Double_t clusterTime   = 0.;
191                     Int_t   clusterMult   = 0;
192                     Double_t clusterPeakAmplitude = 0.;
193                     Int_t its,peakpos     = -1;
194                     Double_t n, baseline;
195                     GetResp()->GetNoiseParam(n,baseline);
196                     for(its=tstart; its<=tstop; its++) {
197                         fadc=(float)Map()->GetSignal(idx,its);
198                         if(fadc>baseline) fadc -= baseline;
199                         else fadc = 0.;
200                         clusterCharge += fadc;
201                         // as a matter of fact we should take the peak
202                         // pos before FFT
203                         // to get the list of tracks !!!
204                         if(fadc > clusterPeakAmplitude) {
205                             clusterPeakAmplitude = fadc;
206                             //peakpos=Map()->GetHitIndex(idx,its);
207                             Int_t shift = (int)(fTimeCorr/fTimeStep);
208                             if(its>shift && its<(fMaxNofSamples-shift))
209                                 peakpos  = Map()->GetHitIndex(idx,its+shift);
210                             else peakpos = Map()->GetHitIndex(idx,its);
211                             if(peakpos<0) peakpos =Map()->GetHitIndex(idx,its);
212                         } // end if
213                         clusterTime += fadc*its;
214                         if(fadc > 0) clusterMult++;
215                         if(its == tstop) {
216                             clusterTime /= (clusterCharge/fTimeStep);   // ns
217                             if(clusterTime>fTimeCorr) clusterTime -=fTimeCorr;
218                             //ns
219                         } // end if
220                     } // end for its
221
222                     Double_t clusteranodePath = (clusterAnode - fNofAnodes/2)*
223                                                  anodePitch;
224                     Double_t clusterDriftPath = clusterTime*fDriftSpeed;
225                     clusterDriftPath = fSddLength-clusterDriftPath;
226                     if(clusterCharge <= 0.) break;
227                     AliITSRawClusterSDD clust(j+1,//i
228                                               clusterAnode,clusterTime,//ff
229                                               clusterCharge, //f
230                                               clusterPeakAmplitude, //f
231                                               peakpos, //i
232                                               0.,0.,clusterDriftPath,//fff
233                                               clusteranodePath, //f
234                                               clusterMult, //i
235                                               0,0,0,0,0,0,0);//7*i
236                     fITS->AddCluster(1,&clust);
237                     it = tstop;
238                 } // ilcl
239                 it++;
240             } // while (samples)
241         } // anodes
242     } // detectors (2)
243
244     for(i=0;i<fNofAnodes;i++) delete[] dfadc[i];
245     delete [] dfadc;
246
247     return;
248 }
249 //______________________________________________________________________
250 void AliITSClusterFinderSDD::Find1DClustersE(){
251     // find 1D clusters
252     // retrieve the parameters 
253     Int_t fNofMaps = GetSeg()->Npz();
254     Int_t fMaxNofSamples = GetSeg()->Npx();
255     Int_t fNofAnodes = fNofMaps/2;
256     Int_t dummy=0;
257     Double_t fTimeStep = GetSeg()->Dpx( dummy );
258     Double_t fSddLength = GetSeg()->Dx();
259     Double_t fDriftSpeed = GetResp()->DriftSpeed();
260     Double_t anodePitch = GetSeg()->Dpz( dummy );
261     Double_t n, baseline;
262     GetResp()->GetNoiseParam( n, baseline );
263     // map the signal
264     Map()->ClearMap();
265     Map()->SetThreshold( fCutAmplitude );
266     Map()->FillMap();
267     
268     Int_t nClu = 0;
269     //        cout << "Search  cluster... "<< endl;
270     for( Int_t j=0; j<2; j++ ){
271         for( Int_t k=0; k<fNofAnodes; k++ ){
272             Int_t idx = j*fNofAnodes+k;
273             Bool_t on = kFALSE;
274             Int_t start = 0;
275             Int_t nTsteps = 0;
276             Double_t fmax = 0.;
277             Int_t lmax = 0;
278             Double_t charge = 0.;
279             Double_t time = 0.;
280             Double_t anode = k+0.5;
281             Int_t peakpos = -1;
282             for( Int_t l=0; l<fMaxNofSamples; l++ ){
283                 Double_t fadc = (Double_t)Map()->GetSignal( idx, l );
284                 if( fadc > 0.0 ){
285                     if( on == kFALSE && l<fMaxNofSamples-4){
286                         // star RawCluster (reset var.)
287                         Double_t fadc1 = (Double_t)Map()->GetSignal( idx, l+1 );
288                         if( fadc1 < fadc ) continue;
289                         start = l;
290                         fmax = 0.;
291                         lmax = 0;
292                         time = 0.;
293                         charge = 0.; 
294                         on = kTRUE; 
295                         nTsteps = 0;
296                     } // end if on...
297                     nTsteps++ ;
298                     if( fadc > baseline ) fadc -= baseline;
299                     else fadc=0.;
300                     charge += fadc;
301                     time += fadc*l;
302                     if( fadc > fmax ){ 
303                         fmax = fadc; 
304                         lmax = l; 
305                         Int_t shift = (Int_t)(fTimeCorr/fTimeStep + 0.5);
306                         if( l > shift && l < (fMaxNofSamples-shift) )  
307                             peakpos = Map()->GetHitIndex( idx, l+shift );
308                         else
309                             peakpos = Map()->GetHitIndex( idx, l );
310                         if( peakpos < 0) peakpos = Map()->GetHitIndex(idx,l);
311                     } // end if fadc
312                 }else{ // end fadc>0
313                     if( on == kTRUE ){        
314                         if( nTsteps > 2 ){
315                             //  min # of timesteps for a RawCluster
316                             // Found a RawCluster...
317                             Int_t stop = l-1;
318                             time /= (charge/fTimeStep);   // ns
319                                 // time = lmax*fTimeStep;   // ns
320                             if( time > fTimeCorr ) time -= fTimeCorr;   // ns
321                             Double_t anodePath =(anode-fNofAnodes/2)*anodePitch;
322                             Double_t driftPath = time*fDriftSpeed;
323                             driftPath = fSddLength-driftPath;
324                             AliITSRawClusterSDD clust(j+1,anode,time,charge,
325                                                       fmax, peakpos,0.,0.,
326                                                       driftPath,anodePath,
327                                                       nTsteps,start,stop,
328                                                       start, stop, 1, k, k );
329                             fITS->AddCluster( 1, &clust );
330                             if(AliDebugLevel()>=5) clust.PrintInfo();
331                             nClu++;
332                         } // end if nTsteps
333                         on = kFALSE;
334                     } // end if on==kTRUE
335                 } // end if fadc>0
336             } // samples
337         } // anodes
338     } // wings
339     AliDebug(3,Form("# Rawclusters %d",nClu));         
340     return; 
341 }
342 //_______________________________________________________________________
343 Int_t AliITSClusterFinderSDD::SearchPeak(Double_t *spect,Int_t xdim,Int_t zdim,
344                                          Int_t *peakX, Int_t *peakZ, 
345                                          Double_t *peakAmp, Double_t minpeak ){
346     // search peaks on a 2D cluster
347     Int_t npeak = 0;    // # peaks
348     Int_t i,j;
349     // search peaks
350     for( Int_t z=1; z<zdim-1; z++ ){
351         for( Int_t x=1; x<xdim-2; x++ ){
352             Double_t sxz = spect[x*zdim+z];
353             Double_t sxz1 = spect[(x+1)*zdim+z];
354             Double_t sxz2 = spect[(x-1)*zdim+z];
355             // search a local max. in s[x,z]
356             if( sxz < minpeak || sxz1 <= 0 || sxz2 <= 0 ) continue;
357             if( sxz >= spect[(x+1)*zdim+z  ] && sxz >= spect[(x-1)*zdim+z  ] &&
358                 sxz >= spect[x*zdim    +z+1] && sxz >= spect[x*zdim    +z-1] &&
359                 sxz >= spect[(x+1)*zdim+z+1] && sxz >= spect[(x+1)*zdim+z-1] &&
360                 sxz >= spect[(x-1)*zdim+z+1] && sxz >= spect[(x-1)*zdim+z-1] ){
361                 // peak found
362                 peakX[npeak] = x;
363                 peakZ[npeak] = z;
364                 peakAmp[npeak] = sxz;
365                 npeak++;
366             } // end if ....
367         } // end for x
368     } // end for z
369     // search groups of peaks with same amplitude.
370     Int_t *flag = new Int_t[npeak];
371     for( i=0; i<npeak; i++ ) flag[i] = 0;
372     for( i=0; i<npeak; i++ ){
373         for( j=0; j<npeak; j++ ){
374             if( i==j) continue;
375             if( flag[j] > 0 ) continue;
376             if( peakAmp[i] == peakAmp[j] && 
377                 TMath::Abs(peakX[i]-peakX[j])<=1 && 
378                 TMath::Abs(peakZ[i]-peakZ[j])<=1 ){
379                 if( flag[i] == 0) flag[i] = i+1;
380                 flag[j] = flag[i];
381             } // end if ...
382         } // end for j
383     } // end for i
384     // make average of peak groups        
385     for( i=0; i<npeak; i++ ){
386         Int_t nFlag = 1;
387         if( flag[i] <= 0 ) continue;
388         for( j=0; j<npeak; j++ ){
389             if( i==j ) continue;
390             if( flag[j] != flag[i] ) continue;
391             peakX[i] += peakX[j];
392             peakZ[i] += peakZ[j];
393             nFlag++;
394             npeak--;
395             for( Int_t k=j; k<npeak; k++ ){
396                 peakX[k] = peakX[k+1];
397                 peakZ[k] = peakZ[k+1];
398                 peakAmp[k] = peakAmp[k+1];
399                 flag[k] = flag[k+1];
400             } // end for k        
401             j--;
402         } // end for j
403         if( nFlag > 1 ){
404             peakX[i] /= nFlag;
405             peakZ[i] /= nFlag;
406         } // end fi nFlag
407     } // end for i
408     delete [] flag;
409     return( npeak );
410 }
411 //______________________________________________________________________
412 void AliITSClusterFinderSDD::PeakFunc( Int_t xdim, Int_t zdim, Double_t *par,
413                                        Double_t *spe, Double_t *integral){
414     // function used to fit the clusters
415     // par -> parameters..
416     // par[0]  number of peaks.
417     // for each peak i=1, ..., par[0]
418     //                 par[i] = Ampl.
419     //                 par[i+1] = xpos
420     //                 par[i+2] = zpos
421     //                 par[i+3] = tau
422     //                 par[i+4] = sigma.
423     Int_t electronics = GetResp()->Electronics(); // 1 = PASCAL, 2 = OLA
424     const Int_t knParam = 5;
425     Int_t npeak = (Int_t)par[0];
426
427     memset( spe, 0, sizeof( Double_t )*zdim*xdim );
428
429     Int_t k = 1;
430     for( Int_t i=0; i<npeak; i++ ){
431         if( integral != 0 ) integral[i] = 0.;
432         Double_t sigmaA2 = par[k+4]*par[k+4]*2.;
433         Double_t t2 = par[k+3];   // PASCAL
434         if( electronics == 2 ) { t2 *= t2; t2 *= 2; } // OLA
435         for( Int_t z=0; z<zdim; z++ ){
436             for( Int_t x=0; x<xdim; x++ ){
437                 Double_t z2 = (z-par[k+2])*(z-par[k+2])/sigmaA2;
438                 Double_t x2 = 0.;
439                 Double_t signal = 0.;
440                 if( electronics == 1 ){ // PASCAL
441                     x2 = (x-par[k+1]+t2)/t2;
442                     signal = (x2>0.) ? par[k]*x2*exp(-x2+1.-z2) :0.0; // RCCR2
443                 //  signal =(x2>0.) ? par[k]*x2*x2*exp(-2*x2+2.-z2 ):0.0;//RCCR
444                 }else if( electronics == 2 ) { // OLA
445                     x2 = (x-par[k+1])*(x-par[k+1])/t2;
446                     signal = par[k]  * exp( -x2 - z2 );
447                 } else {
448                     Warning("PeakFunc","Wrong SDD Electronics = %d",
449                             electronics);
450                     // exit( 1 );
451                 } // end if electronicx
452                 spe[x*zdim+z] += signal;
453                 if( integral != 0 ) integral[i] += signal;
454             } // end for x
455         } // end for z
456         k += knParam;
457     } // end for i
458     return;
459 }
460 //__________________________________________________________________________
461 Double_t AliITSClusterFinderSDD::ChiSqr( Int_t xdim, Int_t zdim, Double_t *spe,
462                                         Double_t *speFit ) const{
463     // EVALUATES UNNORMALIZED CHI-SQUARED
464     Double_t chi2 = 0.;
465     for( Int_t z=0; z<zdim; z++ ){
466         for( Int_t x=1; x<xdim-1; x++ ){
467             Int_t index = x*zdim+z;
468             Double_t tmp = spe[index] - speFit[index];
469             chi2 += tmp*tmp;
470         } // end for x
471     } // end for z
472     return( chi2 );
473 }
474 //_______________________________________________________________________
475 void AliITSClusterFinderSDD::Minim( Int_t xdim, Int_t zdim, Double_t *param,
476                                     Double_t *prm0,Double_t *steprm,
477                                     Double_t *chisqr,Double_t *spe,
478                                     Double_t *speFit ){
479     // 
480     Int_t   k, nnn, mmm, i;
481     Double_t p1, delta, d1, chisq1, p2, chisq2, t, p3, chisq3, a, b, p0, chisqt;
482     const Int_t knParam = 5;
483     Int_t npeak = (Int_t)param[0];
484     for( k=1; k<(npeak*knParam+1); k++ ) prm0[k] = param[k];
485     for( k=1; k<(npeak*knParam+1); k++ ){
486         p1 = param[k];
487         delta = steprm[k];
488         d1 = delta;
489         // ENSURE THAT STEP SIZE IS SENSIBLY LARGER THAN MACHINE ROUND OFF
490         if( fabs( p1 ) > 1.0E-6 ) 
491             if ( fabs( delta/p1 ) < 1.0E-4 ) delta = p1/1000;
492             else  delta = (Double_t)1.0E-4;
493         //  EVALUATE CHI-SQUARED AT FIRST TWO SEARCH POINTS
494         PeakFunc( xdim, zdim, param, speFit );
495         chisq1 = ChiSqr( xdim, zdim, spe, speFit );
496         p2 = p1+delta;
497         param[k] = p2;
498         PeakFunc( xdim, zdim, param, speFit );
499         chisq2 = ChiSqr( xdim, zdim, spe, speFit );
500         if( chisq1 < chisq2 ){
501             // REVERSE DIRECTION OF SEARCH IF CHI-SQUARED IS INCREASING
502             delta = -delta;
503             t = p1;
504             p1 = p2;
505             p2 = t;
506             t = chisq1;
507             chisq1 = chisq2;
508             chisq2 = t;
509         } // end if
510         i = 1; nnn = 0;
511         do {   // INCREMENT param(K) UNTIL CHI-SQUARED STARTS TO INCREASE
512             nnn++;
513             p3 = p2 + delta;
514             mmm = nnn - (nnn/5)*5;  // multiplo de 5
515             if( mmm == 0 ){
516                 d1 = delta;
517                 // INCREASE STEP SIZE IF STEPPING TOWARDS MINIMUM IS TOO SLOW 
518                 delta *= 5;
519             } // end if
520             param[k] = p3;
521             // Constrain paramiters
522             Int_t kpos = (k-1) % knParam;
523             switch( kpos ){
524             case 0 :
525                 if( param[k] <= 20 ) param[k] = fMinPeak;
526                 break;
527             case 1 :
528                 if( fabs( param[k] - prm0[k] ) > 1.5 ) param[k] = prm0[k];
529                 break;
530             case 2 :
531                 if( fabs( param[k] - prm0[k] ) > 1. ) param[k] = prm0[k];
532                 break;
533             case 3 :
534                 if( param[k] < .5 ) param[k] = .5;        
535                 break;
536             case 4 :
537                 if( param[k] < .288 ) param[k] = .288;// 1/sqrt(12) = 0.288
538                 if( param[k] > zdim*.5 ) param[k] = zdim*.5;
539                 break;
540             }; // end switch
541             PeakFunc( xdim, zdim, param, speFit );
542             chisq3 = ChiSqr( xdim, zdim, spe, speFit );
543             if( chisq3 < chisq2 && nnn < 50 ){
544                 p1 = p2;
545                 p2 = p3;
546                 chisq1 = chisq2;
547                 chisq2 = chisq3;
548             }else i=0;
549         } while( i );
550         // FIND MINIMUM OF PARABOLA DEFINED BY LAST THREE POINTS
551         a = chisq1*(p2-p3)+chisq2*(p3-p1)+chisq3*(p1-p2);
552         b = chisq1*(p2*p2-p3*p3)+chisq2*(p3*p3-p1*p1)+chisq3*(p1*p1-p2*p2);
553         if( a!=0 ) p0 = (Double_t)(0.5*b/a);
554         else p0 = 10000;
555         //--IN CASE OF NEARLY EQUAL CHI-SQUARED AND TOO SMALL STEP SIZE PREVENT
556         //   ERRONEOUS EVALUATION OF PARABOLA MINIMUM
557         //---NEXT TWO LINES CAN BE OMITTED FOR HIGHER PRECISION MACHINES
558         //dp = (Double_t) max (fabs(p3-p2), fabs(p2-p1));
559         //if( fabs( p2-p0 ) > dp ) p0 = p2;
560         param[k] = p0;
561         // Constrain paramiters
562         Int_t kpos = (k-1) % knParam;
563         switch( kpos ){
564         case 0 :
565             if( param[k] <= 20 ) param[k] = fMinPeak;   
566             break;
567         case 1 :
568             if( fabs( param[k] - prm0[k] ) > 1.5 ) param[k] = prm0[k];
569             break;
570         case 2 :
571             if( fabs( param[k] - prm0[k] ) > 1. ) param[k] = prm0[k];
572             break;
573         case 3 :
574             if( param[k] < .5 ) param[k] = .5;        
575             break;
576         case 4 :
577             if( param[k] < .288 ) param[k] = .288;  // 1/sqrt(12) = 0.288
578             if( param[k] > zdim*.5 ) param[k] = zdim*.5;
579             break;
580         }; // end switch
581         PeakFunc( xdim, zdim, param, speFit );
582         chisqt = ChiSqr( xdim, zdim, spe, speFit );
583         // DO NOT ALLOW ERRONEOUS INTERPOLATION
584         if( chisqt <= *chisqr ) *chisqr = chisqt;
585         else param[k] = prm0[k];
586         // OPTIMIZE SEARCH STEP FOR EVENTUAL NEXT CALL OF MINIM
587         steprm[k] = (param[k]-prm0[k])/5;
588         if( steprm[k] >= d1 ) steprm[k] = d1/5;
589     } // end for k
590     // EVALUATE FIT AND CHI-SQUARED FOR OPTIMIZED PARAMETERS
591     PeakFunc( xdim, zdim, param, speFit );
592     *chisqr = ChiSqr( xdim, zdim, spe, speFit );
593     return;
594 }
595 //_________________________________________________________________________
596 Int_t AliITSClusterFinderSDD::NoLinearFit( Int_t xdim, Int_t zdim, 
597                                            Double_t *param, Double_t *spe, 
598                                            Int_t *niter, Double_t *chir ){
599     // fit method from Comput. Phys. Commun 46(1987) 149
600     const Double_t kchilmt = 0.01;  //        relative accuracy           
601     const Int_t   knel = 3;        //        for parabolic minimization  
602     const Int_t   knstop = 50;     //        Max. iteration number          
603     const Int_t   knParam = 5;
604     Int_t npeak = (Int_t)param[0];
605     // RETURN IF NUMBER OF DEGREES OF FREEDOM IS NOT POSITIVE 
606     if( (xdim*zdim - npeak*knParam) <= 0 ) return( -1 );
607     Double_t degFree = (xdim*zdim - npeak*knParam)-1;
608     Int_t   n, k, iterNum = 0;
609     Double_t *prm0 = new Double_t[npeak*knParam+1];
610     Double_t *step = new Double_t[npeak*knParam+1];
611     Double_t *schi = new Double_t[npeak*knParam+1]; 
612     Double_t *sprm[3];
613     sprm[0] = new Double_t[npeak*knParam+1];
614     sprm[1] = new Double_t[npeak*knParam+1];
615     sprm[2] = new Double_t[npeak*knParam+1];
616     Double_t  chi0, chi1, reldif, a, b, prmin, dp;
617     Double_t *speFit = new Double_t[ xdim*zdim ];
618     PeakFunc( xdim, zdim, param, speFit );
619     chi0 = ChiSqr( xdim, zdim, spe, speFit );
620     chi1 = chi0;
621     for( k=1; k<(npeak*knParam+1); k++) prm0[k] = param[k];
622         for( k=1 ; k<(npeak*knParam+1); k+=knParam ){
623             step[k] = param[k] / 20.0 ;
624             step[k+1] = param[k+1] / 50.0;
625             step[k+2] = param[k+2] / 50.0;                 
626             step[k+3] = param[k+3] / 20.0;                 
627             step[k+4] = param[k+4] / 20.0;                 
628         } // end for k
629     Int_t out = 0;
630     do{
631         iterNum++;
632             chi0 = chi1;
633             Minim( xdim, zdim, param, prm0, step, &chi1, spe, speFit );
634             reldif = ( chi1 > 0 ) ? ((Double_t) fabs( chi1-chi0)/chi1 ) : 0;
635         // EXIT conditions
636         if( reldif < (float) kchilmt ){
637             *chir  = (chi1>0) ? (float) TMath::Sqrt (chi1/degFree) :0;
638             *niter = iterNum;
639             out = 0;
640             break;
641         } // end if
642         if( (reldif < (float)(5*kchilmt)) && (iterNum > knstop) ){
643             *chir = (chi1>0) ?(float) TMath::Sqrt (chi1/degFree):0;
644             *niter = iterNum;
645             out = 0;
646             break;
647         } // end if
648         if( iterNum > 5*knstop ){
649             *chir  = (chi1>0) ?(float) TMath::Sqrt (chi1/degFree):0;
650             *niter = iterNum;
651             out = 1;
652             break;
653         } // end if
654         if( iterNum <= knel ) continue;
655         n = iterNum - (iterNum/knel)*knel; // EXTRAPOLATION LIMIT COUNTER N
656         if( n > 3 || n == 0 ) continue;
657         schi[n-1] = chi1;
658         for( k=1; k<(npeak*knParam+1); k++ ) sprm[n-1][k] = param[k];
659         if( n != 3 ) continue;
660         // -EVALUATE EXTRAPOLATED VALUE OF EACH PARAMETER BY FINDING MINIMUM OF
661         //    PARABOLA DEFINED BY LAST THREE CALLS OF MINIM
662         for( k=1; k<(npeak*knParam+1); k++ ){
663             Double_t tmp0 = sprm[0][k];
664             Double_t tmp1 = sprm[1][k];
665             Double_t tmp2 = sprm[2][k];
666             a  = schi[0]*(tmp1-tmp2) + schi[1]*(tmp2-tmp0);
667             a += (schi[2]*(tmp0-tmp1));
668             b  = schi[0]*(tmp1*tmp1-tmp2*tmp2);
669             b += (schi[1]*(tmp2*tmp2-tmp0*tmp0)+(schi[2]*
670                                              (tmp0*tmp0-tmp1*tmp1)));
671             if ((double)a < 1.0E-6) prmin = 0;
672             else prmin = (float) (0.5*b/a);
673             dp = 5*(tmp2-tmp0);
674             if( fabs(prmin-tmp2) > fabs(dp) ) prmin = tmp2+dp;
675             param[k] = prmin;
676             step[k]  = dp/10; // OPTIMIZE SEARCH STEP
677         } // end for k
678     } while( kTRUE );
679     delete [] prm0;
680     delete [] step;
681     delete [] schi; 
682     delete [] sprm[0];
683     delete [] sprm[1];
684     delete [] sprm[2];
685     delete [] speFit;
686     return( out );
687 }
688
689 //______________________________________________________________________
690 void AliITSClusterFinderSDD::ResolveClusters(){
691     // The function to resolve clusters if the clusters overlapping exists
692     Int_t i;
693     // get number of clusters for this module
694     Int_t nofClusters = NClusters();
695     nofClusters -= fNclusters;
696     Int_t fNofMaps = GetSeg()->Npz();
697     Int_t fNofAnodes = fNofMaps/2;
698     //Int_t fMaxNofSamples = GetSeg()->Npx();
699     Int_t dummy=0;
700     Double_t fTimeStep = GetSeg()->Dpx( dummy );
701     Double_t fSddLength = GetSeg()->Dx();
702     Double_t fDriftSpeed = GetResp()->DriftSpeed();
703     Double_t anodePitch = GetSeg()->Dpz( dummy );
704     Double_t n, baseline;
705     GetResp()->GetNoiseParam( n, baseline );
706     Int_t electronics = GetResp()->Electronics(); // 1 = PASCAL, 2 = OLA
707
708     for( Int_t j=0; j<nofClusters; j++ ){ 
709         // get cluster information
710         AliITSRawClusterSDD *clusterJ=(AliITSRawClusterSDD*) Cluster(j);
711         Int_t astart = clusterJ->Astart();
712         Int_t astop = clusterJ->Astop();
713         Int_t tstart = clusterJ->Tstartf();
714         Int_t tstop = clusterJ->Tstopf();
715         Int_t wing = (Int_t)clusterJ->W();
716         if( wing == 2 ){
717             astart += fNofAnodes; 
718             astop  += fNofAnodes;
719         } // end if 
720         Int_t xdim = tstop-tstart+3;
721         Int_t zdim = astop-astart+3;
722         if( xdim > 50 || zdim > 30 ) { 
723             Warning("ResolveClusters","xdim: %d , zdim: %d ",xdim,zdim);
724             continue;
725         }
726         Double_t *sp = new Double_t[ xdim*zdim+1 ];
727         memset( sp, 0, sizeof(Double_t)*(xdim*zdim+1) );
728         
729         // make a local map from cluster region
730         for( Int_t ianode=astart; ianode<=astop; ianode++ ){
731             for( Int_t itime=tstart; itime<=tstop; itime++ ){
732                 Double_t fadc = Map()->GetSignal( ianode, itime );
733                 if( fadc > baseline ) fadc -= (Double_t)baseline;
734                 else fadc = 0.;
735                 Int_t index = (itime-tstart+1)*zdim+(ianode-astart+1);
736                 sp[index] = fadc;
737             } // time loop
738         } // anode loop
739         
740         // search peaks on cluster
741         const Int_t kNp = 150;
742         Int_t peakX1[kNp];
743         Int_t peakZ1[kNp];
744         Double_t peakAmp1[kNp];
745         Int_t npeak = SearchPeak(sp,xdim,zdim,peakX1,peakZ1,peakAmp1,fMinPeak);
746
747         // if multiple peaks, split cluster
748         if( npeak >= 1 ){
749             //        cout << "npeak " << npeak << endl;
750             //        clusterJ->PrintInfo();
751             Double_t *par = new Double_t[npeak*5+1];
752             par[0] = (Double_t)npeak;                
753             // Initial parameters in cell dimentions
754             Int_t k1 = 1;
755             for( i=0; i<npeak; i++ ){
756                 par[k1] = peakAmp1[i];
757                 par[k1+1] = peakX1[i]; // local time pos. [timebin]
758                 par[k1+2] = peakZ1[i]; // local anode pos. [anodepitch]
759                 if( electronics == 1 ) par[k1+3] = 2.; // PASCAL
760                 else if(electronics==2) par[k1+3] = 0.7;//tau [timebin] OLA 
761                 par[k1+4] = .4;    // sigma        [anodepich]
762                 k1 += 5;
763             } // end for i                        
764             Int_t niter;
765             Double_t chir;                        
766             NoLinearFit( xdim, zdim, par, sp, &niter, &chir );
767             Double_t peakX[kNp];
768             Double_t peakZ[kNp];
769             Double_t sigma[kNp];
770             Double_t tau[kNp];
771             Double_t peakAmp[kNp];
772             Double_t integral[kNp];
773             //get integrals => charge for each peak
774             PeakFunc( xdim, zdim, par, sp, integral );
775             k1 = 1;
776             for( i=0; i<npeak; i++ ){
777                 peakAmp[i] = par[k1];
778                 peakX[i]   = par[k1+1];
779                 peakZ[i]   = par[k1+2];
780                 tau[i]     = par[k1+3];
781                 sigma[i]   = par[k1+4];
782                 k1+=5;
783             } // end for i
784             // calculate parameter for new clusters
785             for( i=0; i<npeak; i++ ){
786                 AliITSRawClusterSDD clusterI( *clusterJ );
787             
788                 Int_t newAnode = peakZ1[i]-1 + astart;
789
790             //    Int_t newiTime = peakX1[i]-1 + tstart;
791             //    Int_t shift = (Int_t)(fTimeCorr/fTimeStep + 0.5);
792             //    if( newiTime > shift && newiTime < (fMaxNofSamples-shift) ) 
793             //        shift = 0;
794             //    Int_t peakpos = Map()->GetHitIndex(newAnode,newiTime+shift );
795             //    clusterI.SetPeakPos( peakpos );
796             
797                 clusterI.SetPeakAmpl( peakAmp1[i] );
798                 Double_t newAnodef = peakZ[i] - 0.5 + astart;
799                 Double_t newiTimef = peakX[i] - 1 + tstart;
800                 if( wing == 2 ) newAnodef -= fNofAnodes; 
801                 Double_t anodePath = (newAnodef - fNofAnodes/2)*anodePitch;
802                 newiTimef *= fTimeStep;
803                 if( newiTimef > fTimeCorr ) newiTimef -= fTimeCorr;
804                 if( electronics == 1 ){
805                 //    newiTimef *= 0.999438;    // PASCAL
806                 //    newiTimef += (6./fDriftSpeed - newiTimef/3000.);
807                 }else if( electronics == 2 )
808                     newiTimef *= 0.99714;    // OLA
809                     
810                 Int_t timeBin = (Int_t)(newiTimef/fTimeStep+0.5);    
811                 Int_t peakpos = Map()->GetHitIndex( newAnode, timeBin );
812                 if( peakpos < 0 ) { 
813                     for( Int_t ii=0; ii<3; ii++ ) {
814                         peakpos = Map()->GetHitIndex( newAnode, timeBin+ii );
815                         if( peakpos > 0 ) break;
816                         peakpos = Map()->GetHitIndex( newAnode, timeBin-ii );
817                         if( peakpos > 0 ) break;
818                     }
819                 }
820                 
821                 if( peakpos < 0 ) { 
822                     //Warning("ResolveClusters",
823                     //        "Digit not found for cluster");
824                     //if(AliDebugLevel()>=3) clusterI.PrintInfo(); 
825                    continue;
826                 }
827                 clusterI.SetPeakPos( peakpos );    
828                 Double_t driftPath = fSddLength - newiTimef * fDriftSpeed;
829                 Double_t sign = ( wing == 1 ) ? -1. : 1.;
830                 clusterI.SetX( driftPath*sign * 0.0001 );        
831                 clusterI.SetZ( anodePath * 0.0001 );
832                 clusterI.SetAnode( newAnodef );
833                 clusterI.SetTime( newiTimef );
834                 clusterI.SetAsigma( sigma[i]*anodePitch );
835                 clusterI.SetTsigma( tau[i]*fTimeStep );
836                 clusterI.SetQ( integral[i] );
837                 
838                 fITS->AddCluster( 1, &clusterI );
839             } // end for i
840             Clusters()->RemoveAt( j );
841             delete [] par;
842         } else {  // something odd
843             Warning( "ResolveClusters",
844                      "--- Peak not found!!!!  minpeak=%d ,cluster peak= %f"
845                      " , module= %d",
846                      fMinPeak, clusterJ->PeakAmpl(),GetModule()); 
847             clusterJ->PrintInfo();
848             Warning( "ResolveClusters"," xdim= %d zdim= %d", xdim-2, zdim-2 );
849         }
850         delete [] sp;
851     } // cluster loop
852     Clusters()->Compress();
853 //    Map()->ClearMap(); 
854 }
855 //________________________________________________________________________
856 void  AliITSClusterFinderSDD::GroupClusters(){
857     // group clusters
858     Int_t dummy=0;
859     Double_t fTimeStep = GetSeg()->Dpx(dummy);
860     // get number of clusters for this module
861     Int_t nofClusters = NClusters();
862     nofClusters -= fNclusters;
863     AliITSRawClusterSDD *clusterI;
864     AliITSRawClusterSDD *clusterJ;
865     Int_t *label = new Int_t [nofClusters];
866     Int_t i,j;
867     for(i=0; i<nofClusters; i++) label[i] = 0;
868     for(i=0; i<nofClusters; i++) { 
869         if(label[i] != 0) continue;
870         for(j=i+1; j<nofClusters; j++) { 
871             if(label[j] != 0) continue;
872             clusterI = (AliITSRawClusterSDD*) Cluster(i);
873             clusterJ = (AliITSRawClusterSDD*) Cluster(j);
874             // 1.3 good
875             if(clusterI->T() < fTimeStep*60) fDAnode = 4.2;  // TB 3.2  
876             if(clusterI->T() < fTimeStep*10) fDAnode = 1.5;  // TB 1.
877             Bool_t pair = clusterI->Brother(clusterJ,fDAnode,fDTime);
878             if(!pair) continue;
879             if(AliDebugLevel()>=4){
880                 clusterI->PrintInfo();
881                 clusterJ->PrintInfo();
882             } // end if AliDebugLevel
883             clusterI->Add(clusterJ);
884             label[j] = 1;
885             Clusters()->RemoveAt(j);
886             j=i; // <- Ernesto
887         } // J clusters  
888         label[i] = 1;
889     } // I clusters
890     Clusters()->Compress();
891
892     delete [] label;
893     return;
894 }
895 //________________________________________________________________________
896 void AliITSClusterFinderSDD::SelectClusters(){
897     // get number of clusters for this module
898     Int_t nofClusters = NClusters();
899
900     nofClusters -= fNclusters;
901     Int_t i;
902     for(i=0; i<nofClusters; i++) { 
903         AliITSRawClusterSDD *clusterI =(AliITSRawClusterSDD*) Cluster(i);
904         Int_t rmflg = 0;
905         Double_t wy = 0.;
906         if(clusterI->Anodes() != 0.) {
907             wy = ((Double_t) clusterI->Samples())/clusterI->Anodes();
908         } // end if
909         Int_t amp = (Int_t) clusterI->PeakAmpl();
910         Int_t cha = (Int_t) clusterI->Q();
911         if(amp < fMinPeak) rmflg = 1;  
912         if(cha < fMinCharge) rmflg = 1;
913         if(wy < fMinNCells) rmflg = 1;
914         //if(wy > fMaxNCells) rmflg = 1;
915         if(rmflg) Clusters()->RemoveAt(i);
916     } // I clusters
917     Clusters()->Compress();
918     return;
919 }
920
921 //______________________________________________________________________
922 void AliITSClusterFinderSDD::GetRecPoints(){
923     // get rec points
924   
925     // get number of clusters for this module
926     Int_t nofClusters = NClusters();
927     nofClusters -= fNclusters;
928     const Double_t kconvGeV = 1.e-6; // GeV -> KeV
929     const Double_t kconv = 1.0e-4; 
930     const Double_t kRMSx = 38.0*kconv; // microns->cm ITS TDR Table 1.3
931     const Double_t kRMSz = 28.0*kconv; // microns->cm ITS TDR Table 1.3
932     Int_t i;
933     Int_t ix, iz, idx=-1;
934     AliITSdigitSDD *dig=0;
935     Int_t ndigits=NDigits();
936     for(i=0; i<nofClusters; i++) { 
937         AliITSRawClusterSDD *clusterI = (AliITSRawClusterSDD*)Cluster(i);
938         if(!clusterI) Error("SDD: GetRecPoints","i clusterI ",i,clusterI);
939         if(clusterI) idx=clusterI->PeakPos();
940         if(idx>ndigits) Error("SDD: GetRecPoints","idx ndigits",idx,ndigits);
941         // try peak neighbours - to be done 
942         if(idx&&idx<= ndigits) dig =(AliITSdigitSDD*)GetDigit(idx);
943         if(!dig) {
944             // try cog
945             GetSeg()->GetPadIxz(clusterI->X(),clusterI->Z(),ix,iz);
946             dig = (AliITSdigitSDD*)Map()->GetHit(iz-1,ix-1);
947             // if null try neighbours
948             if (!dig) dig = (AliITSdigitSDD*)Map()->GetHit(iz-1,ix); 
949             if (!dig) dig = (AliITSdigitSDD*)Map()->GetHit(iz-1,ix+1); 
950             if (!dig) printf("SDD: cannot assign the track number!\n");
951         } //  end if !dig
952         AliITSRecPoint rnew;
953         rnew.SetX(clusterI->X());
954         rnew.SetZ(clusterI->Z());
955         rnew.SetQ(clusterI->Q());   // in KeV - should be ADC
956         rnew.SetdEdX(kconvGeV*clusterI->Q());
957         rnew.SetSigmaX2(kRMSx*kRMSx);
958         rnew.SetSigmaZ2(kRMSz*kRMSz);
959
960         if(dig) rnew.fTracks[0]=dig->GetTrack(0);
961         if(dig) rnew.fTracks[1]=dig->GetTrack(1);
962         if(dig) rnew.fTracks[2]=dig->GetTrack(2);
963
964         fITS->AddRecPoint(rnew);
965     } // I clusters
966 //    Map()->ClearMap();
967 }
968 //______________________________________________________________________
969 void AliITSClusterFinderSDD::FindRawClusters(Int_t mod){
970     // find raw clusters
971     
972     SetModule(mod);
973     Find1DClustersE();
974     GroupClusters();
975     SelectClusters();
976     ResolveClusters();
977     GetRecPoints();
978 }
979 //_______________________________________________________________________
980 void AliITSClusterFinderSDD::PrintStatus() const{
981     // Print SDD cluster finder Parameters
982
983     cout << "**************************************************" << endl;
984     cout << " Silicon Drift Detector Cluster Finder Parameters " << endl;
985     cout << "**************************************************" << endl;
986     cout << "Number of Clusters: " << fNclusters << endl;
987     cout << "Anode Tolerance: " << fDAnode << endl;
988     cout << "Time  Tolerance: " << fDTime << endl;
989     cout << "Time  correction (electronics): " << fTimeCorr << endl;
990     cout << "Cut Amplitude (threshold): " << fCutAmplitude << endl;
991     cout << "Minimum Amplitude: " << fMinPeak << endl;
992     cout << "Minimum Charge: " << fMinCharge << endl;
993     cout << "Minimum number of cells/clusters: " << fMinNCells << endl;
994     cout << "Maximum number of cells/clusters: " << fMaxNCells << endl;
995     cout << "**************************************************" << endl;
996 }