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