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