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