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