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