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