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