]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSClusterFinderSDD.cxx
Modify the SetCutAmplitude() function
[u/mrichter/AliRoot.git] / ITS / AliITSClusterFinderSDD.cxx
1
2 /**************************************************************************
3  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  *                                                                        *
5  * Author: The ALICE Off-line Project.                                    *
6  * Contributors are mentioned in the code where appropriate.              *
7  *                                                                        *
8  * Permission to use, copy, modify and distribute this software and its   *
9  * documentation strictly for non-commercial purposes is hereby granted   *
10  * without fee, provided that the above copyright notice appears in all   *
11  * copies and that both the copyright notice and this permission notice   *
12  * appear in the supporting documentation. The authors make no claims     *
13  * about the suitability of this software for any purpose. It is          *
14  * provided "as is" without express or implied warranty.                  *
15  **************************************************************************/
16
17 #include <iostream.h>
18
19 #include <TFile.h>
20 #include <TMath.h>
21 #include <math.h>
22
23 #include "AliITSClusterFinderSDD.h"
24 #include "AliITSMapA1.h"
25 #include "AliITS.h"
26 #include "AliITSdigit.h"
27 #include "AliITSRawCluster.h"
28 #include "AliITSRecPoint.h"
29 #include "AliITSsegmentation.h"
30 #include "AliITSresponseSDD.h"
31 #include "AliRun.h"
32
33
34
35 ClassImp(AliITSClusterFinderSDD)
36
37 //----------------------------------------------------------
38 AliITSClusterFinderSDD::AliITSClusterFinderSDD
39 (AliITSsegmentation *seg, AliITSresponse *response, TClonesArray *digits, TClonesArray *recp)   
40 {
41   // standard constructor
42
43     fSegmentation=seg;
44     fResponse=response;
45     fDigits=digits;
46     fClusters=recp;
47     fNclusters= fClusters->GetEntriesFast();
48     SetCutAmplitude();
49     SetDAnode();
50     SetDTime();
51     SetMinPeak();
52     SetMinNCells();
53     SetMaxNCells();
54     SetTimeCorr();
55     SetMinCharge();
56     fMap=new AliITSMapA1(fSegmentation,fDigits,fCutAmplitude);
57
58 }
59
60 //_____________________________________________________________________________
61 AliITSClusterFinderSDD::AliITSClusterFinderSDD()
62 {
63   // default constructor
64     fSegmentation=0;
65     fResponse=0;
66     fDigits=0;
67     fClusters=0;
68     fNclusters=0;
69     fMap=0;
70     fCutAmplitude=0;
71     SetDAnode();
72     SetDTime();
73     SetMinPeak();
74     SetMinNCells();
75     SetMaxNCells();
76     SetTimeCorr();
77     SetMinCharge();
78
79 }
80
81 //_____________________________________________________________________________
82 AliITSClusterFinderSDD::~AliITSClusterFinderSDD()
83 {
84     // destructor
85
86     if(fMap) delete fMap;
87
88 }
89
90 //_____________________________________________________________________________
91
92 void AliITSClusterFinderSDD::SetCutAmplitude(Float_t nsigma)
93 {
94   // set the signal threshold for cluster finder
95
96    Float_t baseline,noise;
97    fResponse->GetNoiseParam(noise,baseline);
98    Float_t noise_after_el = ((AliITSresponseSDD*)fResponse)->GetNoiseAfterElectronics();
99    fCutAmplitude=(Int_t)(baseline + nsigma*noise_after_el);
100 }
101 //_____________________________________________________________________________
102
103 void AliITSClusterFinderSDD::Find1DClusters()
104 {
105   // find 1D clusters
106
107   AliITS *iTS=(AliITS*)gAlice->GetModule("ITS");
108   
109   // retrieve the parameters 
110   Int_t fNofMaps = fSegmentation->Npz();
111   Int_t fMaxNofSamples = fSegmentation->Npx();
112   Int_t fNofAnodes = fNofMaps/2;
113   Int_t dummy=0;
114   Float_t fTimeStep = fSegmentation->Dpx(dummy);
115   Float_t fSddLength = fSegmentation->Dx();
116   Float_t fDriftSpeed = fResponse->DriftSpeed();
117   
118   Float_t anodePitch = fSegmentation->Dpz(dummy);
119   // map the signal
120   fMap->SetThreshold(fCutAmplitude);
121
122   fMap->FillMap();
123   
124   Float_t noise;
125   Float_t baseline;
126   fResponse->GetNoiseParam(noise,baseline);
127   
128   Int_t nofFoundClusters = 0;
129   Int_t i;
130   Float_t **dfadc = new Float_t*[fNofAnodes];
131   for(i=0;i<fNofAnodes;i++) dfadc[i] = new Float_t[fMaxNofSamples];
132   Float_t fadc = 0.;
133   Float_t fadc1 = 0.;
134   Float_t fadc2 = 0.;
135   Int_t j,k,idx,l,m;
136   for(j=0;j<2;j++) {
137     for(k=0;k<fNofAnodes;k++) {
138       idx = j*fNofAnodes+k;
139       // signal (fadc) & derivative (dfadc)
140       dfadc[k][255]=0.;
141       for(l=0; l<fMaxNofSamples; l++) {
142         fadc2=(Float_t)fMap->GetSignal(idx,l);
143         if(l>0) fadc1=(Float_t)fMap->GetSignal(idx,l-1);
144         if(l>0) dfadc[k][l-1] = fadc2-fadc1;
145       } // samples
146     } // anodes
147     
148     for(k=0;k<fNofAnodes;k++) {
149       //cout << "Anode: " << k+1 << ", Wing: " << j+1 << endl;
150       idx = j*fNofAnodes+k;
151       
152       Int_t imax = 0;
153       Int_t imaxd = 0;
154       Int_t it=0;
155       while(it <= fMaxNofSamples-3) {
156         
157         imax = it;
158         imaxd = it;
159         // maximum of signal      
160         
161         Float_t fadcmax = 0.;
162         Float_t dfadcmax = 0.;
163         Int_t lthrmina = 1;
164         Int_t lthrmint = 3;
165         
166         Int_t lthra = 1;
167         Int_t lthrt = 0;
168         
169         for(m=0;m<20;m++) {
170           Int_t id = it+m;
171           if(id>=fMaxNofSamples) break;
172           fadc=(float)fMap->GetSignal(idx,id);
173           if(fadc > fadcmax) { fadcmax = fadc; imax = id;}
174           if(fadc > (float)fCutAmplitude) { 
175             lthrt++; 
176           }
177           
178           if(dfadc[k][id] > dfadcmax) {
179             dfadcmax = dfadc[k][id];
180             imaxd = id;
181           }
182         }
183         it = imaxd;
184         
185         if(fMap->TestHit(idx,imax) == kEmpty) {it++; continue;}
186         
187         // cluster charge
188         Int_t tstart = it-2;
189         if(tstart < 0) tstart = 0;
190         
191         Bool_t ilcl = 0;
192         if(lthrt >= lthrmint && lthra >= lthrmina) ilcl = 1;
193
194         if(ilcl) {
195           nofFoundClusters++;
196           
197           Int_t tstop = tstart;
198           Float_t dfadcmin = 10000.;
199           Int_t ij;
200           for(ij=0; ij<20; ij++) {
201             if(tstart+ij > 255) { tstop = 255; break; }
202             fadc=(float)fMap->GetSignal(idx,tstart+ij);
203             if((dfadc[k][tstart+ij] < dfadcmin) && (fadc > fCutAmplitude)) {
204               tstop = tstart+ij+5;
205               if(tstop > 255) tstop = 255;
206               dfadcmin = dfadc[k][it+ij];
207             }
208           }
209           
210           Float_t clusterCharge = 0.;
211           Float_t clusterAnode = k+0.5;
212           Float_t clusterTime = 0.;
213           Float_t clusterMult = 0.;
214           Float_t clusterPeakAmplitude = 0.;
215           Int_t its,peakpos=-1;
216           Float_t n, baseline;
217           fResponse->GetNoiseParam(n,baseline);
218           for(its=tstart; its<=tstop; its++) {
219             fadc=(float)fMap->GetSignal(idx,its);
220             if(fadc>baseline)
221               fadc-=baseline;
222             else
223               fadc=0.;
224             clusterCharge += fadc;
225             // as a matter of fact we should take the peak pos before FFT
226             // to get the list of tracks !!!
227             if(fadc > clusterPeakAmplitude) {
228               clusterPeakAmplitude = fadc;
229               //peakpos=fMap->GetHitIndex(idx,its);
230               Int_t shift=(int)(fTimeCorr/fTimeStep);
231               if(its>shift && its<(fMaxNofSamples-shift)) peakpos=fMap->GetHitIndex(idx,its+shift);
232               else peakpos=fMap->GetHitIndex(idx,its);
233               if(peakpos<0) peakpos=fMap->GetHitIndex(idx,its);
234             }
235             clusterTime += fadc*its;
236             if(fadc > 0) clusterMult++;
237             if(its == tstop) {
238               clusterTime /= (clusterCharge/fTimeStep);   // ns
239               if(clusterTime > fTimeCorr) clusterTime -= fTimeCorr;   // ns
240             }
241           }
242           
243           Float_t clusteranodePath = (clusterAnode - fNofAnodes/2)*anodePitch;
244           Float_t clusterDriftPath = clusterTime*fDriftSpeed;
245           clusterDriftPath = fSddLength-clusterDriftPath;
246
247           if(clusterCharge <= 0.) break;
248           AliITSRawClusterSDD clust(j+1,clusterAnode,clusterTime,clusterCharge,clusterPeakAmplitude,peakpos,0.,0.,clusterDriftPath,clusteranodePath,clusterMult,0,0,0,0,0,0,0);
249           iTS->AddCluster(1,&clust);
250           it = tstop;
251         } // ilcl
252         
253         it++;
254         
255       } // while (samples)
256     } // anodes
257   } // detectors (2)
258   
259   //fMap->ClearMap(); 
260   
261   for(i=0;i<fNofAnodes;i++) delete[] dfadc[i];
262   delete [] dfadc;
263   
264   return;
265   
266 }
267
268 //_____________________________________________________________________________
269
270 void AliITSClusterFinderSDD::Find1DClustersE()
271 {
272     // find 1D clusters
273
274         AliITS *iTS=(AliITS*)gAlice->GetModule("ITS");
275
276         // retrieve the parameters 
277         Int_t fNofMaps = fSegmentation->Npz();
278         Int_t fMaxNofSamples = fSegmentation->Npx();
279         Int_t fNofAnodes = fNofMaps/2;
280         Int_t dummy=0;
281         Float_t fTimeStep = fSegmentation->Dpx( dummy );
282         Float_t fSddLength = fSegmentation->Dx();
283         Float_t fDriftSpeed = fResponse->DriftSpeed();
284         Float_t anodePitch = fSegmentation->Dpz( dummy );
285         Float_t n, baseline;
286         fResponse->GetNoiseParam( n, baseline );
287
288         // map the signal
289         fMap->SetThreshold( fCutAmplitude );
290         fMap->FillMap();
291         Int_t nClu = 0;
292         
293 //      cout << "Search  cluster... "<< endl;
294         for( Int_t j=0; j<2; j++ ) 
295         {
296         for( Int_t k=0; k<fNofAnodes; k++ ) 
297                 {
298                 Int_t idx = j*fNofAnodes+k;
299                         
300                         Bool_t on = kFALSE;
301                         Int_t start = 0;
302                         Int_t nTsteps = 0;
303                         Float_t fmax = 0.;
304                         Int_t lmax = 0;
305                         Float_t charge = 0.;
306                         Float_t time = 0.;
307                         Float_t anode = k+0.5;
308                         Int_t peakpos = -1;
309
310                 for( Int_t l=0; l<fMaxNofSamples; l++ ) 
311                         {
312                                 Float_t fadc = (Float_t)fMap->GetSignal( idx, l );
313                                 if( fadc > 0.0 )
314                                 {
315                                         if( on == kFALSE && l<fMaxNofSamples-4)  // star RawCluster (reset var.)
316                                         { 
317                                                 Float_t fadc1 = (Float_t)fMap->GetSignal( idx, l+1 );
318                                                 if( fadc1 < fadc ) continue;
319                                                 start = l;
320                                                 fmax = 0.;
321                                                 lmax = 0;
322                                                 time = 0.;
323                                                 charge = 0.; 
324                                                 on = kTRUE; 
325                                                 nTsteps = 0;
326                                         }
327                                         
328                                         nTsteps++ ;
329                                         if( fadc > baseline )
330                                         fadc -= baseline;
331                                 else
332                                         fadc=0.;
333
334                                 charge += fadc;
335                                         time += fadc*l;
336
337                                         if( fadc > fmax ) 
338                                         { 
339                                                 fmax = fadc; 
340                                                 lmax = l; 
341                                                 Int_t shift = (Int_t)(fTimeCorr/fTimeStep + 0.5);
342                             if( l > shift && l < (fMaxNofSamples-shift) )  
343                                                         peakpos = fMap->GetHitIndex( idx, l+shift );
344                                         else
345                                                         peakpos = fMap->GetHitIndex( idx, l );
346                                             if( peakpos < 0 ) peakpos = fMap->GetHitIndex( idx, l );
347                                         }
348                                 }
349                                 else
350                                 {
351                                         if( on == kTRUE )
352                                         {       
353                                                 if( nTsteps > 2 ) //  min # of timesteps for a RawCluster
354                                                 {
355                                                         // Found a RawCluster...
356                                                         Int_t stop = l-1;
357                                                         time /= (charge/fTimeStep);   // ns
358                                 //                      time = lmax*fTimeStep;   // ns
359                                                         if( time > fTimeCorr ) time -= fTimeCorr;   // ns
360                                                         Float_t anodePath = (anode - fNofAnodes/2)*anodePitch;
361                                                         Float_t driftPath = time*fDriftSpeed;
362                                                         driftPath = fSddLength-driftPath;
363                                                         AliITSRawClusterSDD clust( j+1, anode, time, charge,
364                                     fmax, peakpos, 0., 0., driftPath, anodePath, nTsteps
365                                                         , start, stop, start, stop, 1, k, k );
366                                                         iTS->AddCluster( 1, &clust );
367                                                 //      clust.PrintInfo();
368                                                         nClu++;
369                                                 }
370                                                 on = kFALSE;
371                                         }
372                                 }
373                 } // samples
374         } // anodes
375         } // wings
376
377 //      cout << "# Rawclusters " << nClu << endl;       
378         return; 
379
380 }
381
382 //_____________________________________________________________________________
383
384 Int_t AliITSClusterFinderSDD::SearchPeak( Float_t *spect, Int_t xdim, Int_t zdim, 
385                                                                                   Int_t *peakX, Int_t *peakZ, Float_t *peakAmp, Float_t minpeak )
386 {
387         // search peaks on a 2D cluster
388         Int_t npeak = 0;    // # peaks
389     Int_t i,j;
390
391         // search peaks
392         for( Int_t z=1; z<zdim-1; z++ )
393         {
394                 for( Int_t x=2; x<xdim-3; x++ )
395                 {
396                         Float_t sxz = spect[x*zdim+z];
397                         Float_t sxz1 = spect[(x+1)*zdim+z];
398                         Float_t sxz2 = spect[(x-1)*zdim+z];
399
400                         // search a local max. in s[x,z]
401                         if( sxz < minpeak || sxz1 <= 0 || sxz2 <= 0 ) continue;
402                         if( sxz >= spect[(x+1)*zdim+z  ] && sxz >= spect[(x-1)*zdim+z  ] &&
403                                 sxz >= spect[x*zdim    +z+1] && sxz >= spect[x*zdim    +z-1] &&
404                                 sxz >= spect[(x+1)*zdim+z+1] && sxz >= spect[(x+1)*zdim+z-1] &&
405                                 sxz >= spect[(x-1)*zdim+z+1] && sxz >= spect[(x-1)*zdim+z-1] )
406                         {
407                                 // peak found
408                                 peakX[npeak] = x;
409                                 peakZ[npeak] = z;
410                                 peakAmp[npeak] = sxz;
411                                 npeak++;
412                         }
413                 }
414         }                       
415
416         // search groups of peaks with same amplitude.
417         Int_t *flag = new Int_t[npeak];
418         for( i=0; i<npeak; i++ ) flag[i] = 0;
419         for( i=0; i<npeak; i++ )
420         {
421                 for( j=0; j<npeak; j++ )
422                 {
423                         if( i==j) continue;
424                         if( flag[j] > 0 ) continue;
425                         if( peakAmp[i] == peakAmp[j] && TMath::Abs(peakX[i]-peakX[j])<=1 && TMath::Abs(peakZ[i]-peakZ[j])<=1 )
426                         {
427                                 if( flag[i] == 0) flag[i] = i+1;
428                                 flag[j] = flag[i];
429                         }
430                 }
431         }
432
433         // make average of peak groups  
434         for( i=0; i<npeak; i++ )
435         {
436                 Int_t nFlag = 1;
437                 if( flag[i] <= 0 ) continue;
438                 for( j=0; j<npeak; j++ )
439                 {
440                         if( i==j ) continue;
441                         if( flag[j] != flag[i] ) continue;
442                         peakX[i] += peakX[j];
443                         peakZ[i] += peakZ[j];
444                         nFlag++;
445                         npeak--;
446                         for( Int_t k=j; k<npeak; k++ )
447                         {
448                                 peakX[k] = peakX[k+1];
449                                 peakZ[k] = peakZ[k+1];
450                                 peakAmp[k] = peakAmp[k+1];
451                                 flag[k] = flag[k+1];
452                         }       
453                         j--;
454                 }
455
456                 if( nFlag > 1 )
457                 {
458                         peakX[i] /= nFlag;
459                         peakZ[i] /= nFlag;
460                 }
461         }
462         
463         delete [] flag;
464         return( npeak );
465 }
466
467
468 void AliITSClusterFinderSDD::PeakFunc( Int_t xdim, Int_t zdim, Float_t *par, Float_t *spe, Float_t
469 *integral) 
470 {
471     // function used to fit the clusters
472     // par -> paramiters..
473     // par[0]  number of peaks.
474     // for each peak i=1, ..., par[0]
475     //          par[i] = Ampl.
476     //          par[i+1] = xpos
477     //          par[i+2] = zpos
478     //          par[i+3] = tau
479     //          par[i+4] = sigma.
480
481     Int_t electronics = fResponse->Electronics(); // 1 = PASCAL, 2 = OLA
482     const Int_t knParam = 5;
483     Int_t npeak = (Int_t)par[0];
484   
485     memset( spe, 0, sizeof( Float_t )*zdim*xdim );
486   
487     Int_t k = 1;
488     for( Int_t i=0; i<npeak; i++ )
489     {
490         if( integral != 0 ) integral[i] = 0.;
491         Float_t sigmaA2 = par[k+4]*par[k+4]*2.;
492         Float_t T2 = par[k+3];   // PASCAL
493         if( electronics == 2 ) { T2 *= T2; T2 *= 2; } // OLA
494         for( Int_t z=0; z<zdim; z++ )
495         {
496             for( Int_t x=0; x<xdim; x++ )
497             {
498                 Float_t z2 = (z-par[k+2])*(z-par[k+2])/sigmaA2;
499                 Float_t x2 = 0.;
500                 Float_t signal = 0.;
501                 if( electronics == 1 ) // PASCAL
502                 {
503                     x2 = (x-par[k+1]+T2)/T2;
504                     signal = (x2 > 0.) ? par[k] * x2 * exp( -x2+1. - z2 ) : 0.0;
505                 //  signal = (x2 > 0.) ? par[k] * x2*x2 * exp( -2*x2+2. - z2 ) : 0.0; // RCCR
506                 }
507                 else
508                     if( electronics == 2 ) // OLA
509                     {
510                         x2 = (x-par[k+1])*(x-par[k+1])/T2;
511                         signal = par[k]  * exp( -x2 - z2 );
512                     }
513                     else
514                     {
515                         cout << "Wrong SDD Electronics =" << electronics << endl;
516                         // exit( 1 );
517                     }
518                 spe[x*zdim+z] += signal;
519                 if( integral != 0 ) integral[i] += signal;
520             }
521         }
522         k += knParam;
523     }
524     return;
525 }
526
527
528 Float_t AliITSClusterFinderSDD::ChiSqr( Int_t xdim, Int_t zdim, Float_t *spe, Float_t *speFit )
529 {
530         // EVALUATES UNNORMALIZED CHI-SQUARED
531         
532         Float_t chi2 = 0.;
533         for( Int_t z=0; z<zdim; z++ )
534         {
535                 for( Int_t x=1; x<xdim-1; x++ )
536                 {
537                         Int_t index = x*zdim+z;
538                         Float_t tmp = spe[index] - speFit[index];
539                         chi2 += tmp*tmp;
540                 }
541         }       
542         return( chi2 );
543 }
544
545
546 void AliITSClusterFinderSDD::Minim( Int_t xdim, Int_t zdim, Float_t *param, Float_t *prm0, Float_t *steprm, Float_t *chisqr, 
547                 Float_t *spe, Float_t *speFit )
548 {
549         // 
550         Int_t   k, nnn, mmm, i;
551         Float_t p1, delta, d1, chisq1, p2, chisq2, t, p3, chisq3, a, b, p0, chisqt;
552         
553         const Int_t knParam = 5;
554         Int_t npeak = (Int_t)param[0];
555         for( k=1; k<(npeak*knParam+1); k++ ) prm0[k] = param[k];
556
557         for( k=1; k<(npeak*knParam+1); k++ ) 
558         {
559                 p1 = param[k];
560                 delta = steprm[k];
561                 d1 = delta;
562
563                 // ENSURE THAT STEP SIZE IS SENSIBLY LARGER THAN MACHINE ROUND OFF
564                 if( fabs( p1 ) > 1.0E-6 ) 
565                         if ( fabs( delta/p1 ) < 1.0E-4 ) delta = p1/1000;
566                 else  delta = (Float_t)1.0E-4;
567
568                 //  EVALUATE CHI-SQUARED AT FIRST TWO SEARCH POINTS
569                 PeakFunc( xdim, zdim, param, speFit );
570                 chisq1 = ChiSqr( xdim, zdim, spe, speFit );
571
572                 p2 = p1+delta;
573                 param[k] = p2;
574
575                 PeakFunc( xdim, zdim, param, speFit );
576                 chisq2 = ChiSqr( xdim, zdim, spe, speFit );
577
578                 if( chisq1 < chisq2 ) 
579                 {
580         // REVERSE DIRECTION OF SEARCH IF CHI-SQUARED IS INCREASING
581                 delta = -delta;
582                 t = p1;
583                 p1 = p2;
584                 p2 = t;
585                 t = chisq1;
586                 chisq1 = chisq2;
587                 chisq2 = t;
588         }
589
590                 i = 1; nnn = 0;
591                 do {   // INCREMENT param(K) UNTIL CHI-SQUARED STARTS TO INCREASE
592                 nnn++;
593                         p3 = p2 + delta;
594                 mmm = nnn - (nnn/5)*5;  // multiplo de 5
595                 if( mmm == 0 ) 
596                         {
597                         d1 = delta;
598                         // INCREASE STEP SIZE IF STEPPING TOWARDS MINIMUM IS TOO SLOW 
599                         delta *= 5;
600                 }
601                 param[k] = p3;
602
603                         // Constrain paramiters
604                         Int_t kpos = (k-1) % knParam;
605                         switch( kpos ) 
606                         {
607                                 case 0 :
608                                         if( param[k] <= 20 ) param[k] = fMinPeak;
609                                 case 1 :
610                                         if( fabs( param[k] - prm0[k] ) > 1.5 ) param[k] = prm0[k];
611                                 case 2 :
612                                         if( fabs( param[k] - prm0[k] ) > 1. ) param[k] = prm0[k];
613                                 case 3 :
614                                         if( param[k] < .5 ) param[k] = .5;      
615                                 case 4 :
616                                         if( param[k] < .288 ) param[k] = .288;  // 1/sqrt(12) = 0.288
617                         };
618         
619                         PeakFunc( xdim, zdim, param, speFit );
620                         chisq3 = ChiSqr( xdim, zdim, spe, speFit );
621                 
622                 if( chisq3 < chisq2 && nnn < 50 ) 
623                         {
624                                 p1 = p2;
625                         p2 = p3;
626                         chisq1 = chisq2;
627                         chisq2 = chisq3;
628                 }
629                 else i=0;
630
631         } while( i );
632
633                 // FIND MINIMUM OF PARABOLA DEFINED BY LAST THREE POINTS
634                 a = chisq1*(p2-p3)+chisq2*(p3-p1)+chisq3*(p1-p2);
635                 b = chisq1*(p2*p2-p3*p3)+chisq2*(p3*p3-p1*p1)+chisq3*(p1*p1-p2*p2);
636                 if( a!=0 ) p0 = (Float_t)(0.5*b/a);
637                   else p0 = 10000;
638
639                 //---IN CASE OF NEARLY EQUAL CHI-SQUARED AND TOO SMALL STEP SIZE PREVENT
640                 //   ERRONEOUS EVALUATION OF PARABOLA MINIMUM
641                 //---NEXT TWO LINES CAN BE OMITTED FOR HIGHER PRECISION MACHINES
642
643                 //dp = (Float_t) max (fabs(p3-p2), fabs(p2-p1));
644                 //if( fabs( p2-p0 ) > dp ) p0 = p2;
645                 param[k] = p0;
646                 
647                 // Constrain paramiters
648                 Int_t kpos = (k-1) % knParam;
649                 switch( kpos ) 
650                 {
651                         case 0 :
652                                 if( param[k] <= 20 ) param[k] = fMinPeak;   
653                         case 1 :
654                                 if( fabs( param[k] - prm0[k] ) > 1.5 ) param[k] = prm0[k];
655                         case 2 :
656                                 if( fabs( param[k] - prm0[k] ) > 1. ) param[k] = prm0[k];
657                         case 3 :
658                                 if( param[k] < .5 ) param[k] = .5;      
659                         case 4 :
660                                 if( param[k] < .288 ) param[k] = .288;  // 1/sqrt(12) = 0.288   
661                 };
662         
663                 PeakFunc( xdim, zdim, param, speFit );
664                 chisqt = ChiSqr( xdim, zdim, spe, speFit );
665
666                 // DO NOT ALLOW ERRONEOUS INTERPOLATION
667                 if( chisqt <= *chisqr ) 
668                         *chisqr = chisqt;
669                 else  
670                         param[k] = prm0[k];
671
672                 // OPTIMIZE SEARCH STEP FOR EVENTUAL NEXT CALL OF MINIM
673                 steprm[k] = (param[k]-prm0[k])/5;
674                 if( steprm[k] >= d1 ) steprm[k] = d1/5;
675         }
676
677         // EVALUATE FIT AND CHI-SQUARED FOR OPTIMIZED PARAMETERS
678         PeakFunc( xdim, zdim, param, speFit );
679         *chisqr = ChiSqr( xdim, zdim, spe, speFit );
680         return;
681 }
682
683
684 Int_t AliITSClusterFinderSDD::NoLinearFit( Int_t xdim, Int_t zdim, Float_t *param, Float_t *spe, Int_t *niter, Float_t *chir )
685 {
686         // fit method from Comput. Phys. Commun 46(1987) 149
687         const Float_t kchilmt = 0.01;              //   relative accuracy         
688         const Int_t   knel = 3;                            //   for parabolic minimization  
689         const Int_t   knstop = 50;                         //   Max. iteration number     
690         const Int_t   knParam = 5;
691
692         Int_t npeak = (Int_t)param[0];
693         
694         // RETURN IF NUMBER OF DEGREES OF FREEDOM IS NOT POSITIVE 
695         if( (xdim*zdim - npeak*knParam) <= 0 ) return( -1 );
696         Float_t degFree = (xdim*zdim - npeak*knParam)-1;
697
698         Int_t   n, k, iterNum = 0;
699         Float_t *prm0 = new Float_t[npeak*knParam+1];
700         Float_t *step = new Float_t[npeak*knParam+1];
701         Float_t *schi = new Float_t[npeak*knParam+1]; 
702         Float_t *sprm[3];
703         sprm[0] = new Float_t[npeak*knParam+1];
704         sprm[1] = new Float_t[npeak*knParam+1];
705         sprm[2] = new Float_t[npeak*knParam+1];
706         
707         Float_t  chi0, chi1, reldif, a, b, prmin, dp;
708         
709         Float_t *speFit = new Float_t[ xdim*zdim ];
710         PeakFunc( xdim, zdim, param, speFit );
711         chi0 = ChiSqr( xdim, zdim, spe, speFit );
712         chi1 = chi0;
713
714
715         for( k=1; k<(npeak*knParam+1); k++) prm0[k] = param[k];
716
717         for( k=1 ; k<(npeak*knParam+1); k+=knParam ) 
718         {
719         step[k] = param[k] / 20.0 ;             
720         step[k+1] = param[k+1] / 50.0;
721                 step[k+2] = param[k+2] / 50.0;           
722         step[k+3] = param[k+3] / 20.0;           
723         step[k+4] = param[k+4] / 20.0;           
724     }
725
726         Int_t out = 0;
727         do 
728         {
729             iterNum++;
730         chi0 = chi1;
731                 
732         Minim( xdim, zdim, param, prm0, step, &chi1, spe, speFit );
733         reldif = ( chi1 > 0 ) ? ((Float_t) fabs( chi1-chi0)/chi1 ) : 0;
734
735             // EXIT conditions
736                 if( reldif < (float) kchilmt )  
737                 {
738                         *chir  = (chi1>0) ? (float) TMath::Sqrt (chi1/degFree) :0;
739                 *niter = iterNum;
740                         out = 0;
741                         break;
742         }
743
744         if( (reldif < (float)(5*kchilmt)) && (iterNum > knstop) ) 
745                 {
746                 *chir = (chi1>0) ?(float) TMath::Sqrt (chi1/degFree):0;
747                 *niter = iterNum;
748                         out = 0;
749                         break;
750         }
751
752         if( iterNum > 5*knstop ) 
753                 {
754                 *chir  = (chi1>0) ?(float) TMath::Sqrt (chi1/degFree):0;
755                 *niter = iterNum;
756                         out = 1;
757                         break;
758         }
759
760         if( iterNum <= knel ) continue;
761
762         n = iterNum - (iterNum/knel)*knel; // EXTRAPOLATION LIMIT COUNTER N
763         if( n > 3 || n == 0 ) continue;
764         schi[n-1] = chi1;
765         for( k=1; k<(npeak*knParam+1); k++ ) sprm[n-1][k] = param[k];
766         if( n != 3 ) continue;
767
768                 //   -EVALUATE EXTRAPOLATED VALUE OF EACH PARAMETER BY FINDING MINIMUM OF
769                 //    PARABOLA DEFINED BY LAST THREE CALLS OF MINIM
770
771         for( k=1; k<(npeak*knParam+1); k++ ) 
772                 {
773                 Float_t tmp0 = sprm[0][k];
774                 Float_t tmp1 = sprm[1][k];
775                 Float_t tmp2 = sprm[2][k];
776                 a  = schi[0]*(tmp1-tmp2) + schi[1]*(tmp2-tmp0);
777                 a += (schi[2]*(tmp0-tmp1));
778                 b  = schi[0]*(tmp1*tmp1-tmp2*tmp2);
779                 b += (schi[1]*(tmp2*tmp2-tmp0*tmp0)+(schi[2]*(tmp0*tmp0-tmp1*tmp1)));
780                 if ((double)a < 1.0E-6) prmin = 0;
781                            else prmin = (float) (0.5*b/a);
782                 dp = 5*(tmp2-tmp0);
783
784                 if (fabs(prmin-tmp2) > fabs(dp)) prmin = tmp2+dp;
785                 param[k] = prmin;
786                 step[k]  = dp/10; // OPTIMIZE SEARCH STEP
787         }
788
789         } while( kTRUE );
790
791         delete [] prm0;
792         delete [] step;
793         delete [] schi; 
794         delete [] sprm[0];
795         delete [] sprm[1];
796         delete [] sprm[2];
797         delete [] speFit;
798
799         return( out );
800 }
801 //_____________________________________________________________________________
802 void AliITSClusterFinderSDD::ResolveClustersE()
803 {
804         // The function to resolve clusters if the clusters overlapping exists
805
806     Int_t i;
807
808         AliITS *iTS = (AliITS*)gAlice->GetModule( "ITS" );
809         // get number of clusters for this module
810         Int_t nofClusters = fClusters->GetEntriesFast();
811         nofClusters -= fNclusters;
812
813         Int_t fNofMaps = fSegmentation->Npz();
814         Int_t fNofAnodes = fNofMaps/2;
815         Int_t fMaxNofSamples = fSegmentation->Npx();
816         Int_t dummy=0;
817         Double_t fTimeStep = fSegmentation->Dpx( dummy );
818         Double_t fSddLength = fSegmentation->Dx();
819         Double_t fDriftSpeed = fResponse->DriftSpeed();
820         Double_t anodePitch = fSegmentation->Dpz( dummy );
821         Float_t n, baseline;
822         fResponse->GetNoiseParam( n, baseline );
823         Int_t electronics = fResponse->Electronics(); // 1 = PASCAL, 2 = OLA
824         
825         // fill Map of signals
826         fMap->FillMap(); 
827
828         for( Int_t j=0; j<nofClusters; j++ ) 
829         { 
830                 // get cluster information
831                 AliITSRawClusterSDD *clusterJ = (AliITSRawClusterSDD*) fClusters->At( j );
832                 Int_t astart = clusterJ->Astart();
833                 Int_t astop = clusterJ->Astop();
834                 Int_t tstart = clusterJ->Tstartf();
835                 Int_t tstop = clusterJ->Tstopf();
836                 Int_t wing = (Int_t)clusterJ->W();
837                 if( wing == 2 ) 
838                 {
839                         astart += fNofAnodes; 
840                         astop  += fNofAnodes;
841                 } 
842                 Int_t xdim = tstop-tstart+3;
843                 Int_t zdim = astop-astart+3;
844                 Float_t *sp = new Float_t[ xdim*zdim+1 ];
845                 memset( sp, 0, sizeof(Float_t)*(xdim*zdim+1) );
846
847                 // make a local map from cluster region
848                 for( Int_t ianode=astart; ianode<=astop; ianode++ )
849                 {
850                         for( Int_t itime=tstart; itime<=tstop; itime++ )
851                         {
852                                 Float_t fadc = fMap->GetSignal( ianode, itime );
853                                 if( fadc > baseline ) fadc -= (Double_t)baseline;
854                                 else fadc = 0.;
855                                 Int_t index = (itime-tstart+1)*zdim+(ianode-astart+1);
856                                 sp[index] = fadc;
857                         } // time loop
858                 } // anode loop
859                 
860                 // search peaks on cluster
861                 const Int_t kNp = 150;
862                 Int_t peakX1[kNp];
863                 Int_t peakZ1[kNp];
864                 Float_t peakAmp1[kNp];
865                 Int_t npeak = SearchPeak( sp, xdim, zdim, peakX1, peakZ1, peakAmp1, fMinPeak );
866
867                 // if multiple peaks, split cluster
868                 if( npeak >= 1 )
869                 {
870                 //      cout << "npeak " << npeak << endl;
871                 //      clusterJ->PrintInfo();
872                         
873                         Float_t *par = new Float_t[npeak*5+1];
874                         par[0] = (Float_t)npeak;                
875                         
876                         // Initial paramiters in cell dimentions
877                         Int_t k1 = 1;
878                         for( i=0; i<npeak; i++ ) 
879                         {
880                                 par[k1] = peakAmp1[i];
881                                 par[k1+1] = peakX1[i]; // local time pos. [timebin]
882                                 par[k1+2] = peakZ1[i]; // local anode pos. [anodepitch]
883                                 if( electronics == 1 ) 
884                                 par[k1+3] = 2.; // PASCAL
885                                 else if( electronics == 2 ) 
886                                 par[k1+3] = 0.7; // tau [timebin]  OLA 
887                                 par[k1+4] = .4;    // sigma     [anodepich]
888                                 k1+=5;
889                         }                       
890                         Int_t niter;
891                         Float_t chir;                   
892                         NoLinearFit( xdim, zdim, par, sp, &niter, &chir );
893
894                         Float_t peakX[kNp];
895                         Float_t peakZ[kNp];
896                         Float_t sigma[kNp];
897                         Float_t tau[kNp];
898                         Float_t peakAmp[kNp];
899                         Float_t integral[kNp];
900                         
901                         //get integrals => charge for each peak
902                         PeakFunc( xdim, zdim, par, sp, integral );
903                         
904                         k1 = 1;
905                         for( i=0; i<npeak; i++ ) 
906                         {
907                                 peakAmp[i] = par[k1];
908                                 peakX[i] = par[k1+1];
909                                 peakZ[i] = par[k1+2];
910                                 tau[i] = par[k1+3];
911                                 sigma[i] = par[k1+4];
912                                 k1+=5;
913                         }
914                         
915                         // calculate paramiter for new clusters
916                         for( i=0; i<npeak; i++ )
917                         {
918                                 AliITSRawClusterSDD clusterI( *clusterJ );
919                                 Int_t newAnode = peakZ1[i]-1 + astart;
920                                 Int_t newiTime = peakX1[i]-1 + tstart;
921
922                                 Int_t shift = (Int_t)(fTimeCorr/fTimeStep + 0.5);
923                                 if( newiTime > shift && newiTime < (fMaxNofSamples-shift) )  shift = 0;
924                                 Int_t peakpos = fMap->GetHitIndex( newAnode, newiTime+shift );
925                                 clusterI.SetPeakPos( peakpos );
926                                 clusterI.SetPeakAmpl( peakAmp1[i] );
927
928                                 Float_t newAnodef = peakZ[i] - 0.5 + astart;
929                                 Float_t newiTimef = peakX[i] - 1 + tstart;                              
930                                 if( wing == 2 ) newAnodef -= fNofAnodes; 
931                                 Float_t anodePath = (newAnodef - fNofAnodes/2)*anodePitch;
932                                 newiTimef *= fTimeStep;
933                                 if( newiTimef > fTimeCorr ) newiTimef -= fTimeCorr;
934                                 if( electronics == 1 )
935                                 {
936                                     newiTimef *= 0.999438;    // PASCAL
937                                     newiTimef += (6./fDriftSpeed - newiTimef/3000.);
938                                 }
939                                 else if( electronics == 2 )
940                                     newiTimef *= 0.99714;    // OLA
941
942                                 Float_t driftPath = fSddLength - newiTimef * fDriftSpeed;
943                                 Float_t sign = ( wing == 1 ) ? -1. : 1.;
944                                 clusterI.SetX( driftPath*sign * 0.0001 );       
945                                 clusterI.SetZ( anodePath * 0.0001 );
946                                 clusterI.SetAnode( newAnodef );
947                                 clusterI.SetTime( newiTimef );
948                                 clusterI.SetAsigma( sigma[i]*anodePitch );
949                                 clusterI.SetTsigma( tau[i]*fTimeStep );
950                                 clusterI.SetQ( integral[i] );
951                                 
952                         //      clusterI.PrintInfo();
953                                 iTS->AddCluster( 1, &clusterI );
954                         }
955                         fClusters->RemoveAt( j );
956                         delete [] par;
957                 }
958                 else cout <<" --- Peak not found!!!!  minpeak=" << fMinPeak<< 
959                             " cluster peak=" << clusterJ->PeakAmpl() << endl << endl;
960                 
961                 delete [] sp;
962         } // cluster loop
963
964         fClusters->Compress();
965         fMap->ClearMap(); 
966 }
967
968
969 //_____________________________________________________________________________
970 void  AliITSClusterFinderSDD::GroupClusters()
971 {
972   // group clusters
973   Int_t dummy=0;
974   Float_t fTimeStep = fSegmentation->Dpx(dummy);
975
976
977   // get number of clusters for this module
978   Int_t nofClusters = fClusters->GetEntriesFast();
979   nofClusters -= fNclusters;
980
981   AliITSRawClusterSDD *clusterI;
982   AliITSRawClusterSDD *clusterJ;
983
984   Int_t *label = new Int_t [nofClusters];
985   Int_t i,j;
986   for(i=0; i<nofClusters; i++) label[i] = 0;
987   for(i=0; i<nofClusters; i++) { 
988     if(label[i] != 0) continue;
989     for(j=i+1; j<nofClusters; j++) { 
990       if(label[j] != 0) continue;
991       clusterI = (AliITSRawClusterSDD*) fClusters->At(i);
992       clusterJ = (AliITSRawClusterSDD*) fClusters->At(j);
993       // 1.3 good
994       if(clusterI->T() < fTimeStep*60) fDAnode = 4.2;  // TB 3.2  
995       if(clusterI->T() < fTimeStep*10) fDAnode = 1.5;  // TB 1.
996       Bool_t pair = clusterI->Brother(clusterJ,fDAnode,fDTime);
997       if(!pair) continue;
998       //      clusterI->PrintInfo();
999       //      clusterJ->PrintInfo();
1000       clusterI->Add(clusterJ);
1001       label[j] = 1;
1002       fClusters->RemoveAt(j);
1003       j=i; // <- Ernesto
1004     } // J clusters  
1005     label[i] = 1;
1006   } // I clusters
1007   fClusters->Compress();
1008   
1009   delete [] label;
1010   return;
1011
1012 }
1013
1014 //_____________________________________________________________________________
1015
1016 void AliITSClusterFinderSDD::SelectClusters()
1017 {
1018   // get number of clusters for this module
1019   Int_t nofClusters = fClusters->GetEntriesFast();
1020   nofClusters -= fNclusters;
1021
1022   Int_t i;
1023   for(i=0; i<nofClusters; i++) { 
1024     AliITSRawClusterSDD *clusterI = (AliITSRawClusterSDD*) fClusters->At(i);
1025     Int_t rmflg = 0;
1026     Float_t wy = 0.;
1027     if(clusterI->Anodes() != 0.) {
1028       wy = ((Float_t) clusterI->Samples())/clusterI->Anodes();
1029     }
1030     Int_t amp = (Int_t) clusterI->PeakAmpl();
1031     Int_t cha = (Int_t) clusterI->Q();
1032     if(amp < fMinPeak) rmflg = 1;  
1033     if(cha < fMinCharge) rmflg = 1;
1034     if(wy < fMinNCells) rmflg = 1;
1035     //if(wy > fMaxNCells) rmflg = 1;
1036     if(rmflg) fClusters->RemoveAt(i);
1037   } // I clusters
1038   fClusters->Compress();
1039   return;
1040
1041 }
1042
1043 //_____________________________________________________________________________
1044
1045 void AliITSClusterFinderSDD::ResolveClusters()
1046 {
1047 /*
1048 // The function to resolve clusters if the clusters overlapping exists
1049
1050    AliITS *iTS=(AliITS*)gAlice->GetModule("ITS");
1051
1052   // get number of clusters for this module
1053    Int_t nofClusters = fClusters->GetEntriesFast();
1054    nofClusters -= fNclusters;
1055    //   cout<<"Resolve Cl: nofClusters, fNclusters ="<<nofClusters<<","<<fNclusters<<endl;
1056
1057    Int_t fNofMaps = fSegmentation->Npz();
1058    Int_t fNofAnodes = fNofMaps/2;
1059    Int_t dummy=0;
1060    Double_t fTimeStep = fSegmentation->Dpx(dummy);
1061    Double_t fSddLength = fSegmentation->Dx();
1062    Double_t fDriftSpeed = fResponse->DriftSpeed();
1063    Double_t anodePitch = fSegmentation->Dpz(dummy);
1064    Float_t n, baseline;
1065    fResponse->GetNoiseParam(n,baseline);
1066    Float_t dzz_1A = anodePitch * anodePitch / 12;
1067
1068   // fill Map of signals
1069     fMap->FillMap(); 
1070
1071   Int_t j,i,ii,ianode,anode,itime;
1072   Int_t wing,astart,astop,tstart,tstop,nanode;
1073   Double_t fadc,ClusterTime;
1074   Double_t q[400],x[400],z[400]; // digit charges and coordinates
1075
1076   for(j=0; j<nofClusters; j++) { 
1077
1078     AliITSRawClusterSDD *clusterJ = (AliITSRawClusterSDD*) fClusters->At(j);
1079
1080     Int_t ndigits = 0;
1081     astart=clusterJ->Astart();
1082     astop=clusterJ->Astop();
1083     tstart=clusterJ->Tstartf();
1084     tstop=clusterJ->Tstopf();
1085     nanode=clusterJ->Anodes();  // <- Ernesto
1086     wing=(Int_t)clusterJ->W();
1087     if(wing == 2) {
1088         astart += fNofAnodes; 
1089         astop  += fNofAnodes;
1090     } 
1091
1092 //          cout<<"astart,astop,tstart,tstop ="<<astart<<","<<astop<<","<<tstart<<","<<tstop<<endl;
1093
1094             // clear the digit arrays
1095    for(ii=0; ii<400; ii++) { 
1096        q[ii] = 0.; 
1097        x[ii] = 0.;
1098        z[ii] = 0.;
1099    }
1100     
1101    for(ianode=astart; ianode<=astop; ianode++) { 
1102     for(itime=tstart; itime<=tstop; itime++) { 
1103           fadc=fMap->GetSignal(ianode,itime);
1104           if(fadc>baseline) {
1105             fadc-=(Double_t)baseline;
1106               q[ndigits] = fadc*(fTimeStep/160);  // KeV
1107               anode = ianode;
1108               if(wing == 2) anode -= fNofAnodes;
1109               z[ndigits] = (anode + 0.5 - fNofAnodes/2)*anodePitch;
1110               ClusterTime = itime*fTimeStep;
1111               if(ClusterTime > fTimeCorr) ClusterTime -= fTimeCorr;   // ns
1112               x[ndigits] = fSddLength - ClusterTime*fDriftSpeed;
1113               if(wing == 1) x[ndigits] *= (-1);
1114 //              cout<<"ianode,itime,fadc ="<<ianode<<","<<itime<<","<<fadc<<endl;
1115 //            cout<<"wing,anode,ndigits,charge ="<<wing<<","<<anode<<","<<ndigits<<","<<q[ndigits]<<endl;
1116               ndigits++;
1117               continue;
1118           }
1119               fadc=0;
1120               //              cout<<"fadc=0, ndigits ="<<ndigits<<endl;
1121     } // time loop
1122    } // anode loop
1123    //     cout<<"for new cluster ndigits ="<<ndigits<<endl;
1124
1125
1126    // Fit cluster to resolve for two separate ones --------------------
1127
1128    Double_t qq=0., xm=0., zm=0., xx=0., zz=0., xz=0.;
1129    Double_t dxx=0., dzz=0., dxz=0.;
1130    Double_t scl = 0., tmp, tga, elps = -1.;
1131    Double_t xfit[2], zfit[2], qfit[2];
1132    Double_t pitchz = anodePitch*1.e-4;             // cm
1133    Double_t pitchx = fTimeStep*fDriftSpeed*1.e-4;  // cm
1134    Double_t sigma2;
1135    Int_t nfhits;
1136    Int_t nbins = ndigits;
1137    Int_t separate = 0;
1138
1139    // now, all lengths are in microns
1140
1141    for (ii=0; ii<nbins; ii++) {
1142        qq += q[ii];
1143        xm += x[ii]*q[ii];
1144        zm += z[ii]*q[ii];
1145        xx += x[ii]*x[ii]*q[ii];
1146        zz += z[ii]*z[ii]*q[ii];
1147        xz += x[ii]*z[ii]*q[ii];
1148    }
1149
1150    xm /= qq;
1151    zm /= qq;
1152    xx /= qq;
1153    zz /= qq;
1154    xz /= qq;
1155
1156    dxx = xx - xm*xm;
1157    dzz = zz - zm*zm;
1158    dxz = xz - xm*zm;
1159    
1160    // shrink the cluster in the time direction proportionaly to the 
1161    // dxx/dzz, which lineary depends from the drift path
1162  
1163         // new  Ernesto........  
1164    if( nanode == 1 )
1165    {
1166            dzz = dzz_1A; // for one anode cluster dzz = anode**2/12
1167            scl = TMath::Sqrt( 7.2/(-0.57*xm*1.e-3+71.8) );
1168    }
1169    if( nanode == 2 )
1170    {
1171            scl = TMath::Sqrt( (-0.18*xm*1.e-3+21.3)/(-0.57*xm*1.e-3+71.8) );
1172    }
1173    
1174    if( nanode == 3 )       
1175    {
1176            scl = TMath::Sqrt( (-0.5*xm*1.e-3+34.5)/(-0.57*xm*1.e-3+71.8) );
1177    }
1178
1179    if( nanode > 3 )        
1180    {
1181            scl = TMath::Sqrt( (1.3*xm*1.e-3+49.)/(-0.57*xm*1.e-3+71.8) );
1182    }
1183
1184    //   cout<<"1 microns: zm,dzz,xm,dxx,dxz,qq ="<<zm<<","<<dzz<<","<<xm<<","<<dxx<<","<<dxz<<","<<qq<<endl;
1185
1186  //  old Boris.........
1187  //  tmp=29730. - 585.*fabs(xm/1000.); 
1188  //  scl=TMath::Sqrt(tmp/130000.);
1189    
1190    xm *= scl;
1191    xx *= scl*scl;
1192    xz *= scl;
1193
1194    
1195    dxx = xx - xm*xm;
1196 //   dzz = zz - zm*zm;
1197    dxz = xz - xm*zm;
1198
1199    //   cout<<"microns: zm,dzz,xm,dxx,xz,dxz,qq ="<<zm<<","<<dzz<<","<<xm<<","<<dxx<<","<<xz<<","<<dxz<<","<<qq<<endl;
1200
1201 //   if(dzz < 7200.) dzz = 7200.; // for one anode cluster dzz = anode**2/12
1202   
1203    if (dxx < 0.) dxx=0.;
1204
1205    // the data if no cluster overlapping (the coordunates are in cm) 
1206    nfhits = 1;
1207    xfit[0] = xm*1.e-4;
1208    zfit[0] = zm*1.e-4;
1209    qfit[0] = qq;
1210
1211 //   if(nbins < 7) cout<<"**** nbins ="<<nbins<<endl;
1212   
1213    if (nbins >= 7) {
1214       if (dxz==0.) tga=0.;
1215       else {
1216          tmp=0.5*(dzz-dxx)/dxz;
1217          tga = (dxz<0.) ? tmp-TMath::Sqrt(tmp*tmp+1) : tmp+TMath::Sqrt(tmp*tmp+1);
1218       }
1219       elps=(tga*tga*dxx-2*tga*dxz+dzz)/(dxx+2*tga*dxz+tga*tga*dzz);
1220
1221       // change from microns to cm
1222       xm *= 1.e-4; 
1223       zm *= 1.e-4; 
1224       zz *= 1.e-8;
1225       xx *= 1.e-8;
1226       xz *= 1.e-8;
1227       dxz *= 1.e-8;
1228       dxx *= 1.e-8;
1229       dzz *= 1.e-8;
1230
1231    //   cout<<"cm: zm,dzz,xm,dxx,xz,dxz,qq ="<<zm<<","<<dzz<<","<<xm<<","<<dxx<<","<<xz<<","<<dxz<<","<<qq<<endl;
1232
1233      for (i=0; i<nbins; i++) {     
1234        x[i] = x[i] *= scl;
1235        x[i] = x[i] *= 1.e-4;
1236        z[i] = z[i] *= 1.e-4;
1237      }
1238
1239      //     cout<<"!!! elps ="<<elps<<endl;
1240
1241      if (elps < 0.3) { // try to separate hits 
1242          separate = 1;
1243          tmp=atan(tga);
1244          Double_t cosa=cos(tmp),sina=sin(tmp);
1245          Double_t a1=0., x1=0., xxx=0.;
1246          for (i=0; i<nbins; i++) {
1247              tmp=x[i]*cosa + z[i]*sina;
1248              if (q[i] > a1) {
1249                 a1=q[i];
1250                 x1=tmp;
1251              }
1252              xxx += tmp*tmp*tmp*q[i];
1253          }
1254          xxx /= qq;
1255          Double_t z12=-sina*xm + cosa*zm;
1256          sigma2=(sina*sina*xx-2*cosa*sina*xz+cosa*cosa*zz) - z12*z12;
1257          xm=cosa*xm + sina*zm;
1258          xx=cosa*cosa*xx + 2*cosa*sina*xz + sina*sina*zz;
1259          Double_t x2=(xx - xm*x1 - sigma2)/(xm - x1);
1260          Double_t r=a1*2*TMath::ACos(-1.)*sigma2/(qq*pitchx*pitchz);
1261          for (i=0; i<33; i++) { // solve a system of equations
1262              Double_t x1_old=x1, x2_old=x2, r_old=r;
1263              Double_t c11=x1-x2;
1264              Double_t c12=r;
1265              Double_t c13=1-r;
1266              Double_t c21=x1*x1 - x2*x2;
1267              Double_t c22=2*r*x1;
1268              Double_t c23=2*(1-r)*x2;
1269              Double_t c31=3*sigma2*(x1-x2) + x1*x1*x1 - x2*x2*x2;
1270              Double_t c32=3*r*(sigma2 + x1*x1);
1271              Double_t c33=3*(1-r)*(sigma2 + x2*x2);
1272              Double_t f1=-(r*x1 + (1-r)*x2 - xm);
1273              Double_t f2=-(r*(sigma2 + x1*x1) + (1-r)*(sigma2 + x2*x2) - xx);
1274              Double_t f3=-(r*x1*(3*sigma2+x1*x1) + (1-r)*x2*(3*sigma2+x2*x2)-xxx);
1275              Double_t d=c11*c22*c33 + c21*c32*c13 + c12*c23*c31 - c31*c22*c13 - c21*c12*c33 - c32*c23*c11;
1276              if (d==0.) {
1277                 cout<<"*********** d=0 ***********\n";
1278                 break;
1279              }
1280              Double_t dr=f1*c22*c33 + f2*c32*c13 + c12*c23*f3 -
1281                        f3*c22*c13 - f2*c12*c33 - c32*c23*f1;
1282              Double_t d1=c11*f2*c33 + c21*f3*c13 + f1*c23*c31 -
1283                        c31*f2*c13 - c21*f1*c33 - f3*c23*c11;
1284              Double_t d2=c11*c22*f3 + c21*c32*f1 + c12*f2*c31 -
1285                        c31*c22*f1 - c21*c12*f3 - c32*f2*c11;
1286              r  += dr/d;
1287              x1 += d1/d;
1288              x2 += d2/d;
1289
1290              if (fabs(x1-x1_old) > 0.0001) continue;
1291              if (fabs(x2-x2_old) > 0.0001) continue;
1292              if (fabs(r-r_old)/5 > 0.001) continue;
1293
1294              a1=r*qq*pitchx*pitchz/(2*TMath::ACos(-1.)*sigma2);
1295              Double_t a2=a1*(1-r)/r;
1296              qfit[0]=a1; xfit[0]=x1*cosa - z12*sina; zfit[0]=x1*sina + z12*cosa;
1297              qfit[1]=a2; xfit[1]=x2*cosa - z12*sina; zfit[1]=x2*sina + z12*cosa;
1298              nfhits=2;
1299              break; // Ok !
1300          }
1301          if (i==33) cerr<<"No more iterations ! "<<endl;
1302     } // end of attempt to separate overlapped clusters
1303    } // end of nbins cut 
1304
1305    if(elps < 0.) cout<<" elps=-1 ="<<elps<<endl;
1306    if(elps >0. && elps< 0.3 && nfhits == 1) cout<<" small elps, nfh=1 ="<<elps<<","<<nfhits<<endl;
1307    if(nfhits == 2) cout<<" nfhits=2 ="<<nfhits<<endl;
1308
1309    for (i=0; i<nfhits; i++) {
1310        xfit[i] *= (1.e+4/scl);
1311        if(wing == 1) xfit[i] *= (-1);
1312        zfit[i] *= 1.e+4;
1313        //       cout<<" ---------  i,xfiti,zfiti,qfiti ="<<i<<","<<xfit[i]<<","<<zfit[i]<<","<<qfit[i]<<endl;
1314    }
1315
1316     Int_t ncl = nfhits;
1317     if(nfhits == 1 && separate == 1) {
1318       cout<<"!!!!! no separate"<<endl;
1319       ncl = -2;
1320     } 
1321
1322    if(nfhits == 2) {
1323      cout << "Split cluster: " << endl;
1324      clusterJ->PrintInfo();
1325      cout << " in: " << endl;
1326      for (i=0; i<nfhits; i++) {
1327
1328         // AliITSRawClusterSDD *clust = new AliITSRawClusterSDD(wing,-1,-1,(Float_t)qfit[i],ncl,0,0,(Float_t)xfit[i],(Float_t)zfit[i],0,0,0,0,tstart,tstop,astart,astop);
1329         //      AliITSRawClusterSDD *clust = new AliITSRawClusterSDD(wing,-1,-1,(Float_t)qfit[i],0,0,0,(Float_t)xfit[i],(Float_t)zfit[i],0,0,0,0,tstart,tstop,astart,astop,ncl);
1330
1331         // ???????????
1332        // if(wing == 1) xfit[i] *= (-1);
1333               Float_t Anode = (zfit[i]/anodePitch+fNofAnodes/2-0.5);
1334               Float_t Time = (fSddLength - xfit[i])/fDriftSpeed;
1335         Float_t clusterPeakAmplitude = clusterJ->PeakAmpl();
1336         Float_t peakpos = clusterJ->PeakPos();
1337
1338         Float_t clusteranodePath = (Anode - fNofAnodes/2)*anodePitch;
1339         Float_t clusterDriftPath = Time*fDriftSpeed;
1340         clusterDriftPath = fSddLength-clusterDriftPath;
1341
1342         AliITSRawClusterSDD *clust = new AliITSRawClusterSDD(wing,Anode,Time,qfit[i],
1343                                     clusterPeakAmplitude,peakpos,0.,0.,clusterDriftPath,clusteranodePath,clusterJ->Samples()/2
1344                                     ,tstart,tstop,0,0,0,astart,astop);
1345         clust->PrintInfo();
1346         iTS->AddCluster(1,clust);
1347         //      cout<<"new cluster added: tstart,tstop,astart,astop,x,ncl ="<<tstart<<","<<tstop<<","<<astart<<","<<astop<<","<<xfit[i]<<","<<ncl<<endl;
1348         delete clust;
1349      }// nfhits loop
1350      fClusters->RemoveAt(j);
1351
1352    } // if nfhits = 2
1353   } // cluster loop
1354
1355   fClusters->Compress();
1356   fMap->ClearMap(); 
1357 */          
1358   return;
1359 }
1360
1361
1362 //_____________________________________________________________________________
1363
1364 void AliITSClusterFinderSDD::GetRecPoints()
1365 {
1366   // get rec points
1367
1368   AliITS *iTS=(AliITS*)gAlice->GetModule("ITS");
1369
1370   // get number of clusters for this module
1371   Int_t nofClusters = fClusters->GetEntriesFast();
1372   nofClusters -= fNclusters;
1373
1374   const Float_t kconvGeV = 1.e-6; // GeV -> KeV
1375   const Float_t kconv = 1.0e-4; 
1376   const Float_t kRMSx = 38.0*kconv; // microns->cm ITS TDR Table 1.3
1377   const Float_t kRMSz = 28.0*kconv; // microns->cm ITS TDR Table 1.3
1378
1379
1380   Int_t i;
1381   Int_t ix, iz, idx=-1;
1382   AliITSdigitSDD *dig=0;
1383   Int_t ndigits=fDigits->GetEntriesFast();
1384   for(i=0; i<nofClusters; i++) { 
1385     AliITSRawClusterSDD *clusterI = (AliITSRawClusterSDD*)fClusters->At(i);
1386     if(!clusterI) Error("SDD: GetRecPoints","i clusterI ",i,clusterI);
1387     if(clusterI) idx=clusterI->PeakPos();
1388     if(idx>ndigits) Error("SDD: GetRecPoints","idx ndigits",idx,ndigits);
1389     // try peak neighbours - to be done 
1390     if(idx && idx <= ndigits) dig = (AliITSdigitSDD*)fDigits->UncheckedAt(idx);
1391     if(!dig) {
1392         // try cog
1393         fSegmentation->GetPadIxz(clusterI->X(),clusterI->Z(),ix,iz);
1394         dig = (AliITSdigitSDD*)fMap->GetHit(iz-1,ix-1);
1395         // if null try neighbours
1396         if (!dig) dig = (AliITSdigitSDD*)fMap->GetHit(iz-1,ix); 
1397         if (!dig) dig = (AliITSdigitSDD*)fMap->GetHit(iz-1,ix+1); 
1398         if (!dig) printf("SDD: cannot assign the track number!\n");
1399     }
1400
1401     AliITSRecPoint rnew;
1402     rnew.SetX(clusterI->X());
1403     rnew.SetZ(clusterI->Z());
1404     rnew.SetQ(clusterI->Q());   // in KeV - should be ADC
1405     rnew.SetdEdX(kconvGeV*clusterI->Q());
1406     rnew.SetSigmaX2(kRMSx*kRMSx);
1407     rnew.SetSigmaZ2(kRMSz*kRMSz);
1408     if(dig) rnew.fTracks[0]=dig->fTracks[0];
1409     if(dig) rnew.fTracks[1]=dig->fTracks[1];
1410     if(dig) rnew.fTracks[2]=dig->fTracks[2];
1411     //printf("SDD: i %d track1 track2 track3 %d %d %d x y %f %f\n",i,rnew.fTracks[0],rnew.fTracks[1],rnew.fTracks[2],clusterI->X(),clusterI->Z());
1412     iTS->AddRecPoint(rnew);
1413   } // I clusters
1414
1415   fMap->ClearMap();
1416 }
1417
1418 //_____________________________________________________________________________
1419
1420 void AliITSClusterFinderSDD::FindRawClusters(Int_t mod)
1421 {
1422   // find raw clusters
1423     Find1DClustersE();
1424     GroupClusters();
1425     SelectClusters();
1426     ResolveClustersE();
1427     GetRecPoints();
1428 }
1429 //_____________________________________________________________________________
1430
1431 void AliITSClusterFinderSDD::Print()
1432 {
1433   // Print SDD cluster finder Parameters
1434
1435    cout << "**************************************************" << endl;
1436    cout << " Silicon Drift Detector Cluster Finder Parameters " << endl;
1437    cout << "**************************************************" << endl;
1438    cout << "Number of Clusters: " << fNclusters << endl;
1439    cout << "Anode Tolerance: " << fDAnode << endl;
1440    cout << "Time  Tolerance: " << fDTime << endl;
1441    cout << "Time  correction (electronics): " << fTimeCorr << endl;
1442    cout << "Cut Amplitude (threshold): " << fCutAmplitude << endl;
1443    cout << "Minimum Amplitude: " << fMinPeak << endl;
1444    cout << "Minimum Charge: " << fMinCharge << endl;
1445    cout << "Minimum number of cells/clusters: " << fMinNCells << endl;
1446    cout << "Maximum number of cells/clusters: " << fMaxNCells << endl;
1447    cout << "**************************************************" << endl;
1448 }