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