]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/UPGRADE/AliITSUpgradeClusterFinder.cxx
- fixing warnings/coverity
[u/mrichter/AliRoot.git] / ITS / UPGRADE / AliITSUpgradeClusterFinder.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 ////////////////////////////////////////////////////////////////////////////
19 //                                                                        //
20 // Base Class used to find                                                //
21 // the reconstructed points for ITS Upgrade                               // 
22 //                                                                        //
23 ////////////////////////////////////////////////////////////////////////////
24
25 #include "AliITSUpgradeClusterFinder.h"
26 #include "AliITSsegmentationUpgrade.h"
27 #include "AliITSRecPointU.h"
28 #include "AliITSDigitUpgrade.h"
29 #include "AliITSRawStreamSPD.h"
30 #include "AliLog.h"
31 #include <string.h>
32 #include <TObjString.h>
33 #include <TClonesArray.h>
34 #include <TTree.h>
35 #include <TMath.h>
36
37 AliITSUpgradeClusterFinder::AliITSUpgradeClusterFinder() : 
38   TObject(),
39   fNhitsLeft(0),
40   fOldModule(-1),
41   fClusterTypeFlag(kTRUE),
42   fFindClusterType(kFALSE),
43   fClusterTypeOrigCol(0),
44   fClusterTypeOrigRow(0),
45   fColSum(0),fRowSum(0),
46   fCharge(0),
47   fClusterWidthMaxCol(0),
48   fClusterWidthMinCol(0),
49   fClusterWidthMaxRow(0),
50   fClusterWidthMinRow(0),
51   fChargeArray(0x0),
52   fRecPoints(0x0),
53   fNlayers()
54
55   //
56   // Default constructor
57   //
58   AliITSsegmentationUpgrade *s = new AliITSsegmentationUpgrade();
59   fNlayers = s->GetNLayers();
60   delete s;
61   
62   fChargeArray = new TObjArray();
63   fChargeArray->SetOwner(kTRUE);
64   fRecPoints = new TClonesArray("AliITSRecPointU",3000);
65   fTmpLabel[0]=-5;
66   fTmpLabel[1]=-5;
67   fTmpLabel[2]=-5;
68   for(int il=0; il<10;il++) fLabels[il]=-5;
69   for(int k=0; k<999999; k++){
70     fHitCol[k]=999;
71     fHitRow[k]=999;
72   }
73 }
74
75 //___________________________________________________________________________________
76 AliITSUpgradeClusterFinder::~AliITSUpgradeClusterFinder() {
77   if(fChargeArray) delete fChargeArray;
78   if(fRecPoints) delete fRecPoints; 
79 }
80 //___________________________________________________________________________________
81 void AliITSUpgradeClusterFinder::StartEvent() {
82   NewEvent();
83 }
84 //___________________________________________________________________________________
85 Int_t AliITSUpgradeClusterFinder::ProcessHit(Int_t layer , UInt_t col, UInt_t row, UShort_t charge, Int_t label[3]) {
86   //
87   // Adds one pixel to the cluster 
88   //
89   
90   if (layer>=fNlayers) {
91     printf("ERROR: UpgradeClusterFinder::ProcessHit: Out of bounds: layer ,col,row, charge, label(0,1,2)= %d,%d,%d,%d,(%d,%d,%d)\n",layer ,col,row,charge,label[0],label[1],label[2]); 
92     return 1;
93   }
94
95   // do we need to perform clustering on previous module?
96   if (fOldModule!=-1 && (Int_t)layer!=fOldModule) {
97     fChargeArray->AddLast(new TObjString(Form("%i %i %i %i %i %i",col,row,charge,label[0],label[1],label[2])));
98     DoModuleClustering(fOldModule,charge);
99     NewModule();
100   }
101   // add hit
102   fChargeArray->AddLast(new TObjString(Form("%i %i %i %i %i %i",col,row,charge,label[0],label[1],label[2])));
103
104   fOldModule=layer;
105   fHitCol[fNhitsLeft]=col;
106   fHitRow[fNhitsLeft]=row;
107   fHits[col][row]=kTRUE;
108   fTmpLabel[0]=label[0];
109   fTmpLabel[1]=label[1];
110   fTmpLabel[2]=label[2];
111   fNhitsLeft++;
112   return 0;
113 }
114 //___________________________________________________________________________________
115 void AliITSUpgradeClusterFinder::FinishEvent() {
116   if (fNhitsLeft>0) {
117     DoModuleClustering(fOldModule,fCharge);
118   }
119 }
120 //___________________________________________________________________________________
121 UInt_t AliITSUpgradeClusterFinder::GetClusterCount(Int_t layer) const {
122   //
123   // number of clusters in the layer
124   // 
125   if (layer>=fNlayers ) {
126     printf("ERROR: UpgradeClusterFinder::GetClusterCount: Module out of bounds: layer = %d\n",layer);
127     return 0;
128   }
129   return fClusterList[layer].GetNrEntries();
130 }
131 //___________________________________________________________________________________
132 Float_t AliITSUpgradeClusterFinder::GetClusterMeanCol(Int_t layer, UInt_t index) {
133   //
134   // cluster position in terms of colums : mean column ID
135   //
136   if (layer>=fNlayers) {
137     printf("ERROR: UpgradeClusterFinder::GetClusterMeanCol: Module out of bounds: layer = %d\n",layer);
138     return 0;
139   }
140   return fClusterList[layer].GetColIndex(index);
141 }
142 //___________________________________________________________________________________
143 Float_t AliITSUpgradeClusterFinder::GetClusterMeanRow(Int_t layer, UInt_t index) {
144   //
145   // cluster position in terms of rows : mean row ID
146   //
147   if (layer>=fNlayers) {
148     printf("ERROR: UpgradeClusterFinder::GetClusterMeanRow: Module out of bounds: layer = %d\n",layer);
149     return 0;
150   }
151   return fClusterList[layer].GetRowIndex(index);
152 }
153 //___________________________________________________________________________________
154 UInt_t AliITSUpgradeClusterFinder::GetClusterSize(Int_t layer, UInt_t index) {
155   //
156   // total number of pixels of the cluster 
157   //
158   if (layer>=fNlayers) {
159     printf("ERROR: UpgradeClusterFinder::GetClusterSize: Module out of bounds: layer = %d\n",layer);
160     return 0;
161   }
162   return fClusterList[layer].GetSizeIndex(index);
163 }
164 //___________________________________________________________________________________
165 UInt_t AliITSUpgradeClusterFinder::GetClusterWidthZ(Int_t layer, UInt_t index) {
166   //
167   // # pixels of the cluster in Z direction
168   //
169   
170   if (layer>=fNlayers) {
171     printf("ERROR: UpgradeClusterFinder::GetClusterWidthZ: Module out of bounds: layer = %d\n",layer);
172     return 0;
173   }
174   return fClusterList[layer].GetWidthZIndex(index);
175 }
176 //___________________________________________________________________________________
177 UInt_t AliITSUpgradeClusterFinder::GetClusterWidthPhi(Int_t layer, UInt_t index) {
178   //
179   // # pixels of the cluster in phi direction (XY plane)
180   //
181
182   if (layer>=fNlayers) {
183     printf("ERROR: UpgradeClusterFinder::GetClusterWidthPhi: Module out of bounds: layer = %d\n",layer);
184     return 0;
185   }
186   return fClusterList[layer].GetWidthPhiIndex(index);
187 }
188 //___________________________________________________________________________________
189 UInt_t AliITSUpgradeClusterFinder::GetClusterType(Int_t layer, UInt_t index) {
190   //
191   // cluster shape
192   //
193
194   if (layer>=fNlayers) {
195     printf("ERROR: UpgradeClusterFinder::GetClusterType: Module out of bounds: layer = %d\n",layer);
196     return 0;
197   }
198   return fClusterList[layer].GetTypeIndex(index);
199 }
200 //___________________________________________________________________________________
201 void AliITSUpgradeClusterFinder::PrintAllClusters() {
202   //
203   // printout of the cluster information
204   // 
205   
206   for (Int_t layer=0; layer<fNlayers; layer++) {
207     PrintClusters(layer);
208   }
209 }
210 //___________________________________________________________________________________
211 void AliITSUpgradeClusterFinder::PrintClusters(Int_t layer) {
212   //
213   // printout of the cluster information
214   //   
215   
216   if (layer>=fNlayers) {
217     printf("ERROR: UpgradeClusterFinder::PrintClusters: Out of bounds: layer = %d\n",layer);
218     return;
219   }
220   if(fClusterList[layer].GetNrEntries()==0) {
221     printf("no cluster list entries. Exiting... \n");
222     return;
223   }
224   for (UInt_t c=0; c<fClusterList[layer].GetNrEntries(); c++) {
225     printf("layer  %d , z,y=%f,%f , size=%d , type=%d labels=%d %d %d (label printout to be implemented...)\n",layer,fClusterList[layer].GetColIndex(c),fClusterList[layer].GetRowIndex(c),fClusterList[layer].GetSizeIndex(c),fClusterList[layer].GetTypeIndex(c),999,999,999);  
226   }  
227 }
228 //___________________________________________________________________________________
229 void AliITSUpgradeClusterFinder::NewEvent() {
230   //
231   // Cleaning up and preparation for the clustering procedure
232   //
233   
234   for (Int_t i=0; i<fNlayers; i++) {
235     fClusterList[i].Clear();
236   }
237   NewModule();
238   fOldModule = -1;
239 }
240 //___________________________________________________________________________________
241 void AliITSUpgradeClusterFinder::NewModule() {
242   //
243   // Initializations
244   //
245   fNhitsLeft=0;
246   memset(fHits,0,999*999*sizeof(Bool_t));
247   fChargeArray->Clear();
248 }
249 //___________________________________________________________________________________
250 Int_t AliITSUpgradeClusterFinder::DoModuleClustering(Int_t Layer, UShort_t charge) {
251   //
252   // Clustering and cluster-list container filling
253   //
254   UInt_t maxHits=fNhitsLeft;
255   for (UInt_t hit=0; hit<maxHits; hit++) {
256     if (fClusterTypeFlag) fFindClusterType = kTRUE;
257     
258     fClusterWidthMinCol = fHitCol[hit];
259     fClusterWidthMinRow = fHitRow[hit];
260     fClusterWidthMaxCol = fHitCol[hit];
261     fClusterWidthMaxRow = fHitRow[hit];
262  
263     fClusterTypeOrigCol = fHitCol[hit];
264     fClusterTypeOrigRow = fHitRow[hit];
265     memset(fClusterTypeArea,0,kMAXCLUSTERTYPESIDEZ*kMAXCLUSTERTYPESIDEY*sizeof(Bool_t));
266     fColSum=0;
267     fRowSum=0;
268     UInt_t size = FindClusterRecu(fClusterTypeOrigCol,fClusterTypeOrigRow,charge);
269     if(size==1){
270       fCharge=GetPixelCharge(fColSum,fRowSum);
271       AddLabelIndex(fColSum,fRowSum);
272     }
273     if (size>0) {
274       if(size>1) AliDebug(2,Form("DoModuleClustering, size %i , labels :  %i  %i  %i \n",size,fLabels[0],fLabels[1],fLabels[2]));
275       fClusterList[Layer].Insert((Float_t)fColSum/size, (Float_t)fRowSum/size, size, GetClusterWidthZ(), GetClusterWidthPhi(), GetClusterType(size), fCharge,fLabels);
276       fCharge=0;
277       for(Int_t i=0; i<10; i++) fLabels[i]=-5;
278     }
279     if (fNhitsLeft==0) break;
280   }
281   return 0;
282 }
283 //___________________________________________________________________________________
284 UInt_t AliITSUpgradeClusterFinder::FindClusterRecu(Int_t col, Int_t row, UShort_t charge) {
285   //
286   // Pixel selection for the clusters (adjiacent pixels or pixels at the corners) 
287   //
288   if (col<0 || !fHits[col][row]) {
289     return 0;
290   }
291   fHits[col][row]=kFALSE;
292   fColSum+=col;
293   fRowSum+=row;
294   fNhitsLeft--;
295   UInt_t size=1;
296
297   if (col>fClusterWidthMaxCol) fClusterWidthMaxCol = col;
298   if (col<fClusterWidthMinCol) fClusterWidthMinCol = col;
299   if (row>fClusterWidthMaxRow) fClusterWidthMaxRow = row;
300   if (row<fClusterWidthMaxRow) fClusterWidthMinRow = row;
301
302   if (fFindClusterType) {
303     Short_t diffz = fClusterTypeOrigCol - col;
304     Short_t diffy = row - fClusterTypeOrigRow;
305     if (diffz>=kMAXCLUSTERTYPESIDEZ || diffy>=kMAXCLUSTERTYPESIDEY) {
306       fFindClusterType=kFALSE;
307     }
308     else {
309       if (diffz==-1) {
310         ShiftClusterTypeArea(kSHIFTRIGHT);
311         diffz=0;
312       }
313       if (diffy==-1) {
314         ShiftClusterTypeArea(kSHIFTDOWN);
315         diffy=0;
316       }
317       fClusterTypeArea[diffz][diffy] = kTRUE;
318     }
319   }
320   // straight:
321   size+=FindClusterRecu(col  ,row-1,charge);
322   fCharge+=GetPixelCharge(col,row-1);
323   AddLabelIndex(col,row-1);
324
325   size+=FindClusterRecu(col-1,row  ,charge);
326   fCharge+=GetPixelCharge(col-1,row);
327   AddLabelIndex(col-1,row);
328
329   size+=FindClusterRecu(col  ,row+1,charge);
330   fCharge+=GetPixelCharge(col,row+1);
331   AddLabelIndex(col,row+1);
332
333   size+=FindClusterRecu(col+1,row ,charge );
334   fCharge+=GetPixelCharge(col+1,row);
335   AddLabelIndex(col+1,row);
336
337   // diagonal:
338   size+=FindClusterRecu(col-1,row-1,charge);
339   fCharge+=GetPixelCharge(col-1,row-1);
340   AddLabelIndex(col-1,row-1);
341
342   size+=FindClusterRecu(col-1,row+1,charge);
343   fCharge+=GetPixelCharge(col-1,row+1);
344   AddLabelIndex(col-1,row+1);
345
346   size+=FindClusterRecu(col+1,row-1,charge);
347   fCharge+=GetPixelCharge(col+1,row-1);
348   AddLabelIndex(col+1,row-1);
349
350   size+=FindClusterRecu(col+1,row+1,charge);
351   fCharge+=GetPixelCharge(col+1,row+1);
352   AddLabelIndex(col+1,row+1);
353   return size;
354 }
355 //___________________________________________________________________________________
356 void AliITSUpgradeClusterFinder::ShiftClusterTypeArea(UInt_t direction) {
357   //
358   // Cluster checks 
359   //
360   
361   if (direction == kSHIFTRIGHT) {
362     fClusterTypeOrigCol++;
363     Bool_t tmpClusterTypeArea[kMAXCLUSTERTYPESIDEZ][kMAXCLUSTERTYPESIDEY];
364     memset(tmpClusterTypeArea,0,kMAXCLUSTERTYPESIDEZ*kMAXCLUSTERTYPESIDEY*sizeof(Bool_t));
365     for (UInt_t z=0; z<kMAXCLUSTERTYPESIDEZ; z++) {
366       for (UInt_t y=0; y<kMAXCLUSTERTYPESIDEY; y++) {
367         if (fClusterTypeArea[z][y]) {
368           if (z==kMAXCLUSTERTYPESIDEZ-1) {
369             fFindClusterType=kFALSE;
370             return;
371           }
372           else {
373             tmpClusterTypeArea[z+1][y] = kTRUE;
374           }
375         }
376       }
377     }
378     memcpy(fClusterTypeArea,tmpClusterTypeArea,kMAXCLUSTERTYPESIDEZ*kMAXCLUSTERTYPESIDEY*sizeof(Bool_t));
379   }
380   else if (direction == kSHIFTDOWN) {
381     fClusterTypeOrigRow--;
382     Bool_t tmpClusterTypeArea[kMAXCLUSTERTYPESIDEZ][kMAXCLUSTERTYPESIDEY];
383     memset(tmpClusterTypeArea,0,kMAXCLUSTERTYPESIDEZ*kMAXCLUSTERTYPESIDEY*sizeof(Bool_t));
384     for (UInt_t z=0; z<kMAXCLUSTERTYPESIDEZ; z++) {
385       for (UInt_t y=0; y<kMAXCLUSTERTYPESIDEY; y++) {
386         if (fClusterTypeArea[z][y]) {
387           if (y==kMAXCLUSTERTYPESIDEY-1) {
388             fFindClusterType=kFALSE;
389             return;
390           }
391           else {
392             tmpClusterTypeArea[z][y+1] = kTRUE;
393           }
394         }
395       }
396     }
397     memcpy(fClusterTypeArea,tmpClusterTypeArea,kMAXCLUSTERTYPESIDEZ*kMAXCLUSTERTYPESIDEY*sizeof(Bool_t));
398   }
399 }
400 //___________________________________________________________________________________
401 UShort_t AliITSUpgradeClusterFinder::GetCharge(Int_t layer,UInt_t index ) {
402   return fClusterList[layer].GetCharge(index);
403 }
404 //___________________________________________________________________________________
405 Int_t * AliITSUpgradeClusterFinder::GetLabels(Int_t layer,UInt_t index) {
406   return fClusterList[layer].GetLabels(index);
407 }
408
409 //___________________________________________________________________________________
410 UInt_t AliITSUpgradeClusterFinder::GetClusterWidthZ() {
411   return fClusterWidthMaxCol-fClusterWidthMinCol+1;
412 }
413 //___________________________________________________________________________________
414 UInt_t AliITSUpgradeClusterFinder::GetClusterWidthPhi() {
415   return fClusterWidthMaxRow-fClusterWidthMinRow+1;
416 }
417 //___________________________________________________________________________________
418 UInt_t AliITSUpgradeClusterFinder::GetClusterType(UInt_t size) {
419   //
420   // Cluster shape
421   //
422   if (!fFindClusterType || size>4) 
423     return 0;
424
425   // X
426   if (size==1) 
427     return 1;
428
429   // X
430   // X
431   if (size==2 && 
432       fClusterTypeArea[0][1] && 
433       fClusterTypeArea[0][0] ) 
434     return 2;
435
436   // XX
437   if (size==2 && 
438       fClusterTypeArea[1][0] &&
439       fClusterTypeArea[0][0] ) 
440     return 3;
441
442   // X
443   // X
444   // X
445   if (size==3 && 
446       fClusterTypeArea[0][2] && 
447       fClusterTypeArea[0][1] && 
448       fClusterTypeArea[0][0] ) 
449     return 4;
450   
451   // XX
452   // X      and equivalent...
453   if (size==3 && 
454       (
455        (
456         fClusterTypeArea[0][0] && 
457         fClusterTypeArea[0][1] && 
458         (fClusterTypeArea[1][1] || 
459          fClusterTypeArea[1][0])
460         )
461        ||
462        (
463         fClusterTypeArea[1][0] && 
464         fClusterTypeArea[1][1] && 
465         (fClusterTypeArea[0][1] || 
466          fClusterTypeArea[0][0])
467         )
468        )
469       ) 
470     return 5;
471   
472   // XXX
473   if (size==3 && 
474       fClusterTypeArea[2][0] && 
475       fClusterTypeArea[1][0] && 
476       fClusterTypeArea[0][0] )
477     return 6;
478   
479   // X
480   // X
481   // X
482   // X
483   if (size==4 &&
484       fClusterTypeArea[0][3] && 
485       fClusterTypeArea[0][2] &&
486       fClusterTypeArea[0][1] && 
487       fClusterTypeArea[0][0] ) 
488     return 7;
489   
490   // XX
491   // XX
492   if (size==4 && 
493       fClusterTypeArea[1][1] && 
494       fClusterTypeArea[1][0] && 
495       fClusterTypeArea[0][1] &&
496       fClusterTypeArea[0][0] ) 
497     return 8;
498
499   // XX
500   // X
501   // X     and equivalent...
502   if (size==4 &&
503       (
504        (
505         fClusterTypeArea[1][0] && 
506         fClusterTypeArea[1][1] && 
507         fClusterTypeArea[1][2] && 
508         fClusterTypeArea[0][0]
509         )
510        ||
511        (
512         fClusterTypeArea[1][0] && 
513         fClusterTypeArea[1][1] && 
514         fClusterTypeArea[1][2] && 
515         fClusterTypeArea[0][2]
516         )
517        ||
518        (
519         fClusterTypeArea[0][0] && 
520         fClusterTypeArea[0][1] && 
521         fClusterTypeArea[0][2] && 
522         fClusterTypeArea[1][0]
523         )
524        ||
525        (
526         fClusterTypeArea[0][0] && 
527         fClusterTypeArea[0][1] && 
528         fClusterTypeArea[0][2] && 
529         fClusterTypeArea[1][2]
530         )
531        )
532       )
533     return 9;
534
535   // X
536   // XX
537   // X     and equivalent...
538   if (size==4 &&
539       (
540        (
541         fClusterTypeArea[1][0] && 
542         fClusterTypeArea[1][1] && 
543         fClusterTypeArea[1][2] && 
544         fClusterTypeArea[0][1]
545         )
546        ||
547        (
548         fClusterTypeArea[0][0] && 
549         fClusterTypeArea[0][1] && 
550         fClusterTypeArea[0][2] && 
551         fClusterTypeArea[1][1]
552         )
553        )
554       )
555     return 10;
556
557   // X
558   // XXX    and equivalent...
559   if (size==4 &&
560       (
561        (
562         fClusterTypeArea[2][0] && 
563         fClusterTypeArea[2][1] && 
564         fClusterTypeArea[1][1] && 
565         fClusterTypeArea[0][1]
566         )
567        ||
568        (
569         fClusterTypeArea[0][0] && 
570         fClusterTypeArea[2][1] && 
571         fClusterTypeArea[1][1] && 
572         fClusterTypeArea[0][1]
573         )
574        ||
575        (
576         fClusterTypeArea[2][1] && 
577         fClusterTypeArea[2][0] && 
578         fClusterTypeArea[1][0] && 
579         fClusterTypeArea[0][0]
580         )
581        ||
582        (
583         fClusterTypeArea[0][1] && 
584         fClusterTypeArea[2][0] && 
585         fClusterTypeArea[1][0] && 
586         fClusterTypeArea[0][0]
587         )
588        )
589       )
590     return 11;
591
592   //  X
593   // XXX     and equivalent...
594   if (size==4 &&
595       (
596        (
597         fClusterTypeArea[1][0] && 
598         fClusterTypeArea[2][1] && 
599         fClusterTypeArea[1][1] && 
600         fClusterTypeArea[0][1]
601         )
602        ||
603        (
604         fClusterTypeArea[1][1] && 
605         fClusterTypeArea[2][0] && 
606         fClusterTypeArea[1][0] && 
607         fClusterTypeArea[0][0]
608         )
609        )
610       )
611     return 12;
612
613   //  X
614   // X      and equivalent...
615   if (size==2 &&
616       (
617        (
618         fClusterTypeArea[0][0] &&
619         fClusterTypeArea[1][1] 
620         )
621        ||
622        (
623         fClusterTypeArea[1][0] &&
624         fClusterTypeArea[0][1] 
625         )
626        )
627       )
628     return 13;
629
630   //  X
631   //  X
632   // X      and equivalent...
633   if (size==3 &&
634       (
635        (
636         fClusterTypeArea[0][0] &&
637         fClusterTypeArea[0][1] && 
638         fClusterTypeArea[1][2] 
639         )
640        ||
641        (
642         fClusterTypeArea[1][0] &&
643         fClusterTypeArea[1][1] && 
644         fClusterTypeArea[0][2] 
645         )
646        ||
647        (
648         fClusterTypeArea[1][0] &&
649         fClusterTypeArea[0][1] && 
650         fClusterTypeArea[0][2] 
651         )
652        ||
653        (
654         fClusterTypeArea[0][0] &&
655         fClusterTypeArea[1][1] && 
656         fClusterTypeArea[1][2] 
657         )
658        )
659       )
660     return 14;
661
662   //  XX
663   // X      and equivalent...
664   if (size==3 &&
665       (
666        (
667         fClusterTypeArea[0][0] &&
668         fClusterTypeArea[1][0] && 
669         fClusterTypeArea[2][1] 
670         )
671        ||
672        (
673         fClusterTypeArea[0][1] &&
674         fClusterTypeArea[1][1] && 
675         fClusterTypeArea[2][0] 
676         )
677        ||
678        (
679         fClusterTypeArea[0][0] &&
680         fClusterTypeArea[1][1] && 
681         fClusterTypeArea[2][1] 
682         )
683        ||
684        (
685         fClusterTypeArea[0][1] &&
686         fClusterTypeArea[1][0] && 
687         fClusterTypeArea[2][0] 
688         )
689        )
690       )
691     return 15;
692
693
694
695   return 0;
696 }
697
698 //___________________________________________________________________________________________________
699 UInt_t AliITSUpgradeClusterFinder::GetPixelCharge(UInt_t col, UInt_t row){
700   //
701   //...self explaining
702   //
703   Int_t q=0;
704   for(Int_t entry =0; entry < fChargeArray->GetEntries(); entry++){
705     TObjString *s = (TObjString*)fChargeArray->At(entry);
706     TString name = s->GetString();
707     if(!name.Contains(Form("%i %i",col,row)))
708       continue;
709     TObjArray *array = name.Tokenize(" ");
710     TString charge = ((TObjString*)array->At(2))->String();
711     TString rowS, colS;
712     rowS = ((TObjString*)array->At(0))->String();
713     colS = ((TObjString*)array->At(1))->String();
714     delete array;
715     q=charge.Atoi();
716     return q;
717
718   }
719   return q;
720 }
721 //____________________________________________________
722
723 void AliITSUpgradeClusterFinder::AddLabelIndex(UInt_t col, UInt_t row){
724   //
725   // Adding cluster labels
726   //
727
728   for(Int_t entry =0; entry < fChargeArray->GetEntries(); entry++){
729     TObjString *s = (TObjString*)fChargeArray->At(entry);
730     TString name = s->GetString();
731     if(!name.Contains(Form("%i %i",col,row)))
732       continue;
733     TObjArray *array = name.Tokenize(" ");
734     TString index[3];
735     Int_t label[3];
736     for(Int_t i=0; i<3; i++){
737       index[i]=((TObjString*)array->At(3+i))->String();
738       label[i]=index[i].Atoi();
739
740     }
741     SetLabels(label);
742     delete array;
743   }
744 }
745 //____________________________________________________
746
747 void AliITSUpgradeClusterFinder::SetLabels(Int_t label[3]){
748   
749   //Set the MC lables to the cluster
750
751   Bool_t isAssigned[3]={kFALSE,kFALSE,kFALSE};
752
753   for(Int_t i=0; i<10; i++){
754     for(Int_t k=0; k<3; k++){
755       if(fLabels[i]!=label[k] && label[k]>-1 && fLabels[i]<0) {
756         if(!isAssigned[k]) {
757           //  printf("assignign...\n.");
758           //  printf("Assignign  fLabels[%i]=%i  label[%i]=%i \n",i,fLabels[i],k,label[k]);
759           fLabels[i+k]=label[k];
760           isAssigned[k]=kTRUE;
761         }
762       }
763     }
764   }
765 }
766
767 //____________________________________________________
768 void AliITSUpgradeClusterFinder::MakeRecPointBranch(TTree *treeR){
769   //
770   // Creating the branch (see AliITSUpgradeReconstructor::Reconstruct)
771   //
772
773   if(!fRecPoints)fRecPoints = new TClonesArray("AliITSRecPointU",1000);
774   if (treeR) {
775     TBranch *branch = treeR->GetBranch("ITSRecPoints");
776     if (branch) return ;
777     else branch = treeR->Branch("ITSRecPoints",&fRecPoints); 
778   }
779 }
780 //____________________________________________________
781 void AliITSUpgradeClusterFinder::SetRecPointTreeAddress(TTree *treeR){
782   //
783   // Addressing the branch (see AliITSUpgradeReconstructor::Reconstruct)
784   //
785   if(!treeR) return;
786   if(!fRecPoints) fRecPoints = new TClonesArray("AliITSRecPointU",1000);
787
788   TBranch *branch;
789   branch = treeR->GetBranch("ITSRecPoints");
790   if (branch) {
791     branch->SetAddress(&fRecPoints);
792   } else AliError("No RecPoint branch available");
793
794 }
795 //____________________________________________________
796 void AliITSUpgradeClusterFinder::DigitsToRecPoints(const TObjArray *digList) {
797   //
798   // the clusterization is performed here
799   //
800   AliITSsegmentationUpgrade *segmentation = new AliITSsegmentationUpgrade(); 
801   AliITSRecPointU  recpnt;
802   Int_t nClusters =0;
803   TClonesArray &lrecp = *fRecPoints;
804
805   for(Int_t ilayer=0; ilayer < fNlayers ;ilayer ++){
806     TClonesArray *pArrDig= (TClonesArray*)digList->At(ilayer);
807     StartEvent();
808     AliDebug(1,Form("layer %i : # digits %i",ilayer,pArrDig->GetEntries()));
809     for(Int_t ientr =0; ientr < pArrDig->GetEntries() ; ientr++){
810       AliITSDigitUpgrade *dig = (AliITSDigitUpgrade*)pArrDig->At(ientr);
811       Int_t colz=dig->GetzPixelNumber();
812       Int_t rowx=dig->GetxPixelNumber();
813       Double_t hitcharge= (dig->GetNelectrons());
814
815       ProcessHit(ilayer,colz, rowx,(Short_t)hitcharge,dig->GetTracks());
816     }//ientr
817     FinishEvent();
818     
819     for(UInt_t nClu = 0; nClu <  GetClusterCount(ilayer); nClu++){
820       UShort_t charge = GetCharge(ilayer, nClu);
821       recpnt.SetQ(charge);
822       recpnt.SetLayer(ilayer);
823       Int_t *lab=GetLabels(ilayer,nClu);
824       for(Int_t l=0; l<3; l++) {recpnt.SetLabel(lab[l],l);}
825
826       Bool_t check2;
827       Double_t xcheck2=0.;
828       Double_t ycheck2=0.;
829       Double_t zcheck2=0.;
830       Double_t xzl2[2]={0.,0.};
831       Double_t xPixC2,zPixC2=0.;
832
833       xPixC2 = GetClusterMeanRow(ilayer, nClu);
834       zPixC2 = GetClusterMeanCol(ilayer, nClu);
835       xzl2[0] = xPixC2*(segmentation->GetCellSizeX(ilayer))+0.5*(segmentation-> GetCellSizeX(ilayer));
836       xzl2[1] = zPixC2*(segmentation->GetCellSizeZ(ilayer))+0.5*(segmentation->GetCellSizeZ(ilayer))-(segmentation->GetHalfLength(ilayer));
837       check2 = segmentation->DetToGlobal(ilayer,xzl2[0], xzl2[1],xcheck2,ycheck2,zcheck2);
838       recpnt.SetType(GetClusterType(ilayer,nClu ));
839       recpnt.SetLocalCoord(xzl2[0],xzl2[1]); //temporary solution (no LocalToTrack Matrix)
840       //from global to tracking system coordinate
841       // global coordinate -> local coordinate getting alpha angle of the recpoint
842       Float_t xclg = xcheck2;//upgrade clusters global coordinate ( ITS official: GetX tracking coordinate)
843       Float_t yclg = ycheck2;
844       Float_t zclg = zcheck2;
845       Double_t phiclu1rad, phiclu1deg;
846       phiclu1rad=TMath::ATan2(yclg,xclg);//cluster phi angle (rad)
847       if (phiclu1rad<0) phiclu1rad+=TMath::TwoPi();//from 0 to 360
848       else if (phiclu1rad>=TMath::TwoPi()) phiclu1rad-=TMath::TwoPi();//
849
850       phiclu1deg=180.*phiclu1rad/TMath::Pi();// in deg
851       Int_t ladder;// virtual segmentation starting from the cluster phi
852
853       ladder=(Int_t)phiclu1deg/18;// in which ladder the cluster is
854       Double_t alpha= (ladder*18.+9.)*TMath::Pi()/180.;//angle at the center of the ladder (rad)
855
856       //alpha rotation
857       Float_t xclu1Tr = xclg*TMath::Cos(alpha)-yclg*TMath::Sin(alpha);
858       Float_t yclu1 = yclg*TMath::Cos(alpha)+xclg*TMath::Sin(alpha);
859       Float_t xclu1 = TMath::Sqrt(xclu1Tr*xclu1Tr+yclu1*yclu1);
860       Float_t zclu1 = zclg;
861       Double_t phiTrk= (phiclu1rad-alpha);// cluster angle in the rotated system (rad)
862
863       yclu1=xclu1*phiTrk; // tracking system coordinate: r*phi
864       recpnt.SetX(0.);
865       recpnt.SetY(yclu1);
866       recpnt.SetZ(zclu1);
867
868       Double_t xsize, zsize;
869       segmentation->GetSegmentation(ilayer,xsize, zsize);
870       recpnt.SetSigmaY2(xsize/TMath::Sqrt(12)*xsize/TMath::Sqrt(12));
871       recpnt.SetSigmaZ2(zsize/TMath::Sqrt(12)*zsize/TMath::Sqrt(12));
872       new(lrecp[nClusters++]) AliITSRecPointU(recpnt);
873       //Int_t idx = fRecPoints->GetEntries();
874       AliDebug(1,Form("recpoint : Nelectrons %f (entry %i)",recpnt.GetQ(),fRecPoints->GetEntries()));
875     }//cluster list entries
876   }//ilayer
877   if(segmentation) delete segmentation;
878 }
879