]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSclustererV2.cxx
Using Fatal instead of Error+exit
[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 "AliITSdigit.h"
23
24 ClassImp(AliITSclustererV2)
25
26 extern AliRun *gAlice;
27
28 AliITSclustererV2::AliITSclustererV2(const AliITSgeom *geom) {
29   //------------------------------------------------------------
30   // Standard constructor
31   //------------------------------------------------------------
32   AliITSgeom *g=(AliITSgeom*)geom;
33
34   fEvent=0;
35   fI=0;
36
37   Int_t mmax=geom->GetIndexMax();
38   if (mmax>2200) {
39      Fatal("AliITSclustererV2","Too many ITS subdetectors !"); 
40   }
41   Int_t m;
42   for (m=0; m<mmax; m++) {
43      Int_t lay,lad,det; g->GetModuleId(m,lay,lad,det);
44      Float_t x,y,z;     g->GetTrans(lay,lad,det,x,y,z); 
45      Double_t rot[9];   g->GetRotMatrix(lay,lad,det,rot);
46      fYshift[m] = x*rot[0] + y*rot[1];
47      fZshift[m] = (Double_t)z;
48      fNdet[m] = (lad-1)*g->GetNdetectors(lay) + (det-1);
49   }
50   fNModules = g->GetIndexMax();
51
52   //SPD geometry  
53   fLastSPD1=g->GetModuleIndex(2,1,1)-1;
54   fNySPD=256; fNzSPD=160;
55   fYpitchSPD=0.0050;
56   fZ1pitchSPD=0.0425; fZ2pitchSPD=0.0625;
57   fHwSPD=0.64; fHlSPD=3.48;
58   fYSPD[0]=0.5*fYpitchSPD;
59   for (m=1; m<fNySPD; m++) fYSPD[m]=fYSPD[m-1]+fYpitchSPD; 
60   fZSPD[0]=fZ1pitchSPD;
61   for (m=1; m<fNzSPD; m++) {
62     Double_t dz=fZ1pitchSPD;
63     if (m==31 || m==32 || m==63  || m==64  || m==95 || m==96 || 
64         m==127 || m==128) dz=fZ2pitchSPD; 
65     fZSPD[m]=fZSPD[m-1]+dz;
66   }
67   for (m=0; m<fNzSPD; m++) {
68     Double_t dz=0.5*fZ1pitchSPD;
69     if (m==31 || m==32 || m==63  || m==64  || m==95 || m==96 || 
70         m==127 || m==128) dz=0.5*fZ2pitchSPD; 
71     fZSPD[m]-=dz;
72   }
73
74   //SDD geometry 
75   fNySDD=256; fNzSDD=256;
76   fYpitchSDD=0.01825;
77   fZpitchSDD=0.02940;
78   fHwSDD=3.5085; fHlSDD=3.7632;
79   fYoffSDD=0.0425;
80
81   //SSD geometry
82   fLastSSD1=g->GetModuleIndex(6,1,1)-1;
83   fYpitchSSD=0.0095;
84   fHwSSD=3.65;
85   fHlSSD=2.00;
86   fTanP=0.0275;
87   fTanN=0.0075;
88 }
89
90 Int_t AliITSclustererV2::Digits2Clusters(TTree *dTree, TTree *cTree) {
91   //------------------------------------------------------------
92   // This function creates ITS clusters
93   //------------------------------------------------------------
94   Int_t ncl=0;
95
96   if (!dTree) {
97     Error("Digits2Clusters","Can't get the tree with digits !");
98     return 1;
99   }
100
101   TClonesArray *digitsSPD=new TClonesArray("AliITSdigitSPD",3000);
102   dTree->SetBranchAddress("ITSDigitsSPD",&digitsSPD);
103   TClonesArray *digitsSDD=new TClonesArray("AliITSdigitSDD",3000);
104   dTree->SetBranchAddress("ITSDigitsSDD",&digitsSDD);
105   TClonesArray *digitsSSD=new TClonesArray("AliITSdigitSSD",3000);
106   dTree->SetBranchAddress("ITSDigitsSSD",&digitsSSD);
107
108   TClonesArray *clusters=new TClonesArray("AliITSclusterV2",1000);
109   TBranch *branch=cTree->GetBranch("Clusters");
110   if (!branch) cTree->Branch("Clusters",&clusters);
111   else branch->SetAddress(&clusters);
112
113   Int_t mmax=(Int_t)dTree->GetEntries();
114
115   for (fI=0; fI<mmax; fI++) {
116     dTree->GetEvent(fI);
117
118     if     (digitsSPD->GetEntriesFast()!=0) 
119                           FindClustersSPD(digitsSPD,clusters);
120     else if(digitsSDD->GetEntriesFast()!=0) 
121                           FindClustersSDD(digitsSDD,clusters);
122     else if(digitsSSD->GetEntriesFast()!=0) 
123                           FindClustersSSD(digitsSSD,clusters);
124
125     ncl+=clusters->GetEntriesFast();
126
127     cTree->Fill();
128
129     digitsSPD->Clear();
130     digitsSDD->Clear();
131     digitsSSD->Clear();
132     clusters->Clear();
133   }
134
135   //cTree->Write();
136
137   delete clusters;
138
139   delete digitsSPD;
140   delete digitsSDD;
141   delete digitsSSD;
142
143   //delete dTree;
144
145   Info("Digits2Clusters","Number of found clusters : %d",ncl);
146
147   return 0;
148 }
149
150 void AliITSclustererV2::Digits2Clusters(AliRawReader* rawReader) {
151   //------------------------------------------------------------
152   // This function creates ITS clusters from raw data
153   //------------------------------------------------------------
154   AliRunLoader* runLoader = AliRunLoader::GetRunLoader();
155   if (!runLoader) {
156     Error("Digits2Clusters", "no run loader found");
157     return;
158   }
159   AliLoader* itsLoader = runLoader->GetLoader("ITSLoader");
160   if (!itsLoader) {
161     Error("Digits2Clusters", "no loader for ITS found");
162     return;
163   }
164   if (!itsLoader->TreeR()) itsLoader->MakeTree("R");
165   TTree* cTree = itsLoader->TreeR();
166
167   TClonesArray *array=new TClonesArray("AliITSclusterV2",1000);
168   cTree->Branch("Clusters",&array);
169   delete array;
170
171   TClonesArray** clusters = new TClonesArray*[fNModules]; 
172   // one TClonesArray per module
173
174   rawReader->Reset();
175   AliITSRawStreamSPD inputSPD(rawReader);
176   FindClustersSPD(&inputSPD, clusters);
177
178   rawReader->Reset();
179   AliITSRawStreamSDD inputSDD(rawReader);
180   FindClustersSDD(&inputSDD, clusters);
181
182   rawReader->Reset();
183   AliITSRawStreamSSD inputSSD(rawReader);
184   FindClustersSSD(&inputSSD, clusters);
185
186   // write all clusters to the tree
187   Int_t nClusters = 0;
188   for (Int_t iModule = 0; iModule < fNModules; iModule++) {
189     TClonesArray* array = clusters[iModule];
190     if (!array) {
191       Error("Digits2Clusters", "data for module %d missing!", iModule);
192       array = new TClonesArray("AliITSclusterV2");
193     }
194     cTree->SetBranchAddress("Clusters", &array);
195     cTree->Fill();
196     nClusters += array->GetEntriesFast();
197     delete array;
198   }
199   itsLoader->WriteRecPoints("OVERWRITE");
200
201   Info("Digits2Clusters", "total number of found clusters in ITS: %d\n", 
202        nClusters);
203 }
204
205 //**** Fast clusters *******************************
206 #include "TParticle.h"
207
208 #include "AliITS.h"
209 #include "AliITSmodule.h"
210 #include "AliITSRecPoint.h"
211 #include "AliITSsimulationFastPoints.h"
212 #include "AliITSRecPoint.h"
213
214 static void CheckLabels(Int_t lab[3]) {
215   //------------------------------------------------------------
216   // Tries to find mother's labels
217   //------------------------------------------------------------
218     Int_t label=lab[0];
219     if (label>=0) {
220        TParticle *part=(TParticle*)gAlice->Particle(label);
221        label=-3;
222        while (part->P() < 0.005) {
223           Int_t m=part->GetFirstMother();
224           if (m<0) {
225              Info("CheckLabels","Primary momentum: %f",part->P()); 
226              break;
227           }
228           if (part->GetStatusCode()>0) {
229              Info("CheckLabels","Primary momentum: %f",part->P()); 
230              break;
231           }
232           label=m;
233           part=(TParticle*)gAlice->Particle(label);
234         }
235         if      (lab[1]<0) lab[1]=label;
236         else if (lab[2]<0) lab[2]=label;
237         else ;//cerr<<"CheckLabels : No empty labels !\n";
238     }
239 }
240
241 void AliITSclustererV2::RecPoints2Clusters
242 (const TClonesArray *points, Int_t idx, TClonesArray *clusters) {
243   //------------------------------------------------------------
244   // Conversion AliITSRecPoint -> AliITSclusterV2 for the ITS 
245   // subdetector indexed by idx 
246   //------------------------------------------------------------
247   TClonesArray &cl=*clusters;
248   Int_t ncl=points->GetEntriesFast();
249   for (Int_t i=0; i<ncl; i++) {
250     AliITSRecPoint *p = (AliITSRecPoint *)points->UncheckedAt(i);
251     Float_t lp[5];
252     lp[0]=-p->GetX()-fYshift[idx]; if (idx<=fLastSPD1) lp[0]*=-1; //SPD1
253     lp[1]=p->GetZ()+fZshift[idx];
254     lp[2]=p->GetSigmaX2();
255     lp[3]=p->GetSigmaZ2();
256     lp[4]=p->GetQ();
257     Int_t lab[4]; 
258     lab[0]=p->GetLabel(0); lab[1]=p->GetLabel(1); lab[2]=p->GetLabel(2);
259     lab[3]=fNdet[idx];
260     CheckLabels(lab);
261     new (cl[i]) AliITSclusterV2(lab,lp);
262   }  
263
264
265 Int_t AliITSclustererV2::Hits2Clusters(TTree *hTree, TTree *cTree) {
266   //------------------------------------------------------------
267   // This function creates ITS clusters
268   //------------------------------------------------------------
269   if (!gAlice) {
270      Error("Hits2Clusters","gAlice==0 !");
271      return 1;
272   }
273
274   AliITS *its  = (AliITS*)gAlice->GetModule("ITS");
275   if (!its) { 
276      Error("Hits2Clusters","Can't find the ITS !"); 
277      return 2; 
278   }
279   AliITSgeom *geom=its->GetITSgeom();
280   Int_t mmax=geom->GetIndexMax();
281
282   its->InitModules(-1,mmax);
283   its->FillModules(hTree,0);
284
285   TClonesArray *clusters=new TClonesArray("AliITSclusterV2",1000);
286   TBranch *branch=cTree->GetBranch("Clusters");
287   if (!branch) cTree->Branch("Clusters",&clusters);
288   else branch->SetAddress(&clusters);
289
290   static TClonesArray *points=its->RecPoints();
291   AliITSsimulationFastPoints sim;
292   Int_t ncl=0;
293   for (Int_t m=0; m<mmax; m++) {
294     AliITSmodule *mod=its->GetModule(m);      
295     sim.CreateFastRecPoints(mod,m,gRandom);      
296
297     RecPoints2Clusters(points, m, clusters);
298     its->ResetRecPoints();
299
300     ncl+=clusters->GetEntriesFast();
301     cTree->Fill();
302     clusters->Clear();
303   }
304
305   Info("Hits2Clusters","Number of found fast clusters : %d",ncl);
306
307   //cTree->Write();
308
309   delete clusters;
310
311   return 0;
312 }
313
314 //***********************************
315
316 #ifndef V1
317
318 void AliITSclustererV2:: 
319 FindCluster(Int_t k,Int_t maxz,AliBin *bins,Int_t &n,Int_t *idx) {
320   //------------------------------------------------------------
321   // returns an array of indices of digits belonging to the cluster
322   // (needed when the segmentation is not regular) 
323   //------------------------------------------------------------
324   if (n<200) idx[n++]=bins[k].GetIndex();
325   bins[k].Use();
326
327   if (bins[k-maxz].IsNotUsed()) FindCluster(k-maxz,maxz,bins,n,idx);
328   if (bins[k-1   ].IsNotUsed()) FindCluster(k-1   ,maxz,bins,n,idx);
329   if (bins[k+maxz].IsNotUsed()) FindCluster(k+maxz,maxz,bins,n,idx);
330   if (bins[k+1   ].IsNotUsed()) FindCluster(k+1   ,maxz,bins,n,idx);
331   /*
332   if (bins[k-maxz-1].IsNotUsed()) FindCluster(k-maxz-1,maxz,bins,n,idx);
333   if (bins[k-maxz+1].IsNotUsed()) FindCluster(k-maxz+1,maxz,bins,n,idx);
334   if (bins[k+maxz-1].IsNotUsed()) FindCluster(k+maxz-1,maxz,bins,n,idx);
335   if (bins[k+maxz+1].IsNotUsed()) FindCluster(k+maxz+1,maxz,bins,n,idx);
336   */
337 }
338
339 void AliITSclustererV2::
340 FindClustersSPD(const TClonesArray *digits, TClonesArray *clusters) {
341   //------------------------------------------------------------
342   // Actual SPD cluster finder
343   //------------------------------------------------------------
344   Int_t kNzBins = fNzSPD + 2;
345   const Int_t kMAXBIN=kNzBins*(fNySPD+2);
346
347   Int_t ndigits=digits->GetEntriesFast();
348   AliBin *bins=new AliBin[kMAXBIN];
349
350   Int_t k;
351   AliITSdigitSPD *d=0;
352   for (k=0; k<ndigits; k++) {
353      d=(AliITSdigitSPD*)digits->UncheckedAt(k);
354      Int_t i=d->GetCoord2()+1;   //y
355      Int_t j=d->GetCoord1()+1;
356      bins[i*kNzBins+j].SetIndex(k);
357      bins[i*kNzBins+j].SetMask(1);
358   }
359    
360   Int_t n=0; TClonesArray &cl=*clusters;
361   for (k=0; k<kMAXBIN; k++) {
362      if (!bins[k].IsNotUsed()) continue;
363      Int_t ni=0, idx[200];
364      FindCluster(k,kNzBins,bins,ni,idx);
365      if (ni==200) {
366         Info("FindClustersSPD","Too big cluster !"); 
367         continue;
368      }
369      Int_t lab[4]; 
370      lab[0]=-2;
371      lab[1]=-2;
372      lab[2]=-2;
373      lab[3]=fNdet[fI];
374
375      d=(AliITSdigitSPD*)digits->UncheckedAt(idx[0]);
376      Int_t ymin=d->GetCoord2(),ymax=ymin;
377      Int_t zmin=d->GetCoord1(),zmax=zmin;
378      Float_t y=0.,z=0.,q=0.;
379      for (Int_t l=0; l<ni; l++) {
380         d=(AliITSdigitSPD*)digits->UncheckedAt(idx[l]);
381
382         if (ymin > d->GetCoord2()) ymin=d->GetCoord2();
383         if (ymax < d->GetCoord2()) ymax=d->GetCoord2();
384         if (zmin > d->GetCoord1()) zmin=d->GetCoord1();
385         if (zmax < d->GetCoord1()) zmax=d->GetCoord1();
386
387         Int_t lab0=(d->GetTracks())[0];      
388         if (lab0>=0) {
389           if (lab[0]<0) {
390              lab[0]=lab0;
391           } else if (lab[1]<0) {
392             if (lab0!=lab[0]) lab[1]=lab0;
393           } else if (lab[2]<0) {
394             if (lab0!=lab[0])
395             if (lab0!=lab[1]) lab[2]=lab0;
396           }
397         }
398         Float_t qq=d->GetSignal();
399         y+=qq*fYSPD[d->GetCoord2()]; z+=qq*fZSPD[d->GetCoord1()]; q+=qq;   
400      }
401      y/=q; z/=q;
402      y-=fHwSPD; z-=fHlSPD;
403
404      Float_t lp[5];
405      lp[0]=-y-fYshift[fI]; if (fI<=fLastSPD1) lp[0]=-lp[0];
406      lp[1]= z+fZshift[fI];
407      lp[2]= fYpitchSPD*fYpitchSPD/12.;
408      lp[3]= fZ1pitchSPD*fZ1pitchSPD/12.;
409      //lp[4]= q;
410      lp[4]= (zmax-zmin+1)*100 + (ymax-ymin+1);
411
412      //CheckLabels(lab);
413      d=(AliITSdigitSPD*)digits->UncheckedAt(idx[0]);
414      new (cl[n]) AliITSclusterV2(lab,lp); n++; 
415   }
416
417   delete bins;
418 }
419
420 void AliITSclustererV2::FindClustersSPD(AliITSRawStream* input, 
421                                         TClonesArray** clusters) 
422 {
423   //------------------------------------------------------------
424   // Actual SPD cluster finder for raw data
425   //------------------------------------------------------------
426
427   Int_t nClustersSPD = 0;
428   Int_t kNzBins = fNzSPD + 2;
429   Int_t kNyBins = fNySPD + 2;
430   Int_t kMaxBin = kNzBins * kNyBins;
431   AliBin* bins = NULL;
432
433   // read raw data input stream
434   while (kTRUE) {
435     Bool_t next = input->Next();
436     if (!next || input->IsNewModule()) {
437       Int_t iModule = input->GetPrevModuleID();
438
439       // when all data from a module was read, search for clusters
440       if (bins) { 
441         clusters[iModule] = new TClonesArray("AliITSclusterV2");
442         Int_t nClusters = 0;
443
444         for (Int_t iBin = 0; iBin < kMaxBin; iBin++) {
445           if (bins[iBin].IsUsed()) continue;
446           Int_t nBins = 0;
447           Int_t idxBins[200];
448           FindCluster(iBin, kNzBins, bins, nBins, idxBins);
449           if (nBins == 200) {
450             Error("FindClustersSPD", "SPD: Too big cluster !\n"); 
451             continue;
452           }
453
454           Int_t label[4]; 
455           label[0] = -2;
456           label[1] = -2;
457           label[2] = -2;
458 //        label[3] = iModule;
459           label[3] = fNdet[iModule];
460
461           Int_t ymin = (idxBins[0] / kNzBins) - 1;
462           Int_t ymax = ymin;
463           Int_t zmin = (idxBins[0] % kNzBins) - 1;
464           Int_t zmax = zmin;
465           Float_t y = 0.;
466           Float_t z = 0.;
467           Float_t q = 0.;
468           for (Int_t idx = 0; idx < nBins; idx++) {
469             Int_t iy = (idxBins[idx] / kNzBins) - 1;
470             Int_t iz = (idxBins[idx] % kNzBins) - 1;
471             if (ymin > iy) ymin = iy;
472             if (ymax < iy) ymax = iy;
473             if (zmin > iz) zmin = iz;
474             if (zmax < iz) zmax = iz;
475
476             Float_t qBin = bins[idxBins[idx]].GetQ();
477             y += qBin * fYSPD[iy]; 
478             z += qBin * fZSPD[iz]; 
479             q += qBin;   
480           }
481           y /= q; 
482           z /= q;
483           y -= fHwSPD; 
484           z -= fHlSPD;
485
486           Float_t hit[5];  // y, z, sigma(y)^2, sigma(z)^2, charge
487           hit[0] = -y-fYshift[iModule]; 
488           if (iModule <= fLastSPD1) hit[0] = -hit[0];
489           hit[1] = z+fZshift[iModule];
490           hit[2] = fYpitchSPD*fYpitchSPD/12.;
491           hit[3] = fZ1pitchSPD*fZ1pitchSPD/12.;
492 //        hit[4] = q;
493           hit[4] = (zmax-zmin+1)*100 + (ymax-ymin+1);
494
495           //CheckLabels(label);
496           new (clusters[iModule]->AddrAt(nClusters)) 
497             AliITSclusterV2(label, hit); 
498           nClusters++;
499         }
500
501         nClustersSPD += nClusters;
502         delete bins;
503       }
504
505       if (!next) break;
506       bins = new AliBin[kMaxBin];
507     }
508
509     // fill the current digit into the bins array
510     Int_t index = (input->GetCoord2()+1) * kNzBins + (input->GetCoord1()+1);
511     bins[index].SetIndex(index);
512     bins[index].SetMask(1);
513     bins[index].SetQ(1);
514   }
515
516   Info("FindClustersSPD", "found clusters in ITS SPD: %d", nClustersSPD);
517 }
518
519
520 Bool_t AliITSclustererV2::IsMaximum(Int_t k,Int_t max,const AliBin *bins) {
521   //------------------------------------------------------------
522   //is this a local maximum ?
523   //------------------------------------------------------------
524   UShort_t q=bins[k].GetQ();
525   if (q==1023) return kFALSE;
526   if (bins[k-max].GetQ() > q) return kFALSE;
527   if (bins[k-1  ].GetQ() > q) return kFALSE; 
528   if (bins[k+max].GetQ() > q) return kFALSE; 
529   if (bins[k+1  ].GetQ() > q) return kFALSE; 
530   if (bins[k-max-1].GetQ() > q) return kFALSE;
531   if (bins[k+max-1].GetQ() > q) return kFALSE; 
532   if (bins[k+max+1].GetQ() > q) return kFALSE; 
533   if (bins[k-max+1].GetQ() > q) return kFALSE;
534   return kTRUE; 
535 }
536
537 void AliITSclustererV2::
538 FindPeaks(Int_t k,Int_t max,AliBin *b,Int_t *idx,UInt_t *msk,Int_t& n) {
539   //------------------------------------------------------------
540   //find local maxima
541   //------------------------------------------------------------
542   if (n<31)
543   if (IsMaximum(k,max,b)) {
544     idx[n]=k; msk[n]=(2<<n);
545     n++;
546   }
547   b[k].SetMask(0);
548   if (b[k-max].GetMask()&1) FindPeaks(k-max,max,b,idx,msk,n);
549   if (b[k-1  ].GetMask()&1) FindPeaks(k-1  ,max,b,idx,msk,n);
550   if (b[k+max].GetMask()&1) FindPeaks(k+max,max,b,idx,msk,n);
551   if (b[k+1  ].GetMask()&1) FindPeaks(k+1  ,max,b,idx,msk,n);
552 }
553
554 void AliITSclustererV2::
555 MarkPeak(Int_t k, Int_t max, AliBin *bins, UInt_t m) {
556   //------------------------------------------------------------
557   //mark this peak
558   //------------------------------------------------------------
559   UShort_t q=bins[k].GetQ();
560
561   bins[k].SetMask(bins[k].GetMask()|m); 
562
563   if (bins[k-max].GetQ() <= q)
564      if ((bins[k-max].GetMask()&m) == 0) MarkPeak(k-max,max,bins,m);
565   if (bins[k-1  ].GetQ() <= q)
566      if ((bins[k-1  ].GetMask()&m) == 0) MarkPeak(k-1  ,max,bins,m);
567   if (bins[k+max].GetQ() <= q)
568      if ((bins[k+max].GetMask()&m) == 0) MarkPeak(k+max,max,bins,m);
569   if (bins[k+1  ].GetQ() <= q)
570      if ((bins[k+1  ].GetMask()&m) == 0) MarkPeak(k+1  ,max,bins,m);
571 }
572
573 void AliITSclustererV2::
574 MakeCluster(Int_t k,Int_t max,AliBin *bins,UInt_t m,AliITSclusterV2 &c) {
575   //------------------------------------------------------------
576   //make cluster using digits of this peak
577   //------------------------------------------------------------
578   Float_t q=(Float_t)bins[k].GetQ();
579   Int_t i=k/max, j=k-i*max;
580
581   c.SetQ(c.GetQ()+q);
582   c.SetY(c.GetY()+i*q); 
583   c.SetZ(c.GetZ()+j*q); 
584   c.SetSigmaY2(c.GetSigmaY2()+i*i*q);
585   c.SetSigmaZ2(c.GetSigmaZ2()+j*j*q);
586
587   bins[k].SetMask(0xFFFFFFFE);
588   
589   if (bins[k-max].GetMask() == m) MakeCluster(k-max,max,bins,m,c);
590   if (bins[k-1  ].GetMask() == m) MakeCluster(k-1  ,max,bins,m,c);
591   if (bins[k+max].GetMask() == m) MakeCluster(k+max,max,bins,m,c);
592   if (bins[k+1  ].GetMask() == m) MakeCluster(k+1  ,max,bins,m,c);
593 }
594
595 void AliITSclustererV2::
596 FindClustersSDD(AliBin* bins[2], Int_t nMaxBin, Int_t nzBins, 
597                 const TClonesArray *digits, TClonesArray *clusters) {
598   //------------------------------------------------------------
599   // Actual SDD cluster finder
600   //------------------------------------------------------------
601   Int_t ncl=0; TClonesArray &cl=*clusters;
602   for (Int_t s=0; s<2; s++)
603     for (Int_t i=0; i<nMaxBin; i++) {
604       if (bins[s][i].IsUsed()) continue;
605       Int_t idx[32]; UInt_t msk[32]; Int_t npeaks=0;
606       FindPeaks(i, nzBins, bins[s], idx, msk, npeaks);
607
608       if (npeaks>30) continue;
609
610       Int_t k,l;
611       for (k=0; k<npeaks-1; k++){//mark adjacent peaks
612         if (idx[k] < 0) continue; //this peak is already removed
613         for (l=k+1; l<npeaks; l++) {
614            if (idx[l] < 0) continue; //this peak is already removed
615            Int_t ki=idx[k]/nzBins, kj=idx[k] - ki*nzBins;
616            Int_t li=idx[l]/nzBins, lj=idx[l] - li*nzBins;
617            Int_t di=TMath::Abs(ki - li);
618            Int_t dj=TMath::Abs(kj - lj);
619            if (di>1 || dj>1) continue;
620            if (bins[s][idx[k]].GetQ() > bins[s][idx[l]].GetQ()) {
621               msk[l]=msk[k];
622               idx[l]*=-1;
623            } else {
624               msk[k]=msk[l];
625               idx[k]*=-1;
626               break;
627            } 
628         }
629       }
630
631       for (k=0; k<npeaks; k++) {
632         MarkPeak(TMath::Abs(idx[k]), nzBins, bins[s], msk[k]);
633       }
634         
635       for (k=0; k<npeaks; k++) {
636          if (idx[k] < 0) continue; //removed peak
637          AliITSclusterV2 c;
638          MakeCluster(idx[k], nzBins, bins[s], msk[k], c);
639
640          //if (c.GetQ() < 200) continue; //noise cluster
641
642          /*
643          Float_t s2 = c.GetSigmaY2()/c.GetQ() - c.GetY()*c.GetY();
644          Float_t w=par->GetPadPitchWidth(sec);
645          c.SetSigmaY2((s2 + 1./12.)*w*w);
646          if (s2 != 0.) {
647            c.SetSigmaY2(c.GetSigmaY2()*0.108);
648            if (sec<par->GetNInnerSector()) c.SetSigmaY2(c.GetSigmaY2()*2.07);
649          }
650
651          s2 = c.GetSigmaZ2()/c.GetQ() - c.GetZ()*c.GetZ();
652          w=par->GetZWidth();
653          c.SetSigmaZ2((s2 + 1./12.)*w*w);
654          if (s2 != 0.) {
655            c.SetSigmaZ2(c.GetSigmaZ2()*0.169);
656            if (sec<par->GetNInnerSector()) c.SetSigmaZ2(c.GetSigmaZ2()*1.77);
657          }
658          */
659
660          c.SetSigmaY2(0.0030*0.0030);
661          c.SetSigmaZ2(0.0020*0.0020);
662          c.SetDetectorIndex(fNdet[fI]);
663
664          Float_t y=c.GetY(),z=c.GetZ(), q=c.GetQ();
665          y/=q; z/=q;
666
667          y=(y-0.5)*fYpitchSDD;
668          y-=fHwSDD;
669          y-=fYoffSDD;  //delay ?
670          if (s) y=-y;
671
672          z=(z-0.5)*fZpitchSDD;
673          z-=fHlSDD;
674
675          y=-y-fYshift[fI];
676          z= z+fZshift[fI];
677          c.SetY(y);
678          c.SetZ(z);
679
680          c.SetQ(q/20.);  //to be consistent with the SSD charges
681
682          if (digits) {
683            AliBin *b=&bins[s][idx[k]];
684            AliITSdigitSDD* d=(AliITSdigitSDD*)digits->UncheckedAt(b->GetIndex());
685            Int_t l0=(d->GetTracks())[0];
686            if (l0<0) {
687              b=&bins[s][idx[k]-1];
688              if (b->GetQ()>0) {
689                d=(AliITSdigitSDD*)digits->UncheckedAt(b->GetIndex());
690                l0=(d->GetTracks())[0];
691              }
692            }
693            if (l0<0) {
694              b=&bins[s][idx[k]+1];
695              if (b->GetQ()>0) {
696                d=(AliITSdigitSDD*)digits->UncheckedAt(b->GetIndex());
697                l0=(d->GetTracks())[0];
698              }
699            }
700            if (l0<0) {
701              b=&bins[s][idx[k]-nzBins];
702              if (b->GetQ()>0) {
703                d=(AliITSdigitSDD*)digits->UncheckedAt(b->GetIndex());
704                l0=(d->GetTracks())[0];
705              }
706            }
707            if (l0<0) {
708              b=&bins[s][idx[k]+nzBins];
709              if (b->GetQ()>0) {
710                d=(AliITSdigitSDD*)digits->UncheckedAt(b->GetIndex());
711                l0=(d->GetTracks())[0];
712              }
713            }
714
715            if (l0<0) {
716              b=&bins[s][idx[k]+nzBins+1];
717              if (b->GetQ()>0) {
718                d=(AliITSdigitSDD*)digits->UncheckedAt(b->GetIndex());
719                l0=(d->GetTracks())[0];
720              }
721            }
722            if (l0<0) {
723              b=&bins[s][idx[k]+nzBins-1];
724              if (b->GetQ()>0) {
725                d=(AliITSdigitSDD*)digits->UncheckedAt(b->GetIndex());
726                l0=(d->GetTracks())[0];
727              }
728            }
729            if (l0<0) {
730              b=&bins[s][idx[k]-nzBins+1];
731              if (b->GetQ()>0) {
732                d=(AliITSdigitSDD*)digits->UncheckedAt(b->GetIndex());
733                l0=(d->GetTracks())[0];
734              }
735            }
736            if (l0<0) {
737              b=&bins[s][idx[k]-nzBins-1];
738              if (b->GetQ()>0) {
739                d=(AliITSdigitSDD*)digits->UncheckedAt(b->GetIndex());
740                l0=(d->GetTracks())[0];
741              }
742            }
743
744            {
745              Int_t lab[3];
746              lab[0]=(d->GetTracks())[0];
747              lab[1]=(d->GetTracks())[1];
748              lab[2]=(d->GetTracks())[2];
749              //CheckLabels(lab);
750              c.SetLabel(lab[0],0);
751              c.SetLabel(lab[1],1);
752              c.SetLabel(lab[2],2);
753            }
754          }
755
756          new (cl[ncl]) AliITSclusterV2(c); ncl++;
757       }
758     }
759 }
760
761 void AliITSclustererV2::
762 FindClustersSDD(const TClonesArray *digits, TClonesArray *clusters) {
763   //------------------------------------------------------------
764   // Actual SDD cluster finder
765   //------------------------------------------------------------
766   Int_t kNzBins = fNzSDD + 2;
767   const Int_t kMAXBIN=kNzBins*(fNySDD+2);
768
769   AliBin *bins[2];
770   bins[0]=new AliBin[kMAXBIN];
771   bins[1]=new AliBin[kMAXBIN];
772
773   AliITSdigitSDD *d=0;
774   Int_t i, ndigits=digits->GetEntriesFast();
775   for (i=0; i<ndigits; i++) {
776      d=(AliITSdigitSDD*)digits->UncheckedAt(i);
777      Int_t y=d->GetCoord2()+1;   //y
778      Int_t z=d->GetCoord1()+1;   //z
779      Int_t q=d->GetSignal();
780      if (z <= fNzSDD) {
781        bins[0][y*kNzBins+z].SetQ(q);
782        bins[0][y*kNzBins+z].SetMask(1);
783        bins[0][y*kNzBins+z].SetIndex(i);
784      } else {
785        z-=fNzSDD; 
786        bins[1][y*kNzBins+z].SetQ(q);
787        bins[1][y*kNzBins+z].SetMask(1);
788        bins[1][y*kNzBins+z].SetIndex(i);
789      }
790   }
791   
792   FindClustersSDD(bins, kMAXBIN, kNzBins, digits, clusters);
793
794   delete[] bins[0];
795   delete[] bins[1];
796 }
797
798 void AliITSclustererV2::FindClustersSDD(AliITSRawStream* input, 
799                                         TClonesArray** clusters) 
800 {
801   //------------------------------------------------------------
802   // Actual SDD cluster finder for raw data
803   //------------------------------------------------------------
804   Int_t nClustersSDD = 0;
805   Int_t kNzBins = fNzSDD + 2;
806   Int_t kMaxBin = kNzBins * (fNySDD+2);
807   AliBin* bins[2] = {NULL, NULL};
808
809   // read raw data input stream
810   while (kTRUE) {
811     Bool_t next = input->Next();
812     if (!next || input->IsNewModule()) {
813       Int_t iModule = input->GetPrevModuleID();
814
815       // when all data from a module was read, search for clusters
816       if (bins[0]) { 
817         clusters[iModule] = new TClonesArray("AliITSclusterV2");
818         fI = iModule;
819         FindClustersSDD(bins, kMaxBin, kNzBins, NULL, clusters[iModule]);
820         Int_t nClusters = clusters[iModule]->GetEntriesFast();
821         nClustersSDD += nClusters;
822         delete[] bins[0];
823         delete[] bins[1];
824       }
825
826       if (!next) break;
827       bins[0] = new AliBin[kMaxBin];
828       bins[1] = new AliBin[kMaxBin];
829     }
830
831     // fill the current digit into the bins array
832     Int_t iz = input->GetCoord1()+1;
833     Int_t side = ((iz <= fNzSDD) ? 0 : 1);
834     iz -= side*fNzSDD;
835     Int_t index = (input->GetCoord2()+1) * kNzBins + iz;
836     bins[side][index].SetQ(input->GetSignal());
837     bins[side][index].SetMask(1);
838     bins[side][index].SetIndex(index);
839   }
840
841   Info("FindClustersSDD", "found clusters in ITS SDD: %d", nClustersSDD);
842 }
843
844 void AliITSclustererV2::
845 FindClustersSSD(Ali1Dcluster* neg, Int_t nn, 
846                 Ali1Dcluster* pos, Int_t np,
847                 TClonesArray *clusters) {
848   //------------------------------------------------------------
849   // Actual SSD cluster finder
850   //------------------------------------------------------------
851   TClonesArray &cl=*clusters;
852
853   Int_t lab[4]={-2,-2,-2,-2};
854   Float_t tanp=fTanP, tann=fTanN;
855   if (fI>fLastSSD1) {tann=fTanP; tanp=fTanN;}
856
857   Int_t idet=fNdet[fI];
858   Int_t ncl=0;
859   for (Int_t i=0; i<np; i++) {
860     //Float_t dq_min=1.e+33;
861     Float_t ybest=1000,zbest=1000,qbest=0;
862     Float_t yp=pos[i].GetY()*fYpitchSSD; 
863     for (Int_t j=0; j<nn; j++) {
864       //if (pos[i].fTracks[0] != neg[j].fTracks[0]) continue;
865       Float_t yn=neg[j].GetY()*fYpitchSSD;
866       Float_t zt=(2*fHlSSD*tanp + yp - yn)/(tann+tanp);
867       Float_t yt=yn + tann*zt;
868       zt-=fHlSSD; yt-=fHwSSD;
869       if (TMath::Abs(yt)<fHwSSD+0.01)
870       if (TMath::Abs(zt)<fHlSSD+0.01) {
871       //if (TMath::Abs(pos[i].GetQ()-neg[j].GetQ())<dq_min) {
872         //dq_min=TMath::Abs(pos[i].GetQ()-neg[j].GetQ());
873         ybest=yt; zbest=zt; 
874         qbest=0.5*(pos[i].GetQ()+neg[j].GetQ());
875
876         lab[0]=pos[i].GetLabel(0);
877         lab[1]=pos[i].GetLabel(1);
878         lab[2]=neg[i].GetLabel(0);
879         lab[3]=(((i<<10) + j)<<10) + idet; // pos|neg|det
880         Float_t lp[5];
881         lp[0]=-ybest-fYshift[fI];
882         lp[1]= zbest+fZshift[fI];
883         lp[2]=0.0025*0.0025;  //SigmaY2
884         lp[3]=0.110*0.110;  //SigmaZ2
885         if (pos[i].GetNd()+neg[j].GetNd() > 4) {
886            lp[2]*=9;
887            lp[3]*=9;
888         }
889         lp[4]=qbest;        //Q
890
891         //CheckLabels(lab);
892         new (cl[ncl]) AliITSclusterV2(lab,lp); ncl++;
893       }
894     }
895     /*
896     if (ybest<100) {
897        lab[3]=idet;
898        Float_t lp[5];
899        lp[0]=-ybest-fYshift[fI];
900        lp[1]= zbest+fZshift[fI];
901        lp[2]=0.002*0.002;  //SigmaY2
902        lp[3]=0.080*0.080;  //SigmaZ2
903        lp[4]=qbest;        //Q
904        //
905        new (cl[ncl]) AliITSclusterV2(lab,lp); ncl++;
906     }
907     */
908   }
909 }
910
911 void AliITSclustererV2::
912 FindClustersSSD(const TClonesArray *digits, TClonesArray *clusters) {
913   //------------------------------------------------------------
914   // Actual SSD cluster finder
915   //------------------------------------------------------------
916   Int_t smax=digits->GetEntriesFast();
917   if (smax==0) return;
918
919   const Int_t MAX=1000;
920   Int_t np=0, nn=0; 
921   Ali1Dcluster pos[MAX], neg[MAX];
922   Float_t y=0., q=0., qmax=0.; 
923   Int_t lab[4]={-2,-2,-2,-2};
924
925   AliITSdigitSSD *d=(AliITSdigitSSD*)digits->UncheckedAt(0);
926   q += d->GetSignal();
927   y += d->GetCoord2()*d->GetSignal();
928   qmax=d->GetSignal();
929   lab[0]=d->GetTrack(0); lab[1]=d->GetTrack(1); lab[2]=d->GetTrack(2);
930   Int_t curr=d->GetCoord2();
931   Int_t flag=d->GetCoord1();
932   Int_t *n=&nn;
933   Ali1Dcluster *c=neg;
934   Int_t nd=1;
935   for (Int_t s=1; s<smax; s++) {
936       d=(AliITSdigitSSD*)digits->UncheckedAt(s);
937       Int_t strip=d->GetCoord2();
938       if ((strip-curr) > 1 || flag!=d->GetCoord1()) {
939          c[*n].SetY(y/q);
940          c[*n].SetQ(q);
941          c[*n].SetNd(nd);
942          c[*n].SetLabels(lab);
943          //Split suspiciously big cluster
944          if (nd>3) {
945             c[*n].SetY(y/q-0.5*nd);
946             c[*n].SetQ(0.5*q);
947             (*n)++;
948             if (*n==MAX) {
949               Error("FindClustersSSD","Too many 1D clusters !");
950               return;
951             }
952             c[*n].SetY(y/q+0.5*nd);
953             c[*n].SetQ(0.5*q);
954             c[*n].SetNd(nd);
955             c[*n].SetLabels(lab);
956          }
957          (*n)++;
958          if (*n==MAX) {
959           Error("FindClustersSSD","Too many 1D clusters !");
960           return;
961          }
962          y=q=qmax=0.;
963          nd=0;
964          lab[0]=lab[1]=lab[2]=-2;
965          if (flag!=d->GetCoord1()) { n=&np; c=pos; }
966       }
967       flag=d->GetCoord1();
968       q += d->GetSignal();
969       y += d->GetCoord2()*d->GetSignal();
970       nd++;
971       if (d->GetSignal()>qmax) {
972          qmax=d->GetSignal();
973          lab[0]=d->GetTrack(0); lab[1]=d->GetTrack(1); lab[2]=d->GetTrack(2);
974       }
975       curr=strip;
976   }
977   c[*n].SetY(y/q);
978   c[*n].SetQ(q);
979   c[*n].SetNd(nd);
980   c[*n].SetLabels(lab);
981   //Split suspiciously big cluster
982   if (nd>3) {
983      c[*n].SetY(y/q-0.5*nd);
984      c[*n].SetQ(0.5*q);
985      (*n)++;
986      if (*n==MAX) {
987         Error("FindClustersSSD","Too many 1D clusters !");
988         return;
989      }
990      c[*n].SetY(y/q+0.5*nd);
991      c[*n].SetQ(0.5*q);
992      c[*n].SetNd(nd);
993      c[*n].SetLabels(lab);
994   }
995   (*n)++;
996   if (*n==MAX) {
997      Error("FindClustersSSD","Too many 1D clusters !");
998      return;
999   }
1000
1001   FindClustersSSD(neg, nn, pos, np, clusters);
1002 }
1003
1004 void AliITSclustererV2::FindClustersSSD(AliITSRawStream* input, 
1005                                         TClonesArray** clusters) 
1006 {
1007   //------------------------------------------------------------
1008   // Actual SSD cluster finder for raw data
1009   //------------------------------------------------------------
1010   Int_t nClustersSSD = 0;
1011   const Int_t MAX = 1000;
1012   Ali1Dcluster clusters1D[2][MAX];
1013   Int_t nClusters[2] = {0, 0};
1014   Float_t q = 0.;
1015   Float_t y = 0.;
1016   Int_t nDigits = 0;
1017   Int_t prevStrip = -1;
1018   Int_t prevFlag = -1;
1019
1020   // read raw data input stream
1021   while (kTRUE) {
1022     Bool_t next = input->Next();
1023
1024     // check if a new cluster starts
1025     Int_t strip = input->GetCoord2();
1026     Int_t flag = input->GetCoord1();
1027     if ((!next || input->IsNewModule() ||
1028          (strip-prevStrip > 1) || (flag != prevFlag)) &&
1029         (nDigits > 0)) {
1030       if (nClusters[prevFlag] == MAX) {
1031         Error("FindClustersSSD", "Too many 1D clusters !");
1032         return;
1033       }
1034       Ali1Dcluster& cluster = clusters1D[prevFlag][nClusters[prevFlag]++];
1035       cluster.SetY(y/q);
1036       cluster.SetQ(q);
1037       cluster.SetNd(nDigits);
1038
1039       //Split suspiciously big cluster
1040       if (nDigits > 3) {
1041         cluster.SetY(y/q - 0.5*nDigits);
1042         cluster.SetQ(0.5*q);
1043         if (nClusters[prevFlag] == MAX) {
1044           Error("FindClustersSSD", "Too many 1D clusters !");
1045           return;
1046         }
1047         Ali1Dcluster& cluster2 = clusters1D[prevFlag][nClusters[prevFlag]++];
1048         cluster2.SetY(y/q + 0.5*nDigits);
1049         cluster2.SetQ(0.5*q);
1050         cluster2.SetNd(nDigits);
1051       }
1052       y = q = 0.;
1053       nDigits = 0;
1054     }
1055
1056     if (!next || input->IsNewModule()) {
1057       Int_t iModule = input->GetPrevModuleID();
1058
1059       // when all data from a module was read, search for clusters
1060       if (prevFlag >= 0) {
1061         clusters[iModule] = new TClonesArray("AliITSclusterV2");
1062         fI = iModule;
1063         FindClustersSSD(&clusters1D[0][0], nClusters[0], 
1064                         &clusters1D[1][0], nClusters[1], clusters[iModule]);
1065         Int_t nClusters = clusters[iModule]->GetEntriesFast();
1066         nClustersSSD += nClusters;
1067       }
1068
1069       if (!next) break;
1070       nClusters[0] = nClusters[1] = 0;
1071       y = q = 0.;
1072       nDigits = 0;
1073     }
1074
1075     // add digit to current cluster
1076     q += input->GetSignal();
1077     y += strip * input->GetSignal();
1078     nDigits++;
1079     prevStrip = strip;
1080     prevFlag = flag;
1081   }
1082
1083   Info("FindClustersSSD", "found clusters in ITS SSD: %d", nClustersSSD);
1084 }
1085
1086 #else   //V1
1087
1088 #include "AliITSDetType.h"
1089
1090 #include "AliITSsegmentationSPD.h"
1091 #include "AliITSClusterFinderSPD.h"
1092
1093 #include "AliITSresponseSDD.h"
1094 #include "AliITSsegmentationSDD.h"
1095 #include "AliITSClusterFinderSDD.h"
1096
1097 #include "AliITSsegmentationSSD.h"
1098 #include "AliITSClusterFinderSSD.h"
1099
1100
1101 void AliITSclustererV2::
1102 FindClustersSPD(const TClonesArray *digits, TClonesArray *clusters) {
1103   //------------------------------------------------------------
1104   // Actual SPD cluster finding based on AliITSClusterFinderSPD
1105   //------------------------------------------------------------
1106   static AliITS *its=(AliITS*)gAlice->GetModule("ITS");
1107   static TClonesArray *points=its->RecPoints();
1108   static AliITSsegmentationSPD *seg=
1109          (AliITSsegmentationSPD *)its->DetType(0)->GetSegmentationModel();
1110   static AliITSClusterFinderSPD cf(seg, (TClonesArray*)digits, points);
1111
1112   cf.FindRawClusters(fI);
1113   RecPoints2Clusters(points, fI, clusters);
1114   its->ResetRecPoints();
1115
1116 }
1117
1118 void AliITSclustererV2::
1119 FindClustersSDD(const TClonesArray *digits, TClonesArray *clusters) {
1120   //------------------------------------------------------------
1121   // Actual SDD cluster finding based on AliITSClusterFinderSDD
1122   //------------------------------------------------------------
1123   static AliITS *its=(AliITS*)gAlice->GetModule("ITS");
1124   static TClonesArray *points=its->RecPoints();
1125   static AliITSresponseSDD *resp=
1126         (AliITSresponseSDD *)its->DetType(1)->GetResponseModel();
1127   static AliITSsegmentationSDD *seg=
1128          (AliITSsegmentationSDD *)its->DetType(1)->GetSegmentationModel();
1129   static AliITSClusterFinderSDD 
1130          cf(seg,resp,(TClonesArray*)digits,its->ClustersAddress(1));
1131
1132   cf.FindRawClusters(fI);
1133   Int_t nc=points->GetEntriesFast();
1134   while (nc--) { //To be consistent with the SSD cluster charges
1135      AliITSRecPoint *p=(AliITSRecPoint*)points->UncheckedAt(nc);
1136      p->SetQ(p->GetQ/12.);
1137   }
1138   RecPoints2Clusters(points, fI, clusters);
1139   its->ResetClusters(1);
1140   its->ResetRecPoints();
1141
1142 }
1143
1144 void AliITSclustererV2::
1145 FindClustersSSD(const TClonesArray *digits, TClonesArray *clusters) {
1146   //------------------------------------------------------------
1147   // Actual SSD cluster finding based on AliITSClusterFinderSSD
1148   //------------------------------------------------------------
1149   static AliITS *its=(AliITS*)gAlice->GetModule("ITS");
1150   static TClonesArray *points=its->RecPoints();
1151   static AliITSsegmentationSSD *seg=
1152          (AliITSsegmentationSSD *)its->DetType(2)->GetSegmentationModel();
1153   static AliITSClusterFinderSSD cf(seg,(TClonesArray*)digits);
1154
1155   cf.FindRawClusters(fI);
1156   RecPoints2Clusters(points, fI, clusters);
1157   its->ResetRecPoints();
1158
1159 }
1160
1161 #endif