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