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