]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSclustererV2.cxx
New version of SPD raw-data reconstruction. The format now correponds to the actual...
[u/mrichter/AliRoot.git] / ITS / AliITSclustererV2.cxx
1 //-------------------------------------------------------------------------
2 //            Implementation of the ITS clusterer V2 class
3 //
4 //          Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
5 //-------------------------------------------------------------------------
6 //uncomment the line below for running with V1 cluster finder classes 
7 //#define V1
8
9 #include "AliRun.h"
10
11 #include "AliITSclustererV2.h"
12 #include "AliITSclusterV2.h"
13 #include "AliITSDetTypeRec.h"
14 #include "AliRawReader.h"
15 #include "AliITSRawStreamSPD.h"
16 #include "AliITSRawStreamSDD.h"
17 #include "AliITSRawStreamSSD.h"
18
19 #include <TFile.h>
20 #include <TTree.h>
21 #include <TClonesArray.h>
22 #include "AliITSgeom.h"
23 #include "AliITSdigitSPD.h"
24 #include "AliITSdigitSDD.h"
25 #include "AliITSdigitSSD.h"
26 #include "AliMC.h"
27
28 ClassImp(AliITSclustererV2)
29
30 extern AliRun *gAlice;
31
32 AliITSclustererV2::AliITSclustererV2():
33 fNModules(0),
34 fEvent(0),
35 fI(0),
36 fLastSPD1(0),
37 fNySPD(0),
38 fNzSPD(0),
39 fYpitchSPD(0),
40 fZ1pitchSPD(0),
41 fZ2pitchSPD(0),
42 fHwSPD(0),
43 fHlSPD(0),
44 fNySDD(0),
45 fNzSDD(0),
46 fYpitchSDD(0),
47 fZpitchSDD(0),
48 fHwSDD(0),
49 fHlSDD(0),
50 fYoffSDD(0),
51 fLastSSD1(0),
52 fYpitchSSD(0),
53 fHwSSD(0),
54 fHlSSD(0),
55 fTanP(0),
56 fTanN(0){
57    //default constructor
58  }
59 AliITSclustererV2::AliITSclustererV2(const AliITSgeom *geom):
60 fNModules(0),
61 fEvent(0),
62 fI(0),
63 fLastSPD1(0),
64 fNySPD(256),
65 fNzSPD(160),
66 fYpitchSPD(0.0050),
67 fZ1pitchSPD(0.0425),
68 fZ2pitchSPD(0.0625),
69 fHwSPD(0.64),
70 fHlSPD(3.48),
71 fNySDD(256),
72 fNzSDD(256),
73 fYpitchSDD(0.01825),
74 fZpitchSDD(0.02940),
75 fHwSDD(3.5085),
76 fHlSDD(3.7632),
77 fYoffSDD(0.0425),
78 fLastSSD1(0),
79 fYpitchSSD(0.0095),
80 fHwSSD(3.65),
81 fHlSSD(2.00),
82 fTanP(0.0275),
83 fTanN(0.0075) {
84   //------------------------------------------------------------
85   // Standard constructor
86   //------------------------------------------------------------
87   AliITSgeom *g=(AliITSgeom*)geom;
88
89   Int_t mmax=geom->GetIndexMax();
90   if (mmax>2200) {
91      Fatal("AliITSclustererV2","Too many ITS subdetectors !"); 
92   }
93   Int_t m;
94   for (m=0; m<mmax; m++) {
95      Int_t lay,lad,det; g->GetModuleId(m,lay,lad,det);
96      Float_t x,y,z;     g->GetTrans(lay,lad,det,x,y,z); 
97      Double_t rot[9];   g->GetRotMatrix(lay,lad,det,rot);
98      Double_t alpha=TMath::ATan2(rot[1],rot[0])+TMath::Pi();
99      Double_t ca=TMath::Cos(alpha), sa=TMath::Sin(alpha);
100      fYshift[m] = x*ca + y*sa;
101      fZshift[m] = (Double_t)z;
102      fNdet[m] = (lad-1)*g->GetNdetectors(lay) + (det-1);
103      fNlayer[m] = lay-1;
104   }
105   fNModules = g->GetIndexMax();
106
107   //SPD geometry  
108   fLastSPD1=g->GetModuleIndex(2,1,1)-1;
109   fYSPD[0]=0.5*fYpitchSPD;
110   for (m=1; m<fNySPD; m++) fYSPD[m]=fYSPD[m-1]+fYpitchSPD; 
111   fZSPD[0]=fZ1pitchSPD;
112   for (m=1; m<fNzSPD; m++) {
113     Double_t dz=fZ1pitchSPD;
114     if (m==31 || m==32 || m==63  || m==64  || m==95 || m==96 || 
115         m==127 || m==128) dz=fZ2pitchSPD; 
116     fZSPD[m]=fZSPD[m-1]+dz;
117   }
118   for (m=0; m<fNzSPD; m++) {
119     Double_t dz=0.5*fZ1pitchSPD;
120     if (m==31 || m==32 || m==63  || m==64  || m==95 || m==96 || 
121         m==127 || m==128) dz=0.5*fZ2pitchSPD; 
122     fZSPD[m]-=dz;
123   }
124
125   //SSD geometry
126   fLastSSD1=g->GetModuleIndex(6,1,1)-1;
127
128 }
129
130
131 Int_t AliITSclustererV2::Digits2Clusters(TTree *dTree, TTree *cTree) {
132   //------------------------------------------------------------
133   // This function creates ITS clusters
134   //------------------------------------------------------------
135   Int_t ncl=0;
136
137   if (!dTree) {
138     Error("Digits2Clusters","Can't get the tree with digits !");
139     return 1;
140   }
141
142   TClonesArray *digitsSPD=new TClonesArray("AliITSdigitSPD",3000);
143   dTree->SetBranchAddress("ITSDigitsSPD",&digitsSPD);
144   TClonesArray *digitsSDD=new TClonesArray("AliITSdigitSDD",3000);
145   dTree->SetBranchAddress("ITSDigitsSDD",&digitsSDD);
146   TClonesArray *digitsSSD=new TClonesArray("AliITSdigitSSD",3000);
147   dTree->SetBranchAddress("ITSDigitsSSD",&digitsSSD);
148
149   TClonesArray *clusters=new TClonesArray("AliITSclusterV2",1000);
150   TBranch *branch=cTree->GetBranch("Clusters");
151   if (!branch) cTree->Branch("Clusters",&clusters);
152   else branch->SetAddress(&clusters);
153
154   Int_t mmax=(Int_t)dTree->GetEntries();
155   if (mmax!=fNModules) {
156     Error("Digits2Clusters","Number of entries != number of modules !");
157     return 1;
158   }
159
160   for (fI=0; fI<mmax; fI++) {
161     dTree->GetEvent(fI);
162
163     if     (digitsSPD->GetEntriesFast()!=0) 
164       FindClustersSPD(digitsSPD,clusters);
165     else 
166       if(digitsSDD->GetEntriesFast()!=0) 
167         FindClustersSDD(digitsSDD,clusters);
168       else if(digitsSSD->GetEntriesFast()!=0) 
169         FindClustersSSD(digitsSSD,clusters);
170     
171     ncl+=clusters->GetEntriesFast();
172
173     cTree->Fill();
174
175     digitsSPD->Clear();
176     digitsSDD->Clear();
177     digitsSSD->Clear();
178     clusters->Clear();
179   }
180
181   //cTree->Write();
182
183   delete clusters;
184
185   delete digitsSPD;
186   delete digitsSDD;
187   delete digitsSSD;
188
189   //delete dTree;
190
191   Info("Digits2Clusters","Number of found clusters : %d",ncl);
192
193   return 0;
194 }
195
196 void AliITSclustererV2::Digits2Clusters(AliRawReader* rawReader) {
197   //------------------------------------------------------------
198   // This function creates ITS clusters from raw data
199   //------------------------------------------------------------
200   AliRunLoader* runLoader = AliRunLoader::GetRunLoader();
201   if (!runLoader) {
202     Error("Digits2Clusters", "no run loader found");
203     return;
204   }
205   runLoader->LoadKinematics();
206   AliLoader* itsLoader = runLoader->GetLoader("ITSLoader");
207   if (!itsLoader) {
208     Error("Digits2Clusters", "no loader for ITS found");
209     return;
210   }
211   if (!itsLoader->TreeR()) itsLoader->MakeTree("R");
212   TTree* cTree = itsLoader->TreeR();
213
214   TClonesArray *array=new TClonesArray("AliITSclusterV2",1000);
215   cTree->Branch("Clusters",&array);
216   delete array;
217
218   TClonesArray** clusters = new TClonesArray*[fNModules]; 
219   for (Int_t iModule = 0; iModule < fNModules; iModule++) {
220     clusters[iModule] = NULL;
221   }
222   // one TClonesArray per module
223
224   rawReader->Reset();
225   AliITSRawStreamSPD inputSPD(rawReader);
226   FindClustersSPD(&inputSPD, clusters);
227
228   rawReader->Reset();
229   AliITSRawStreamSDD inputSDD(rawReader);
230   FindClustersSDD(&inputSDD, clusters);
231
232   rawReader->Reset();
233   AliITSRawStreamSSD inputSSD(rawReader);
234   FindClustersSSD(&inputSSD, clusters);
235
236   // write all clusters to the tree
237   Int_t nClusters = 0;
238   for (Int_t iModule = 0; iModule < fNModules; iModule++) {
239     array = clusters[iModule];
240     if (!array) {
241       Error("Digits2Clusters", "data for module %d missing!", iModule);
242       array = new TClonesArray("AliITSclusterV2");
243     }
244     cTree->SetBranchAddress("Clusters", &array);
245     cTree->Fill();
246     nClusters += array->GetEntriesFast();
247     delete array;
248   }
249   itsLoader->WriteRecPoints("OVERWRITE");
250
251   delete[] clusters;
252
253   Info("Digits2Clusters", "total number of found clusters in ITS: %d\n", 
254        nClusters);
255 }
256
257 //**** Fast clusters *******************************
258 #include "TParticle.h"
259
260 //#include "AliITS.h"
261 #include "AliITSmodule.h"
262 #include "AliITSRecPoint.h"
263 #include "AliITSsimulationFastPoints.h"
264 #include "AliITSRecPoint.h"
265
266
267 static void CheckLabels(Int_t lab[3]) {
268   //------------------------------------------------------------
269   // Tries to find mother's labels
270   //------------------------------------------------------------
271     Int_t label=lab[0];
272     if (label>=0) {
273        TParticle *part=(TParticle*)gAlice->GetMCApp()->Particle(label);
274        label=-3;
275        while (part->P() < 0.005) {
276           Int_t m=part->GetFirstMother();
277           if (m<0) {
278              Info("CheckLabels","Primary momentum: %f",part->P()); 
279              break;
280           }
281           if (part->GetStatusCode()>0) {
282              Info("CheckLabels","Primary momentum: %f",part->P()); 
283              break;
284           }
285           label=m;
286           part=(TParticle*)gAlice->GetMCApp()->Particle(label);
287         }
288         if      (lab[1]<0) lab[1]=label;
289         else if (lab[2]<0) lab[2]=label;
290         else ;//cerr<<"CheckLabels : No empty labels !\n";
291     }
292 }
293
294 /*
295 static void CheckLabels(Int_t lab[3]) {
296   //------------------------------------------------------------
297   // Tries to find mother's labels
298   //------------------------------------------------------------
299
300   if(lab[0]<0 && lab[1]<0 && lab[2]<0) return; // In case of no labels just exit
301
302   Int_t ntracks = gAlice->GetMCApp()->GetNtrack();
303   for (Int_t i=0;i<3;i++){
304     Int_t label = lab[i];
305     if (label>=0 && label<ntracks) {
306       TParticle *part=(TParticle*)gAlice->GetMCApp()->Particle(label);
307       if (part->P() < 0.005) {
308         Int_t m=part->GetFirstMother();
309         if (m<0) {      
310           continue;
311         }
312         if (part->GetStatusCode()>0) {
313           continue;
314         }
315         lab[i]=m;       
316       }
317     }    
318   }
319   
320 }
321 */
322 static void CheckLabels2(Int_t lab[10]) {
323   //------------------------------------------------------------
324   // Tries to find mother's labels
325   //------------------------------------------------------------
326   Int_t nlabels =0; 
327   for (Int_t i=0;i<10;i++) if (lab[i]>=0) nlabels++;
328   if(nlabels == 0) return; // In case of no labels just exit
329
330   Int_t ntracks = gAlice->GetMCApp()->GetNtrack();
331
332   for (Int_t i=0;i<10;i++){
333     Int_t label = lab[i];
334     if (label>=0 && label<ntracks) {
335       TParticle *part=(TParticle*)gAlice->GetMCApp()->Particle(label);
336       if (part->P() < 0.02) {
337           Int_t m=part->GetFirstMother();
338           if (m<0) {    
339             continue;
340           }
341           if (part->GetStatusCode()>0) {
342             continue;
343           }
344           lab[i]=m;       
345       }
346       else
347         if (part->P() < 0.12 && nlabels>3) {
348           lab[i]=-2;
349           nlabels--;
350         } 
351     }
352     else{
353       if ( (label>ntracks||label <0) && nlabels>3) {
354         lab[i]=-2;
355         nlabels--;
356       } 
357     }
358   }  
359   if (nlabels>3){
360     for (Int_t i=0;i<10;i++){
361       if (nlabels>3){
362         Int_t label = lab[i];
363         if (label>=0 && label<ntracks) {
364           TParticle *part=(TParticle*)gAlice->GetMCApp()->Particle(label);
365           if (part->P() < 0.1) {
366             lab[i]=-2;
367             nlabels--;
368           }
369         }
370       }
371     }
372   }
373
374   //compress labels -- if multi-times the same
375   Int_t lab2[10];
376   for (Int_t i=0;i<10;i++) lab2[i]=-2;
377   for (Int_t i=0;i<10  ;i++){
378     if (lab[i]<0) continue;
379     for (Int_t j=0;j<10 &&lab2[j]!=lab[i];j++){
380       if (lab2[j]<0) {
381         lab2[j]= lab[i];
382         break;
383       }
384     }
385   }
386   for (Int_t j=0;j<10;j++) lab[j]=lab2[j];
387   
388 }
389
390 static void AddLabel(Int_t lab[10], Int_t label) {
391
392   if(label<0) return; // In case of no label just exit
393
394   Int_t ntracks = gAlice->GetMCApp()->GetNtrack();
395   if (label>ntracks) return;
396   for (Int_t i=0;i<10;i++){
397     //    if (label<0) break;
398     if (lab[i]==label) break;
399     if (lab[i]<0) {
400       lab[i]= label;
401       break;
402     }
403   }
404 }
405
406 void AliITSclustererV2::RecPoints2Clusters
407 (const TClonesArray *points, Int_t idx, TClonesArray *clusters) {
408   //------------------------------------------------------------
409   // Conversion AliITSRecPoint -> AliITSclusterV2 for the ITS 
410   // subdetector indexed by idx 
411   //------------------------------------------------------------
412   TClonesArray &cl=*clusters;
413   Int_t ncl=points->GetEntriesFast();
414   for (Int_t i=0; i<ncl; i++) {
415     AliITSRecPoint *p = (AliITSRecPoint *)points->UncheckedAt(i);
416     Float_t lp[5];
417     lp[0]=-(-p->GetDetLocalX()+fYshift[idx]); if (idx<=fLastSPD1) lp[0]*=-1; //SPD1
418     lp[1]=  -p->GetZ()+fZshift[idx];
419     lp[2]=p->GetSigmaDetLocX2();
420     lp[3]=p->GetSigmaZ2();
421     lp[4]=p->GetQ()*36./23333.;  //electrons -> ADC
422     Int_t lab[4]; 
423     lab[0]=p->GetLabel(0); lab[1]=p->GetLabel(1); lab[2]=p->GetLabel(2);
424     lab[3]=fNdet[idx];
425     CheckLabels(lab);
426     Int_t dummy[3]={0,0,0};
427     new (cl[i]) AliITSclusterV2(lab,lp, dummy);
428   }  
429
430
431 //***********************************
432
433 #ifndef V1
434
435 void AliITSclustererV2:: 
436 FindCluster(Int_t k,Int_t maxz,AliBin *bins,Int_t &n,Int_t *idx) {
437   //------------------------------------------------------------
438   // returns an array of indices of digits belonging to the cluster
439   // (needed when the segmentation is not regular) 
440   //------------------------------------------------------------
441   if (n<200) idx[n++]=bins[k].GetIndex();
442   bins[k].Use();
443
444   if (bins[k-maxz].IsNotUsed()) FindCluster(k-maxz,maxz,bins,n,idx);
445   if (bins[k-1   ].IsNotUsed()) FindCluster(k-1   ,maxz,bins,n,idx);
446   if (bins[k+maxz].IsNotUsed()) FindCluster(k+maxz,maxz,bins,n,idx);
447   if (bins[k+1   ].IsNotUsed()) FindCluster(k+1   ,maxz,bins,n,idx);
448   /*
449   if (bins[k-maxz-1].IsNotUsed()) FindCluster(k-maxz-1,maxz,bins,n,idx);
450   if (bins[k-maxz+1].IsNotUsed()) FindCluster(k-maxz+1,maxz,bins,n,idx);
451   if (bins[k+maxz-1].IsNotUsed()) FindCluster(k+maxz-1,maxz,bins,n,idx);
452   if (bins[k+maxz+1].IsNotUsed()) FindCluster(k+maxz+1,maxz,bins,n,idx);
453   */
454 }
455
456 void AliITSclustererV2::
457 FindClustersSPD(const TClonesArray *digits, TClonesArray *clusters) {
458   //------------------------------------------------------------
459   // Actual SPD cluster finder
460   //------------------------------------------------------------
461   Int_t kNzBins = fNzSPD + 2;
462   const Int_t kMAXBIN=kNzBins*(fNySPD+2);
463
464   Int_t ndigits=digits->GetEntriesFast();
465   AliBin *bins=new AliBin[kMAXBIN];
466
467   Int_t k;
468   AliITSdigitSPD *d=0;
469   for (k=0; k<ndigits; k++) {
470      d=(AliITSdigitSPD*)digits->UncheckedAt(k);
471      Int_t i=d->GetCoord2()+1;   //y
472      Int_t j=d->GetCoord1()+1;
473      bins[i*kNzBins+j].SetIndex(k);
474      bins[i*kNzBins+j].SetMask(1);
475   }
476    
477   Int_t n=0; TClonesArray &cl=*clusters;
478   for (k=0; k<kMAXBIN; k++) {
479      if (!bins[k].IsNotUsed()) continue;
480      Int_t ni=0, idx[200];
481      FindCluster(k,kNzBins,bins,ni,idx);
482      if (ni==200) {
483         Info("FindClustersSPD","Too big cluster !"); 
484         continue;
485      }
486      Int_t milab[10];
487      for (Int_t ilab=0;ilab<10;ilab++){
488        milab[ilab]=-2;
489      }
490
491      d=(AliITSdigitSPD*)digits->UncheckedAt(idx[0]);
492      Int_t ymin=d->GetCoord2(),ymax=ymin;
493      Int_t zmin=d->GetCoord1(),zmax=zmin;
494
495      for (Int_t l=0; l<ni; l++) {
496         d=(AliITSdigitSPD*)digits->UncheckedAt(idx[l]);
497
498         if (ymin > d->GetCoord2()) ymin=d->GetCoord2();
499         if (ymax < d->GetCoord2()) ymax=d->GetCoord2();
500         if (zmin > d->GetCoord1()) zmin=d->GetCoord1();
501         if (zmax < d->GetCoord1()) zmax=d->GetCoord1();
502         // MI addition - find all labels in cluster
503         for (Int_t dlab=0;dlab<10;dlab++){
504           Int_t digitlab = (d->GetTracks())[dlab];
505           if (digitlab<0) continue;
506           AddLabel(milab,digitlab);       
507         }
508         if (milab[9]>0) CheckLabels2(milab);
509      }
510      CheckLabels2(milab);
511      //
512      //Int_t idy = (fNlayer[fI]==0)? 2:3; 
513      //for (Int_t iz=zmin; iz<=zmax;iz+=2)
514      //Int_t idy = (ymax-ymin)/4.; // max 2 clusters
515      Int_t idy = 0; // max 2 clusters
516      if (fNlayer[fI]==0 &&idy<3) idy=3;
517      if (fNlayer[fI]==1 &&idy<4) idy=4; 
518      Int_t idz =3;
519      for (Int_t iz=zmin; iz<=zmax;iz+=idz)
520        for (Int_t iy=ymin; iy<=ymax;iy+=idy){
521          //
522          Int_t ndigits =0;
523          Float_t y=0.,z=0.,q=0.;         
524          for (Int_t l=0; l<ni; l++) {
525            d=(AliITSdigitSPD*)digits->UncheckedAt(idx[l]);
526            if (zmax-zmin>=idz || ymax-ymin>=idy){
527              if (TMath::Abs( d->GetCoord2()-iy)>0.75*idy) continue;
528              if (TMath::Abs( d->GetCoord1()-iz)>0.75*idz) continue;
529            }
530            ndigits++;
531            Float_t qq=d->GetSignal();
532            y+=qq*fYSPD[d->GetCoord2()]; z+=qq*fZSPD[d->GetCoord1()]; q+=qq;   
533           
534          }     
535          if (ndigits==0) continue;
536          y/=q; z/=q;
537          y-=fHwSPD; z-=fHlSPD;
538          
539          Float_t lp[5];
540          lp[0]=-(-y+fYshift[fI]); if (fI<=fLastSPD1) lp[0]=-lp[0];
541          lp[1]=  -z+fZshift[fI];
542          // Float_t factor=TMath::Max(double(ni-3.),1.5);
543          Float_t factory=TMath::Max(ymax-ymin,1);
544          Float_t factorz=TMath::Max(zmax-zmin,1);
545          factory*= factory;
546          factorz*= factorz;     
547          //lp[2]= (fYpitchSPD*fYpitchSPD/12.)*factory;
548          //lp[3]= (fZ1pitchSPD*fZ1pitchSPD/12.)*factorz;
549          lp[2]= (fYpitchSPD*fYpitchSPD/12.);
550          lp[3]= (fZ1pitchSPD*fZ1pitchSPD/12.);
551          //lp[4]= q;
552          lp[4]= (zmax-zmin+1)*100 + (ymax-ymin+1);
553          
554          milab[3]=fNdet[fI];
555          d=(AliITSdigitSPD*)digits->UncheckedAt(idx[0]);
556          Int_t info[3] = {ymax-ymin+1,zmax-zmin+1,fNlayer[fI]};
557          new (cl[n]) AliITSclusterV2(milab,lp,info); n++;        
558        }
559   }
560   
561   delete [] bins;
562 }
563
564 void AliITSclustererV2::FindClustersSPD(AliITSRawStream* input, 
565                                         TClonesArray** clusters) 
566 {
567   //------------------------------------------------------------
568   // Actual SPD cluster finder for raw data
569   //------------------------------------------------------------
570
571   Int_t nClustersSPD = 0;
572   Int_t kNzBins = fNzSPD + 2;
573   Int_t kNyBins = fNySPD + 2;
574   Int_t kMaxBin = kNzBins * kNyBins;
575   AliBin *binsSPD = new AliBin[kMaxBin];
576   AliBin *binsSPDInit = new AliBin[kMaxBin];
577   AliBin *bins = NULL;
578
579   // read raw data input stream
580   while (kTRUE) {
581     Bool_t next = input->Next();
582     if (!next || input->IsNewModule()) {
583       Int_t iModule = input->GetPrevModuleID();
584
585       // when all data from a module was read, search for clusters
586       if (bins) { 
587         clusters[iModule] = new TClonesArray("AliITSclusterV2");
588         Int_t nClusters = 0;
589
590         for (Int_t iBin = 0; iBin < kMaxBin; iBin++) {
591           if (bins[iBin].IsUsed()) continue;
592           Int_t nBins = 0;
593           Int_t idxBins[200];
594           FindCluster(iBin, kNzBins, bins, nBins, idxBins);
595           if (nBins == 200) {
596             Error("FindClustersSPD", "SPD: Too big cluster !\n"); 
597             continue;
598           }
599
600           Int_t label[4]; 
601           label[0] = -2;
602           label[1] = -2;
603           label[2] = -2;
604 //        label[3] = iModule;
605           label[3] = fNdet[iModule];
606
607           Int_t ymin = (idxBins[0] / kNzBins) - 1;
608           Int_t ymax = ymin;
609           Int_t zmin = (idxBins[0] % kNzBins) - 1;
610           Int_t zmax = zmin;
611           for (Int_t idx = 0; idx < nBins; idx++) {
612             Int_t iy = (idxBins[idx] / kNzBins) - 1;
613             Int_t iz = (idxBins[idx] % kNzBins) - 1;
614             if (ymin > iy) ymin = iy;
615             if (ymax < iy) ymax = iy;
616             if (zmin > iz) zmin = iz;
617             if (zmax < iz) zmax = iz;
618           }
619
620           Int_t idy = 0; // max 2 clusters
621           if ((iModule <= fLastSPD1) &&idy<3) idy=3;
622           if ((iModule > fLastSPD1) &&idy<4) idy=4; 
623           Int_t idz =3;
624           for (Int_t iiz=zmin; iiz<=zmax;iiz+=idz)
625             for (Int_t iiy=ymin; iiy<=ymax;iiy+=idy){
626               //
627               Int_t ndigits =0;
628               Float_t y=0.,z=0.,q=0.;    
629               for (Int_t idx = 0; idx < nBins; idx++) {
630                 Int_t iy = (idxBins[idx] / kNzBins) - 1;
631                 Int_t iz = (idxBins[idx] % kNzBins) - 1;
632                 if (zmax-zmin>=idz || ymax-ymin>=idy){
633                   if (TMath::Abs(iy-iiy)>0.75*idy) continue;
634                   if (TMath::Abs(iz-iiz)>0.75*idz) continue;
635                 }
636                 ndigits++;
637                 Float_t qBin = bins[idxBins[idx]].GetQ();
638                 y += qBin * fYSPD[iy]; 
639                 z += qBin * fZSPD[iz]; 
640                 q += qBin;   
641               }
642               if (ndigits==0) continue;
643               y /= q; 
644               z /= q;
645               y -= fHwSPD; 
646               z -= fHlSPD;
647
648               Float_t hit[5];  // y, z, sigma(y)^2, sigma(z)^2, charge
649               hit[0] = -(-y+fYshift[iModule]); 
650               if (iModule <= fLastSPD1) hit[0] = -hit[0];
651               hit[1] = -z+fZshift[iModule];
652               hit[2] = fYpitchSPD*fYpitchSPD/12.;
653               hit[3] = fZ1pitchSPD*fZ1pitchSPD/12.;
654               //          hit[4] = q;
655               hit[4] = (zmax-zmin+1)*100 + (ymax-ymin+1);
656               //          CheckLabels(label);
657               Int_t info[3]={ymax-ymin+1,zmax-zmin+1,fNlayer[iModule]};
658               new (clusters[iModule]->AddrAt(nClusters)) 
659                 AliITSclusterV2(label, hit,info); 
660               nClusters++;
661             }
662         }
663
664         nClustersSPD += nClusters;
665         bins = NULL;
666       }
667
668       if (!next) break;
669       bins = binsSPD;
670       memcpy(binsSPD,binsSPDInit,sizeof(AliBin)*kMaxBin);
671     }
672
673     // fill the current digit into the bins array
674     Int_t index = (input->GetCoord2()+1) * kNzBins + (input->GetCoord1()+1);
675     bins[index].SetIndex(index);
676     bins[index].SetMask(1);
677     bins[index].SetQ(1);
678   }
679
680   delete [] binsSPDInit;
681   delete [] binsSPD;
682
683   Info("FindClustersSPD", "found clusters in ITS SPD: %d", nClustersSPD);
684 }
685
686
687 Bool_t AliITSclustererV2::IsMaximum(Int_t k,Int_t max,const AliBin *bins) {
688   //------------------------------------------------------------
689   //is this a local maximum ?
690   //------------------------------------------------------------
691   UShort_t q=bins[k].GetQ();
692   if (q==1023) return kFALSE;
693   if (bins[k-max].GetQ() > q) return kFALSE;
694   if (bins[k-1  ].GetQ() > q) return kFALSE; 
695   if (bins[k+max].GetQ() > q) return kFALSE; 
696   if (bins[k+1  ].GetQ() > q) return kFALSE; 
697   if (bins[k-max-1].GetQ() > q) return kFALSE;
698   if (bins[k+max-1].GetQ() > q) return kFALSE; 
699   if (bins[k+max+1].GetQ() > q) return kFALSE; 
700   if (bins[k-max+1].GetQ() > q) return kFALSE;
701   return kTRUE; 
702 }
703
704 void AliITSclustererV2::
705 FindPeaks(Int_t k,Int_t max,AliBin *b,Int_t *idx,UInt_t *msk,Int_t& n) {
706   //------------------------------------------------------------
707   //find local maxima
708   //------------------------------------------------------------
709   if (n<31)
710   if (IsMaximum(k,max,b)) {
711     idx[n]=k; msk[n]=(2<<n);
712     n++;
713   }
714   b[k].SetMask(0);
715   if (b[k-max].GetMask()&1) FindPeaks(k-max,max,b,idx,msk,n);
716   if (b[k-1  ].GetMask()&1) FindPeaks(k-1  ,max,b,idx,msk,n);
717   if (b[k+max].GetMask()&1) FindPeaks(k+max,max,b,idx,msk,n);
718   if (b[k+1  ].GetMask()&1) FindPeaks(k+1  ,max,b,idx,msk,n);
719 }
720
721 void AliITSclustererV2::
722 MarkPeak(Int_t k, Int_t max, AliBin *bins, UInt_t m) {
723   //------------------------------------------------------------
724   //mark this peak
725   //------------------------------------------------------------
726   UShort_t q=bins[k].GetQ();
727
728   bins[k].SetMask(bins[k].GetMask()|m); 
729
730   if (bins[k-max].GetQ() <= q)
731      if ((bins[k-max].GetMask()&m) == 0) MarkPeak(k-max,max,bins,m);
732   if (bins[k-1  ].GetQ() <= q)
733      if ((bins[k-1  ].GetMask()&m) == 0) MarkPeak(k-1  ,max,bins,m);
734   if (bins[k+max].GetQ() <= q)
735      if ((bins[k+max].GetMask()&m) == 0) MarkPeak(k+max,max,bins,m);
736   if (bins[k+1  ].GetQ() <= q)
737      if ((bins[k+1  ].GetMask()&m) == 0) MarkPeak(k+1  ,max,bins,m);
738 }
739
740 void AliITSclustererV2::
741 MakeCluster(Int_t k,Int_t max,AliBin *bins,UInt_t m,AliITSclusterV2 &c) {
742   //------------------------------------------------------------
743   //make cluster using digits of this peak
744   //------------------------------------------------------------
745   Float_t q=(Float_t)bins[k].GetQ();
746   Int_t i=k/max, j=k-i*max;
747
748   c.SetQ(c.GetQ()+q);
749   c.SetY(c.GetY()+i*q); 
750   c.SetZ(c.GetZ()+j*q); 
751   c.SetSigmaY2(c.GetSigmaY2()+i*i*q);
752   c.SetSigmaZ2(c.GetSigmaZ2()+j*j*q);
753
754   bins[k].SetMask(0xFFFFFFFE);
755   
756   if (bins[k-max].GetMask() == m) MakeCluster(k-max,max,bins,m,c);
757   if (bins[k-1  ].GetMask() == m) MakeCluster(k-1  ,max,bins,m,c);
758   if (bins[k+max].GetMask() == m) MakeCluster(k+max,max,bins,m,c);
759   if (bins[k+1  ].GetMask() == m) MakeCluster(k+1  ,max,bins,m,c);
760 }
761
762 void AliITSclustererV2::
763 FindClustersSDD(AliBin* bins[2], Int_t nMaxBin, Int_t nzBins, 
764                 const TClonesArray *digits, TClonesArray *clusters) {
765   //------------------------------------------------------------
766   // Actual SDD cluster finder
767   //------------------------------------------------------------
768   Int_t ncl=0; TClonesArray &cl=*clusters;
769   for (Int_t s=0; s<2; s++)
770     for (Int_t i=0; i<nMaxBin; i++) {
771       if (bins[s][i].IsUsed()) continue;
772       Int_t idx[32]; UInt_t msk[32]; Int_t npeaks=0;
773       FindPeaks(i, nzBins, bins[s], idx, msk, npeaks);
774
775       if (npeaks>30) continue;
776       if (npeaks==0) continue;
777
778       Int_t k,l;
779       for (k=0; k<npeaks-1; k++){//mark adjacent peaks
780         if (idx[k] < 0) continue; //this peak is already removed
781         for (l=k+1; l<npeaks; l++) {
782            if (idx[l] < 0) continue; //this peak is already removed
783            Int_t ki=idx[k]/nzBins, kj=idx[k] - ki*nzBins;
784            Int_t li=idx[l]/nzBins, lj=idx[l] - li*nzBins;
785            Int_t di=TMath::Abs(ki - li);
786            Int_t dj=TMath::Abs(kj - lj);
787            if (di>1 || dj>1) continue;
788            if (bins[s][idx[k]].GetQ() > bins[s][idx[l]].GetQ()) {
789               msk[l]=msk[k];
790               idx[l]*=-1;
791            } else {
792               msk[k]=msk[l];
793               idx[k]*=-1;
794               break;
795            } 
796         }
797       }
798
799       for (k=0; k<npeaks; k++) {
800         MarkPeak(TMath::Abs(idx[k]), nzBins, bins[s], msk[k]);
801       }
802         
803       for (k=0; k<npeaks; k++) {
804          if (idx[k] < 0) continue; //removed peak
805          AliITSclusterV2 c;
806          MakeCluster(idx[k], nzBins, bins[s], msk[k], c);
807          //mi change
808          Int_t milab[10];
809          for (Int_t ilab=0;ilab<10;ilab++){
810            milab[ilab]=-2;
811          }
812          Int_t maxi=0,mini=0,maxj=0,minj=0;
813          //AliBin *bmax=&bins[s][idx[k]];
814          //Float_t max = TMath::Max(TMath::Abs(bmax->GetQ())/5.,3.);
815          Float_t max=3;
816          for (Int_t di=-2; di<=2;di++)
817            for (Int_t dj=-3;dj<=3;dj++){
818              Int_t index = idx[k]+di+dj*nzBins;
819              if (index<0) continue;
820              if (index>=nMaxBin) continue;
821              AliBin *b=&bins[s][index];
822              if (TMath::Abs(b->GetQ())>max){
823                if (di>maxi) maxi=di;
824                if (di<mini) mini=di;
825                if (dj>maxj) maxj=dj;
826                if (dj<minj) minj=dj;
827                //
828                if(digits) {
829                  if (TMath::Abs(di)<2&&TMath::Abs(dj)<2){
830                    AliITSdigitSDD* d=(AliITSdigitSDD*)digits->UncheckedAt(b->GetIndex());
831                    for (Int_t itrack=0;itrack<10;itrack++){
832                      Int_t track = (d->GetTracks())[itrack];
833                      if (track>=0) {
834                        AddLabel(milab, track); 
835                      }
836                    }
837                  }
838                }
839              }
840            }
841          
842          /* 
843             Float_t s2 = c.GetSigmaY2()/c.GetQ() - c.GetY()*c.GetY();
844             Float_t w=par->GetPadPitchWidth(sec);
845             c.SetSigmaY2(s2);
846             if (s2 != 0.) {
847             c.SetSigmaY2(c.GetSigmaY2()*0.108);
848             if (sec<par->GetNInnerSector()) c.SetSigmaY2(c.GetSigmaY2()*2.07);
849             }    
850             s2 = c.GetSigmaZ2()/c.GetQ() - c.GetZ()*c.GetZ();
851             w=par->GetZWidth();
852             c.SetSigmaZ2(s2);
853             
854             if (s2 != 0.) {
855             c.SetSigmaZ2(c.GetSigmaZ2()*0.169);
856             if (sec<par->GetNInnerSector()) c.SetSigmaZ2(c.GetSigmaZ2()*1.77);
857             }
858          */
859
860          c.SetSigmaY2(0.0030*0.0030);
861          c.SetSigmaZ2(0.0020*0.0020);
862          c.SetDetectorIndex(fNdet[fI]);
863
864          Float_t y=c.GetY(),z=c.GetZ(), q=c.GetQ();
865          y/=q; z/=q;
866          //
867          //Float_t s2 = c.GetSigmaY2()/c.GetQ() - y*y;
868          // c.SetSigmaY2(s2);
869          //s2 = c.GetSigmaZ2()/c.GetQ() - z*z;
870          //c.SetSigmaZ2(s2);
871          //
872          y=(y-0.5)*fYpitchSDD;
873          y-=fHwSDD;
874          y-=fYoffSDD;  //delay ?
875          if (s) y=-y;
876
877          z=(z-0.5)*fZpitchSDD;
878          z-=fHlSDD;
879
880          y=-(-y+fYshift[fI]);
881          z=  -z+fZshift[fI];
882          c.SetY(y);
883          c.SetZ(z);
884          c.SetNy(maxj-minj+1);
885          c.SetNz(maxi-mini+1);
886          c.SetType(npeaks);
887          c.SetQ(q/12.7);  //to be consistent with the SSD charges
888
889          if (c.GetQ() < 20.) continue; //noise cluster
890          
891          if (digits) {    
892            //      AliBin *b=&bins[s][idx[k]];
893            //      AliITSdigitSDD* d=(AliITSdigitSDD*)digits->UncheckedAt(b->GetIndex());
894            {
895              //Int_t lab[3];
896              //lab[0]=(d->GetTracks())[0];
897              //lab[1]=(d->GetTracks())[1];
898              //lab[2]=(d->GetTracks())[2];
899              //CheckLabels(lab);
900              CheckLabels2(milab); 
901              c.SetLabel(milab[0],0);
902              c.SetLabel(milab[1],1);
903              c.SetLabel(milab[2],2);
904              c.SetLayer(fNlayer[fI]);
905            }
906          }
907          new (cl[ncl]) AliITSclusterV2(c); ncl++;
908       }
909     }
910 }
911
912 void AliITSclustererV2::
913 FindClustersSDD(const TClonesArray *digits, TClonesArray *clusters) {
914   //------------------------------------------------------------
915   // Actual SDD cluster finder
916   //------------------------------------------------------------
917   Int_t kNzBins = fNzSDD + 2;
918   const Int_t kMAXBIN=kNzBins*(fNySDD+2);
919
920   AliBin *bins[2];
921   bins[0]=new AliBin[kMAXBIN];
922   bins[1]=new AliBin[kMAXBIN];
923
924   AliITSdigitSDD *d=0;
925   Int_t i, ndigits=digits->GetEntriesFast();
926   for (i=0; i<ndigits; i++) {
927      d=(AliITSdigitSDD*)digits->UncheckedAt(i);
928      Int_t y=d->GetCoord2()+1;   //y
929      Int_t z=d->GetCoord1()+1;   //z
930      Int_t q=d->GetSignal();
931      if (q<3) continue;
932
933      if (z <= fNzSDD) {
934        bins[0][y*kNzBins+z].SetQ(q);
935        bins[0][y*kNzBins+z].SetMask(1);
936        bins[0][y*kNzBins+z].SetIndex(i);
937      } else {
938        z-=fNzSDD; 
939        bins[1][y*kNzBins+z].SetQ(q);
940        bins[1][y*kNzBins+z].SetMask(1);
941        bins[1][y*kNzBins+z].SetIndex(i);
942      }
943   }
944   
945   FindClustersSDD(bins, kMAXBIN, kNzBins, digits, clusters);
946
947   delete[] bins[0];
948   delete[] bins[1];
949
950 }
951
952 void AliITSclustererV2::FindClustersSDD(AliITSRawStream* input, 
953                                         TClonesArray** clusters) 
954 {
955   //------------------------------------------------------------
956   // Actual SDD cluster finder for raw data
957   //------------------------------------------------------------
958   Int_t nClustersSDD = 0;
959   Int_t kNzBins = fNzSDD + 2;
960   Int_t kMaxBin = kNzBins * (fNySDD+2);
961   AliBin *binsSDDInit = new AliBin[kMaxBin];
962   AliBin *binsSDD1 = new AliBin[kMaxBin];
963   AliBin *binsSDD2 = new AliBin[kMaxBin];
964   AliBin *bins[2] = {NULL, NULL};
965
966   // read raw data input stream
967   while (kTRUE) {
968     Bool_t next = input->Next();
969     if (!next || input->IsNewModule()) {
970       Int_t iModule = input->GetPrevModuleID();
971
972       // when all data from a module was read, search for clusters
973       if (bins[0]) { 
974         clusters[iModule] = new TClonesArray("AliITSclusterV2");
975         fI = iModule;
976         FindClustersSDD(bins, kMaxBin, kNzBins, NULL, clusters[iModule]);
977         Int_t nClusters = clusters[iModule]->GetEntriesFast();
978         nClustersSDD += nClusters;
979         bins[0] = bins[1] = NULL;
980       }
981
982       if (!next) break;
983       bins[0]=binsSDD1;
984       bins[1]=binsSDD2;
985       memcpy(binsSDD1,binsSDDInit,sizeof(AliBin)*kMaxBin);
986       memcpy(binsSDD2,binsSDDInit,sizeof(AliBin)*kMaxBin);
987     }
988
989     // fill the current digit into the bins array
990     if(input->GetSignal()>=3) {
991       Int_t iz = input->GetCoord1()+1;
992       Int_t side = ((iz <= fNzSDD) ? 0 : 1);
993       iz -= side*fNzSDD;
994       Int_t index = (input->GetCoord2()+1) * kNzBins + iz;
995       bins[side][index].SetQ(input->GetSignal());
996       bins[side][index].SetMask(1);
997       bins[side][index].SetIndex(index);
998     }
999   }
1000
1001   delete[] binsSDD1;
1002   delete[] binsSDD2;
1003   delete[] binsSDDInit;
1004
1005   Info("FindClustersSDD", "found clusters in ITS SDD: %d", nClustersSDD);
1006 }
1007
1008
1009
1010 void AliITSclustererV2::
1011 FindClustersSSD(Ali1Dcluster* neg, Int_t nn, 
1012                 Ali1Dcluster* pos, Int_t np,
1013                 TClonesArray *clusters) {
1014   //------------------------------------------------------------
1015   // Actual SSD cluster finder
1016   //------------------------------------------------------------
1017   TClonesArray &cl=*clusters;
1018   //
1019   Float_t tanp=fTanP, tann=fTanN;
1020   if (fI>fLastSSD1) {tann=fTanP; tanp=fTanN;}
1021   Int_t idet=fNdet[fI];
1022   Int_t ncl=0;
1023   //
1024   Int_t negativepair[30000];
1025   Int_t cnegative[3000];  
1026   Int_t cused1[3000];
1027   Int_t positivepair[30000];
1028   Int_t cpositive[3000];
1029   Int_t cused2[3000];
1030   for (Int_t i=0;i<3000;i++) {cnegative[i]=0; cused1[i]=0;}
1031   for (Int_t i=0;i<3000;i++) {cpositive[i]=0; cused2[i]=0;}
1032   for (Int_t i=0;i<30000;i++){negativepair[i]=positivepair[i]=0;}
1033   static Short_t pairs[1000][1000];
1034   memset(pairs,0,sizeof(Short_t)*1000000);
1035 //   Short_t ** pairs = new Short_t*[1000];
1036 //   for (Int_t i=0; i<1000; i++) {
1037 //     pairs[i] = new Short_t[1000];
1038 //     memset(pairs[i],0,sizeof(Short_t)*1000);
1039 //   }  
1040   //
1041   // find available pairs
1042   //
1043   for (Int_t i=0; i<np; i++) {
1044     Float_t yp=pos[i].GetY()*fYpitchSSD; 
1045     if (pos[i].GetQ()<3) continue;
1046     for (Int_t j=0; j<nn; j++) {
1047       if (neg[j].GetQ()<3) continue;
1048       Float_t yn=neg[j].GetY()*fYpitchSSD;
1049       Float_t zt=(2*fHlSSD*tanp + yp - yn)/(tann+tanp);
1050       Float_t yt=yn + tann*zt;
1051       zt-=fHlSSD; yt-=fHwSSD;
1052       if (TMath::Abs(yt)<fHwSSD+0.01)
1053       if (TMath::Abs(zt)<fHlSSD+0.01*(neg[j].GetNd()+pos[i].GetNd())) {
1054         negativepair[i*10+cnegative[i]] =j;  //index
1055         positivepair[j*10+cpositive[j]] =i;
1056         cnegative[i]++;  //counters
1057         cpositive[j]++; 
1058         pairs[i][j]=100;
1059       }
1060     }
1061   }
1062   //
1063   for (Int_t i=0; i<np; i++) {
1064     Float_t yp=pos[i].GetY()*fYpitchSSD; 
1065     if (pos[i].GetQ()<3) continue;
1066     for (Int_t j=0; j<nn; j++) {
1067       if (neg[j].GetQ()<3) continue;
1068       if (cpositive[j]&&cnegative[i]) continue;
1069       Float_t yn=neg[j].GetY()*fYpitchSSD;
1070       Float_t zt=(2*fHlSSD*tanp + yp - yn)/(tann+tanp);
1071       Float_t yt=yn + tann*zt;
1072       zt-=fHlSSD; yt-=fHwSSD;
1073       if (TMath::Abs(yt)<fHwSSD+0.1)
1074       if (TMath::Abs(zt)<fHlSSD+0.15) {
1075         if (cnegative[i]==0) pos[i].SetNd(100);  // not available pair
1076         if (cpositive[j]==0) neg[j].SetNd(100);  // not available pair
1077         negativepair[i*10+cnegative[i]] =j;  //index
1078         positivepair[j*10+cpositive[j]] =i;
1079         cnegative[i]++;  //counters
1080         cpositive[j]++; 
1081         pairs[i][j]=100;
1082       }
1083     }
1084   }
1085   //
1086   Float_t lp[5];
1087   Int_t milab[10];
1088   Double_t ratio;
1089   
1090   //
1091   // sign gold tracks
1092   //
1093   for (Int_t ip=0;ip<np;ip++){
1094     Float_t ybest=1000,zbest=1000,qbest=0;
1095     //
1096     // select gold clusters
1097     if ( (cnegative[ip]==1) && cpositive[negativepair[10*ip]]==1){ 
1098       Float_t yp=pos[ip].GetY()*fYpitchSSD; 
1099       Int_t j = negativepair[10*ip];      
1100       ratio = (pos[ip].GetQ()-neg[j].GetQ())/(pos[ip].GetQ()+neg[j].GetQ());
1101       //
1102       Float_t yn=neg[j].GetY()*fYpitchSSD;
1103       Float_t zt=(2*fHlSSD*tanp + yp - yn)/(tann+tanp);
1104       Float_t yt=yn + tann*zt;
1105       zt-=fHlSSD; yt-=fHwSSD;
1106       ybest=yt; zbest=zt; 
1107       qbest=0.5*(pos[ip].GetQ()+neg[j].GetQ());
1108       lp[0]=-(-ybest+fYshift[fI]);
1109       lp[1]=  -zbest+fZshift[fI];
1110       lp[2]=0.0025*0.0025;  //SigmaY2
1111       lp[3]=0.110*0.110;  //SigmaZ2
1112       
1113       lp[4]=qbest;        //Q
1114       for (Int_t ilab=0;ilab<10;ilab++) milab[ilab]=-2;
1115       for (Int_t ilab=0;ilab<3;ilab++){
1116         milab[ilab] = pos[ip].GetLabel(ilab);
1117         milab[ilab+3] = neg[j].GetLabel(ilab);
1118       }
1119       //
1120       CheckLabels2(milab);
1121       milab[3]=(((ip<<10) + j)<<10) + idet; // pos|neg|det
1122       Int_t info[3] = {pos[ip].GetNd(),neg[j].GetNd(),fNlayer[fI]};
1123       AliITSclusterV2 * cl2 = new (cl[ncl]) AliITSclusterV2(milab,lp,info);
1124       ncl++;
1125       cl2->SetChargeRatio(ratio);       
1126       cl2->SetType(1);
1127       pairs[ip][j]=1;
1128       if ((pos[ip].GetNd()+neg[j].GetNd())>6){ //multi cluster
1129         cl2->SetType(2);
1130         pairs[ip][j]=2;
1131       }
1132       cused1[ip]++;
1133       cused2[j]++;
1134     }
1135   }
1136     
1137   for (Int_t ip=0;ip<np;ip++){
1138     Float_t ybest=1000,zbest=1000,qbest=0;
1139     //
1140     //
1141     // select "silber" cluster
1142     if ( cnegative[ip]==1 && cpositive[negativepair[10*ip]]==2){
1143       Int_t in  = negativepair[10*ip];
1144       Int_t ip2 = positivepair[10*in];
1145       if (ip2==ip) ip2 =  positivepair[10*in+1];
1146       Float_t pcharge = pos[ip].GetQ()+pos[ip2].GetQ();
1147       if (TMath::Abs(pcharge-neg[in].GetQ())<10){
1148         //
1149         // add first pair
1150         if (pairs[ip][in]==100){  //
1151           Float_t yp=pos[ip].GetY()*fYpitchSSD; 
1152           Float_t yn=neg[in].GetY()*fYpitchSSD;
1153           Float_t zt=(2*fHlSSD*tanp + yp - yn)/(tann+tanp);
1154           Float_t yt=yn + tann*zt;
1155           zt-=fHlSSD; yt-=fHwSSD;
1156           ybest =yt;  zbest=zt; 
1157           qbest =pos[ip].GetQ();
1158           lp[0]=-(-ybest+fYshift[fI]);
1159           lp[1]=  -zbest+fZshift[fI];
1160           lp[2]=0.0025*0.0025;  //SigmaY2
1161           lp[3]=0.110*0.110;  //SigmaZ2
1162           
1163           lp[4]=qbest;        //Q
1164           for (Int_t ilab=0;ilab<10;ilab++) milab[ilab]=-2;
1165           for (Int_t ilab=0;ilab<3;ilab++){
1166             milab[ilab] = pos[ip].GetLabel(ilab);
1167             milab[ilab+3] = neg[in].GetLabel(ilab);
1168           }
1169           //
1170           CheckLabels2(milab);
1171           ratio = (pos[ip].GetQ()-neg[in].GetQ())/(pos[ip].GetQ()+neg[in].GetQ());
1172           milab[3]=(((ip<<10) + in)<<10) + idet; // pos|neg|det
1173           Int_t info[3] = {pos[ip].GetNd(),neg[in].GetNd(),fNlayer[fI]};
1174           AliITSclusterV2 * cl2 = new (cl[ncl]) AliITSclusterV2(milab,lp,info);
1175           ncl++;
1176           cl2->SetChargeRatio(ratio);           
1177           cl2->SetType(5);
1178           pairs[ip][in] = 5;
1179           if ((pos[ip].GetNd()+neg[in].GetNd())>6){ //multi cluster
1180             cl2->SetType(6);
1181             pairs[ip][in] = 6;
1182           }
1183         }
1184         //
1185         // add second pair
1186         
1187         //      if (!(cused1[ip2] || cused2[in])){  //
1188         if (pairs[ip2][in]==100){
1189           Float_t yp=pos[ip2].GetY()*fYpitchSSD;
1190           Float_t yn=neg[in].GetY()*fYpitchSSD;
1191           Float_t zt=(2*fHlSSD*tanp + yp - yn)/(tann+tanp);
1192           Float_t yt=yn + tann*zt;
1193           zt-=fHlSSD; yt-=fHwSSD;
1194           ybest =yt;  zbest=zt; 
1195           qbest =pos[ip2].GetQ();
1196           lp[0]=-(-ybest+fYshift[fI]);
1197           lp[1]=  -zbest+fZshift[fI];
1198           lp[2]=0.0025*0.0025;  //SigmaY2
1199           lp[3]=0.110*0.110;  //SigmaZ2
1200           
1201           lp[4]=qbest;        //Q
1202           for (Int_t ilab=0;ilab<10;ilab++) milab[ilab]=-2;
1203           for (Int_t ilab=0;ilab<3;ilab++){
1204             milab[ilab] = pos[ip2].GetLabel(ilab);
1205             milab[ilab+3] = neg[in].GetLabel(ilab);
1206           }
1207           //
1208           CheckLabels2(milab);
1209           ratio = (pos[ip2].GetQ()-neg[in].GetQ())/(pos[ip2].GetQ()+neg[in].GetQ());
1210           milab[3]=(((ip2<<10) + in)<<10) + idet; // pos|neg|det
1211           Int_t info[3] = {pos[ip2].GetNd(),neg[in].GetNd(),fNlayer[fI]};
1212           AliITSclusterV2 *cl2 = new (cl[ncl]) AliITSclusterV2(milab,lp,info);
1213           ncl++;
1214           cl2->SetChargeRatio(ratio);           
1215           cl2->SetType(5);
1216           pairs[ip2][in] =5;
1217           if ((pos[ip2].GetNd()+neg[in].GetNd())>6){ //multi cluster
1218             cl2->SetType(6);
1219             pairs[ip2][in] =6;
1220           }
1221         }       
1222         cused1[ip]++;
1223         cused1[ip2]++;
1224         cused2[in]++;
1225       }
1226     }    
1227   }
1228   
1229   //  
1230   for (Int_t jn=0;jn<nn;jn++){
1231     if (cused2[jn]) continue;
1232     Float_t ybest=1000,zbest=1000,qbest=0;
1233     // select "silber" cluster
1234     if ( cpositive[jn]==1 && cnegative[positivepair[10*jn]]==2){
1235       Int_t ip  = positivepair[10*jn];
1236       Int_t jn2 = negativepair[10*ip];
1237       if (jn2==jn) jn2 =  negativepair[10*ip+1];
1238       Float_t pcharge = neg[jn].GetQ()+neg[jn2].GetQ();
1239       //
1240       if (TMath::Abs(pcharge-pos[ip].GetQ())<10){
1241         //
1242         // add first pair
1243         //      if (!(cused1[ip]||cused2[jn])){
1244         if (pairs[ip][jn]==100){
1245           Float_t yn=neg[jn].GetY()*fYpitchSSD; 
1246           Float_t yp=pos[ip].GetY()*fYpitchSSD;
1247           Float_t zt=(2*fHlSSD*tanp + yp - yn)/(tann+tanp);
1248           Float_t yt=yn + tann*zt;
1249           zt-=fHlSSD; yt-=fHwSSD;
1250           ybest =yt;  zbest=zt; 
1251           qbest =neg[jn].GetQ();
1252           lp[0]=-(-ybest+fYshift[fI]);
1253           lp[1]=  -zbest+fZshift[fI];
1254           lp[2]=0.0025*0.0025;  //SigmaY2
1255           lp[3]=0.110*0.110;  //SigmaZ2
1256           
1257           lp[4]=qbest;        //Q
1258           for (Int_t ilab=0;ilab<10;ilab++) milab[ilab]=-2;
1259           for (Int_t ilab=0;ilab<3;ilab++){
1260             milab[ilab] = pos[ip].GetLabel(ilab);
1261             milab[ilab+3] = neg[jn].GetLabel(ilab);
1262           }
1263           //
1264           CheckLabels2(milab);
1265           ratio = (pos[ip].GetQ()-neg[jn].GetQ())/(pos[ip].GetQ()+neg[jn].GetQ());
1266           milab[3]=(((ip<<10) + jn)<<10) + idet; // pos|neg|det
1267           Int_t info[3] = {pos[ip].GetNd(),neg[jn].GetNd(),fNlayer[fI]};
1268           AliITSclusterV2 * cl2 = new (cl[ncl]) AliITSclusterV2(milab,lp,info);
1269           ncl++;
1270           cl2->SetChargeRatio(ratio);           
1271           cl2->SetType(7);
1272           pairs[ip][jn] =7;
1273           if ((pos[ip].GetNd()+neg[jn].GetNd())>6){ //multi cluster
1274             cl2->SetType(8);
1275             pairs[ip][jn]=8;
1276           }
1277         }
1278         //
1279         // add second pair
1280         //      if (!(cused1[ip]||cused2[jn2])){
1281         if (pairs[ip][jn2]==100){
1282           Float_t yn=neg[jn2].GetY()*fYpitchSSD; 
1283           Double_t yp=pos[ip].GetY()*fYpitchSSD; 
1284           Double_t zt=(2*fHlSSD*tanp + yp - yn)/(tann+tanp);
1285           Double_t yt=yn + tann*zt;
1286           zt-=fHlSSD; yt-=fHwSSD;
1287           ybest =yt;  zbest=zt; 
1288           qbest =neg[jn2].GetQ();
1289           lp[0]=-(-ybest+fYshift[fI]);
1290           lp[1]=  -zbest+fZshift[fI];
1291           lp[2]=0.0025*0.0025;  //SigmaY2
1292           lp[3]=0.110*0.110;  //SigmaZ2
1293           
1294           lp[4]=qbest;        //Q
1295           for (Int_t ilab=0;ilab<10;ilab++) milab[ilab]=-2;
1296           for (Int_t ilab=0;ilab<3;ilab++){
1297             milab[ilab] = pos[ip].GetLabel(ilab);
1298             milab[ilab+3] = neg[jn2].GetLabel(ilab);
1299           }
1300           //
1301           CheckLabels2(milab);
1302           ratio = (pos[ip].GetQ()-neg[jn2].GetQ())/(pos[ip].GetQ()+neg[jn2].GetQ());
1303           milab[3]=(((ip<<10) + jn2)<<10) + idet; // pos|neg|det
1304           Int_t info[3] = {pos[ip].GetNd(),neg[jn2].GetNd(),fNlayer[fI]};
1305           AliITSclusterV2* cl2 = new (cl[ncl]) AliITSclusterV2(milab,lp,info);
1306           ncl++;
1307           cl2->SetChargeRatio(ratio);           
1308           pairs[ip][jn2]=7;
1309           cl2->SetType(7);
1310           if ((pos[ip].GetNd()+neg[jn2].GetNd())>6){ //multi cluster
1311             cl2->SetType(8);
1312             pairs[ip][jn2]=8;
1313           }
1314         }
1315         cused1[ip]++;
1316         cused2[jn]++;
1317         cused2[jn2]++;
1318       }
1319     }    
1320   }
1321   
1322   for (Int_t ip=0;ip<np;ip++){
1323     Float_t ybest=1000,zbest=1000,qbest=0;
1324     //
1325     // 2x2 clusters
1326     //
1327     if ( (cnegative[ip]<5) && cpositive[negativepair[10*ip]]<5){ 
1328       Float_t minchargediff =4.;
1329       Int_t j=-1;
1330       for (Int_t di=0;di<cnegative[ip];di++){
1331         Int_t   jc = negativepair[ip*10+di];
1332         Float_t chargedif = pos[ip].GetQ()-neg[jc].GetQ();
1333         if (TMath::Abs(chargedif)<minchargediff){
1334           j =jc;
1335           minchargediff = TMath::Abs(chargedif);
1336         }
1337       }
1338       if (j<0) continue;  // not proper cluster      
1339       Int_t count =0;
1340       for (Int_t di=0;di<cnegative[ip];di++){
1341         Int_t   jc = negativepair[ip*10+di];
1342         Float_t chargedif = pos[ip].GetQ()-neg[jc].GetQ();
1343         if (TMath::Abs(chargedif)<minchargediff+3.) count++;
1344       }
1345       if (count>1) continue;  // more than one "proper" cluster for positive
1346       //
1347       count =0;
1348       for (Int_t dj=0;dj<cpositive[j];dj++){
1349         Int_t   ic  = positivepair[j*10+dj];
1350         Float_t chargedif = pos[ic].GetQ()-neg[j].GetQ();
1351         if (TMath::Abs(chargedif)<minchargediff+3.) count++;
1352       }
1353       if (count>1) continue;  // more than one "proper" cluster for negative
1354       
1355       Int_t jp = 0;
1356       
1357       count =0;
1358       for (Int_t dj=0;dj<cnegative[jp];dj++){
1359         Int_t   ic = positivepair[jp*10+dj];
1360         Float_t chargedif = pos[ic].GetQ()-neg[jp].GetQ();
1361         if (TMath::Abs(chargedif)<minchargediff+4.) count++;
1362       }
1363       if (count>1) continue;   
1364       if (pairs[ip][j]<100) continue;
1365       //
1366       //almost gold clusters
1367       Float_t yp=pos[ip].GetY()*fYpitchSSD; 
1368       Float_t yn=neg[j].GetY()*fYpitchSSD;
1369       Float_t zt=(2*fHlSSD*tanp + yp - yn)/(tann+tanp);
1370       Float_t yt=yn + tann*zt;
1371       zt-=fHlSSD; yt-=fHwSSD;
1372       ybest=yt; zbest=zt; 
1373       qbest=0.5*(pos[ip].GetQ()+neg[j].GetQ());
1374       lp[0]=-(-ybest+fYshift[fI]);
1375       lp[1]=  -zbest+fZshift[fI];
1376       lp[2]=0.0025*0.0025;  //SigmaY2
1377       lp[3]=0.110*0.110;  //SigmaZ2     
1378       lp[4]=qbest;        //Q
1379       for (Int_t ilab=0;ilab<10;ilab++) milab[ilab]=-2;
1380       for (Int_t ilab=0;ilab<3;ilab++){
1381         milab[ilab] = pos[ip].GetLabel(ilab);
1382         milab[ilab+3] = neg[j].GetLabel(ilab);
1383       }
1384       //
1385       CheckLabels2(milab);
1386       ratio = (pos[ip].GetQ()-neg[j].GetQ())/(pos[ip].GetQ()+neg[j].GetQ());
1387       milab[3]=(((ip<<10) + j)<<10) + idet; // pos|neg|det
1388       Int_t info[3] = {pos[ip].GetNd(),neg[j].GetNd(),fNlayer[fI]};
1389       AliITSclusterV2 * cl2 = new (cl[ncl]) AliITSclusterV2(milab,lp,info);
1390       ncl++;
1391       cl2->SetChargeRatio(ratio);       
1392       cl2->SetType(10);
1393       pairs[ip][j]=10;
1394       if ((pos[ip].GetNd()+neg[j].GetNd())>6){ //multi cluster
1395         cl2->SetType(11);
1396         pairs[ip][j]=11;
1397       }
1398       cused1[ip]++;
1399       cused2[j]++;      
1400     }
1401
1402   }
1403   
1404   //  
1405   for (Int_t i=0; i<np; i++) {
1406     Float_t ybest=1000,zbest=1000,qbest=0;
1407     Float_t yp=pos[i].GetY()*fYpitchSSD; 
1408     if (pos[i].GetQ()<3) continue;
1409     for (Int_t j=0; j<nn; j++) {
1410     //    for (Int_t di = 0;di<cpositive[i];di++){
1411     //  Int_t j = negativepair[10*i+di];
1412       if (neg[j].GetQ()<3) continue;
1413       if (cused2[j]||cused1[i]) continue;      
1414       if (pairs[i][j]>0 &&pairs[i][j]<100) continue;
1415       ratio = (pos[i].GetQ()-neg[j].GetQ())/(pos[i].GetQ()+neg[j].GetQ());      
1416       Float_t yn=neg[j].GetY()*fYpitchSSD;
1417       Float_t zt=(2*fHlSSD*tanp + yp - yn)/(tann+tanp);
1418       Float_t yt=yn + tann*zt;
1419       zt-=fHlSSD; yt-=fHwSSD;
1420       if (TMath::Abs(yt)<fHwSSD+0.01)
1421       if (TMath::Abs(zt)<fHlSSD+0.01*(neg[j].GetNd()+pos[i].GetNd())) {
1422         ybest=yt; zbest=zt; 
1423         qbest=0.5*(pos[i].GetQ()+neg[j].GetQ());
1424         lp[0]=-(-ybest+fYshift[fI]);
1425         lp[1]=  -zbest+fZshift[fI];
1426         lp[2]=0.0025*0.0025;  //SigmaY2
1427         lp[3]=0.110*0.110;  //SigmaZ2
1428
1429         lp[4]=qbest;        //Q
1430         for (Int_t ilab=0;ilab<10;ilab++) milab[ilab]=-2;
1431         for (Int_t ilab=0;ilab<3;ilab++){
1432           milab[ilab] = pos[i].GetLabel(ilab);
1433           milab[ilab+3] = neg[j].GetLabel(ilab);
1434         }
1435         //
1436         CheckLabels2(milab);
1437         milab[3]=(((i<<10) + j)<<10) + idet; // pos|neg|det
1438         Int_t info[3] = {pos[i].GetNd(),neg[j].GetNd(),fNlayer[fI]};
1439         AliITSclusterV2 * cl2 = new (cl[ncl]) AliITSclusterV2(milab,lp,info); 
1440         ncl++;
1441         cl2->SetChargeRatio(ratio);
1442         cl2->SetType(100+cpositive[j]+cnegative[i]);
1443         //cl2->SetType(0);
1444         /*
1445           if (pairs[i][j]<100){
1446           printf("problem:- %d\n", pairs[i][j]);
1447           }
1448           if (cnegative[i]<2&&cpositive[j]<2){
1449           printf("problem:- %d\n", pairs[i][j]);
1450           }
1451         */
1452       }
1453     }
1454   }
1455
1456 //   for (Int_t i=0; i<1000; i++) delete [] pairs[i];
1457 //   delete [] pairs;
1458
1459 }
1460
1461
1462 void AliITSclustererV2::
1463 FindClustersSSD(const TClonesArray *alldigits, TClonesArray *clusters) {
1464   //------------------------------------------------------------
1465   // Actual SSD cluster finder
1466   //------------------------------------------------------------
1467   Int_t smaxall=alldigits->GetEntriesFast();
1468   if (smaxall==0) return;
1469   TObjArray *digits = new TObjArray;
1470   for (Int_t i=0;i<smaxall; i++){
1471     AliITSdigitSSD *d=(AliITSdigitSSD*)alldigits->UncheckedAt(i);
1472     if (d->GetSignal()<3) continue;
1473     digits->AddLast(d);
1474   }
1475   Int_t smax = digits->GetEntriesFast();
1476   if (smax==0) return;
1477   
1478   const Int_t MAX=1000;
1479   Int_t np=0, nn=0; 
1480   Ali1Dcluster pos[MAX], neg[MAX];
1481   Float_t y=0., q=0., qmax=0.; 
1482   Int_t lab[4]={-2,-2,-2,-2};
1483   
1484   AliITSdigitSSD *d=(AliITSdigitSSD*)digits->UncheckedAt(0);
1485   q += d->GetSignal();
1486   y += d->GetCoord2()*d->GetSignal();
1487   qmax=d->GetSignal();
1488   lab[0]=d->GetTrack(0); lab[1]=d->GetTrack(1); lab[2]=d->GetTrack(2);
1489   Int_t curr=d->GetCoord2();
1490   Int_t flag=d->GetCoord1();
1491   Int_t *n=&nn;
1492   Ali1Dcluster *c=neg;
1493   Int_t nd=1;
1494   Int_t milab[10];
1495   for (Int_t ilab=0;ilab<10;ilab++){
1496     milab[ilab]=-2;
1497   }
1498   milab[0]=d->GetTrack(0); milab[1]=d->GetTrack(1); milab[2]=d->GetTrack(2);
1499
1500   for (Int_t s=1; s<smax; s++) {
1501       d=(AliITSdigitSSD*)digits->UncheckedAt(s);      
1502       Int_t strip=d->GetCoord2();
1503       if ((strip-curr) > 1 || flag!=d->GetCoord1()) {
1504          c[*n].SetY(y/q);
1505          c[*n].SetQ(q);
1506          c[*n].SetNd(nd);
1507          CheckLabels2(milab);
1508          c[*n].SetLabels(milab);
1509          //Split suspiciously big cluster
1510          /*
1511          if (nd>10&&nd<16){
1512            c[*n].SetY(y/q-0.3*nd);
1513            c[*n].SetQ(0.5*q);
1514            (*n)++;
1515            if (*n==MAX) {
1516              Error("FindClustersSSD","Too many 1D clusters !");
1517               return;
1518            }
1519            c[*n].SetY(y/q-0.0*nd);
1520            c[*n].SetQ(0.5*q);
1521            c[*n].SetNd(nd);
1522            (*n)++;
1523            if (*n==MAX) {
1524              Error("FindClustersSSD","Too many 1D clusters !");
1525               return;
1526            }
1527            //
1528            c[*n].SetY(y/q+0.3*nd);
1529            c[*n].SetQ(0.5*q);
1530            c[*n].SetNd(nd);
1531            c[*n].SetLabels(milab);
1532          }
1533          else{
1534          */
1535          if (nd>4&&nd<25) {
1536            c[*n].SetY(y/q-0.25*nd);
1537            c[*n].SetQ(0.5*q);
1538            (*n)++;
1539            if (*n==MAX) {
1540              Error("FindClustersSSD","Too many 1D clusters !");
1541              return;
1542            }
1543            c[*n].SetY(y/q+0.25*nd);
1544            c[*n].SetQ(0.5*q);
1545            c[*n].SetNd(nd);
1546            c[*n].SetLabels(milab);
1547          }       
1548          (*n)++;
1549          if (*n==MAX) {
1550           Error("FindClustersSSD","Too many 1D clusters !");
1551           return;
1552          }
1553          y=q=qmax=0.;
1554          nd=0;
1555          lab[0]=lab[1]=lab[2]=-2;
1556          //
1557          for (Int_t ilab=0;ilab<10;ilab++){
1558            milab[ilab]=-2;
1559          }
1560          //
1561          if (flag!=d->GetCoord1()) { n=&np; c=pos; }
1562       }
1563       flag=d->GetCoord1();
1564       q += d->GetSignal();
1565       y += d->GetCoord2()*d->GetSignal();
1566       nd++;
1567       if (d->GetSignal()>qmax) {
1568          qmax=d->GetSignal();
1569          lab[0]=d->GetTrack(0); lab[1]=d->GetTrack(1); lab[2]=d->GetTrack(2);
1570       }
1571       for (Int_t ilab=0;ilab<10;ilab++) {
1572         if (d->GetTrack(ilab)>=0) AddLabel(milab, (d->GetTrack(ilab))); 
1573       }
1574       curr=strip;
1575   }
1576   c[*n].SetY(y/q);
1577   c[*n].SetQ(q);
1578   c[*n].SetNd(nd);
1579   c[*n].SetLabels(lab);
1580   //Split suspiciously big cluster
1581   if (nd>4 && nd<25) {
1582      c[*n].SetY(y/q-0.25*nd);
1583      c[*n].SetQ(0.5*q);
1584      (*n)++;
1585      if (*n==MAX) {
1586         Error("FindClustersSSD","Too many 1D clusters !");
1587         return;
1588      }
1589      c[*n].SetY(y/q+0.25*nd);
1590      c[*n].SetQ(0.5*q);
1591      c[*n].SetNd(nd);
1592      c[*n].SetLabels(lab);
1593   }
1594   (*n)++;
1595   if (*n==MAX) {
1596      Error("FindClustersSSD","Too many 1D clusters !");
1597      return;
1598   }
1599
1600   FindClustersSSD(neg, nn, pos, np, clusters);
1601 }
1602
1603 void AliITSclustererV2::FindClustersSSD(AliITSRawStream* input, 
1604                                         TClonesArray** clusters) 
1605 {
1606   //------------------------------------------------------------
1607   // Actual SSD cluster finder for raw data
1608   //------------------------------------------------------------
1609   Int_t nClustersSSD = 0;
1610   const Int_t MAX = 1000;
1611   Ali1Dcluster clusters1D[2][MAX];
1612   Int_t nClusters[2] = {0, 0};
1613   Int_t lab[3]={-2,-2,-2};
1614   Float_t q = 0.;
1615   Float_t y = 0.;
1616   Int_t nDigits = 0;
1617   Int_t prevStrip = -1;
1618   Int_t prevFlag = -1;
1619   Int_t prevModule = -1;
1620
1621   // read raw data input stream
1622   while (kTRUE) {
1623     Bool_t next = input->Next();
1624
1625     if(input->GetSignal()<3 && next) continue;
1626     // check if a new cluster starts
1627     Int_t strip = input->GetCoord2();
1628     Int_t flag = input->GetCoord1();
1629     if ((!next || (input->GetModuleID() != prevModule)||
1630          (strip-prevStrip > 1) || (flag != prevFlag)) &&
1631         (nDigits > 0)) {
1632       if (nClusters[prevFlag] == MAX) {
1633         Error("FindClustersSSD", "Too many 1D clusters !");
1634         return;
1635       }
1636       Ali1Dcluster& cluster = clusters1D[prevFlag][nClusters[prevFlag]++];
1637       cluster.SetY(y/q);
1638       cluster.SetQ(q);
1639       cluster.SetNd(nDigits);
1640       cluster.SetLabels(lab);
1641
1642       //Split suspiciously big cluster
1643       if (nDigits > 4&&nDigits < 25) {
1644         cluster.SetY(y/q - 0.25*nDigits);
1645         cluster.SetQ(0.5*q);
1646         if (nClusters[prevFlag] == MAX) {
1647           Error("FindClustersSSD", "Too many 1D clusters !");
1648           return;
1649         }
1650         Ali1Dcluster& cluster2 = clusters1D[prevFlag][nClusters[prevFlag]++];
1651         cluster2.SetY(y/q + 0.25*nDigits);
1652         cluster2.SetQ(0.5*q);
1653         cluster2.SetNd(nDigits);
1654         cluster2.SetLabels(lab);
1655       }
1656       y = q = 0.;
1657       nDigits = 0;
1658     }
1659
1660     if (!next || (input->GetModuleID() != prevModule)) {
1661       Int_t iModule = prevModule;
1662
1663       // when all data from a module was read, search for clusters
1664       if (prevFlag >= 0) {
1665         clusters[iModule] = new TClonesArray("AliITSclusterV2");
1666         fI = iModule;
1667         FindClustersSSD(&clusters1D[0][0], nClusters[0], 
1668                         &clusters1D[1][0], nClusters[1], clusters[iModule]);
1669         Int_t nClusters = clusters[iModule]->GetEntriesFast();
1670         nClustersSSD += nClusters;
1671       }
1672
1673       if (!next) break;
1674       nClusters[0] = nClusters[1] = 0;
1675       y = q = 0.;
1676       nDigits = 0;
1677     }
1678
1679     // add digit to current cluster
1680     q += input->GetSignal();
1681     y += strip * input->GetSignal();
1682     nDigits++;
1683     prevStrip = strip;
1684     prevFlag = flag;
1685     prevModule = input->GetModuleID();
1686
1687   }
1688
1689   Info("FindClustersSSD", "found clusters in ITS SSD: %d", nClustersSSD);
1690 }
1691
1692 #else   //V1
1693
1694 #include "AliITSDetType.h"
1695 #include "AliITS.h"
1696 #include "AliITSsegmentationSPD.h"
1697 #include "AliITSClusterFinderSPD.h"
1698
1699 #include "AliITSresponseSDD.h"
1700 #include "AliITSsegmentationSDD.h"
1701 #include "AliITSClusterFinderSDD.h"
1702
1703 #include "AliITSsegmentationSSD.h"
1704 #include "AliITSClusterFinderSSD.h"
1705
1706
1707 void AliITSclustererV2::
1708 FindClustersSPD(const TClonesArray *digits, TClonesArray *clusters) {
1709   //------------------------------------------------------------
1710   // Actual SPD cluster finding based on AliITSClusterFinderSPD
1711   //------------------------------------------------------------
1712   static AliITS *its=(AliITS*)gAlice->GetModule("ITS");
1713   static TClonesArray *points=its->RecPoints();
1714   static AliITSsegmentationSPD *seg=
1715          (AliITSsegmentationSPD *)its->DetType(0)->GetSegmentationModel();
1716   static AliITSClusterFinderSPD cf(seg, (TClonesArray*)digits, points);
1717
1718   cf.FindRawClusters(fI);
1719   RecPoints2Clusters(points, fI, clusters);
1720   its->ResetRecPoints();
1721
1722 }
1723
1724 void AliITSclustererV2::
1725 FindClustersSDD(const TClonesArray *digits, TClonesArray *clusters) {
1726   //------------------------------------------------------------
1727   // Actual SDD cluster finding based on AliITSClusterFinderSDD
1728   //------------------------------------------------------------
1729   static AliITS *its=(AliITS*)gAlice->GetModule("ITS");
1730   static TClonesArray *points=its->RecPoints();
1731   static AliITSresponseSDD *resp=
1732         (AliITSresponseSDD *)its->DetType(1)->GetResponseModel();
1733   static AliITSsegmentationSDD *seg=
1734          (AliITSsegmentationSDD *)its->DetType(1)->GetSegmentationModel();
1735   static AliITSClusterFinderSDD 
1736          cf(seg,resp,(TClonesArray*)digits,its->ClustersAddress(1));
1737
1738   cf.FindRawClusters(fI);
1739   Int_t nc=points->GetEntriesFast();
1740   while (nc--) { //To be consistent with the SSD cluster charges
1741      AliITSRecPoint *p=(AliITSRecPoint*)points->UncheckedAt(nc);
1742      p->SetQ(p->GetQ()/12.);
1743   }
1744   RecPoints2Clusters(points, fI, clusters);
1745   its->ResetClusters(1);
1746   its->ResetRecPoints();
1747
1748 }
1749
1750 void AliITSclustererV2::
1751 FindClustersSSD(const TClonesArray *digits, TClonesArray *clusters) {
1752   //------------------------------------------------------------
1753   // Actual SSD cluster finding based on AliITSClusterFinderSSD
1754   //------------------------------------------------------------
1755   static AliITS *its=(AliITS*)gAlice->GetModule("ITS");
1756   static TClonesArray *points=its->RecPoints();
1757   static AliITSsegmentationSSD *seg=
1758          (AliITSsegmentationSSD *)its->DetType(2)->GetSegmentationModel();
1759   static AliITSClusterFinderSSD cf(seg,(TClonesArray*)digits);
1760
1761   cf.FindRawClusters(fI);
1762   RecPoints2Clusters(points, fI, clusters);
1763   its->ResetRecPoints();
1764
1765 }
1766
1767 #endif