]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/UPGRADE/AliITSUpgradeClusterFinder.cxx
359b97d6f7515992c76a7f3b2bd78f8fee9be228
[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 "AliITSRawStreamSPD.h"
27 #include "AliLog.h"
28 #include <string.h>
29 #include <TObjString.h>
30
31
32 AliITSUpgradeClusterFinder::AliITSUpgradeClusterFinder() : 
33   fNhitsLeft(0),
34   fOldModule(-1),
35   fClusterTypeFlag(kTRUE),
36   fFindClusterType(kFALSE),
37   fClusterTypeOrigCol(0),
38   fClusterTypeOrigRow(0),
39   fColSum(0),fRowSum(0),
40   fCharge(0),
41   fClusterWidthMaxCol(0),
42   fClusterWidthMinCol(0),
43   fClusterWidthMaxRow(0),
44   fClusterWidthMinRow(0),
45   fChargeArray(0x0)
46
47   fChargeArray = new TObjArray();
48   fTmpLabel[0]=-5;
49   fTmpLabel[1]=-5;
50   fTmpLabel[2]=-5;
51    for(int il=0; il<10;il++) fLabels[il]=-5;
52 }
53
54 //___________________________________________________________________________________
55 AliITSUpgradeClusterFinder::~AliITSUpgradeClusterFinder() {}
56 //___________________________________________________________________________________
57 void AliITSUpgradeClusterFinder::StartEvent() {
58   NewEvent();
59 }
60 //___________________________________________________________________________________
61 Int_t AliITSUpgradeClusterFinder::ProcessHitOnline(Int_t layer,  UInt_t col, UInt_t row, UShort_t charge, Int_t label[3]) {
62   if(layer < 6 ) { 
63     return ProcessHit(layer,col,row, charge, label);
64   }
65   else {
66     printf("ERROR: UpgradeClusterFinder::ProcessHitOnline: Out of bounds: layer,col,row, charge = %d,%d,%d,%d\n",layer ,col,row,charge);
67     return 1;
68   }
69 }
70 //___________________________________________________________________________________
71 Int_t AliITSUpgradeClusterFinder::ProcessHit(Int_t layer , UInt_t col, UInt_t row, UShort_t charge, Int_t label[3]) {
72   
73   if (layer>=6) {
74    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]); 
75    return 1;
76   }
77
78   // do we need to perform clustering on previous module?
79   if (fOldModule!=-1 && (Int_t)layer!=fOldModule) {
80     fChargeArray->AddLast(new TObjString(Form("%i %i %i %i %i %i",col,row,charge,label[0],label[1],label[2])));
81     DoModuleClustering(fOldModule,charge);
82     NewModule();
83   }
84   // add hit
85   fChargeArray->AddLast(new TObjString(Form("%i %i %i %i %i %i",col,row,charge,label[0],label[1],label[2])));
86
87   fOldModule=layer;
88   fHitCol[fNhitsLeft]=col;
89   fHitRow[fNhitsLeft]=row;
90   fHits[col][row]=kTRUE;
91   fTmpLabel[0]=label[0];
92   fTmpLabel[1]=label[1];
93   fTmpLabel[2]=label[2];
94   fNhitsLeft++;
95   return 0;
96 }
97 //___________________________________________________________________________________
98 void AliITSUpgradeClusterFinder::FinishEvent() {
99   if (fNhitsLeft>0) {
100     DoModuleClustering(fOldModule,fCharge);
101   }
102 }
103 //___________________________________________________________________________________
104 UInt_t AliITSUpgradeClusterFinder::GetClusterCount(Int_t layer) {
105   if (layer>=6) {
106     printf("ERROR: UpgradeClusterFinder::GetClusterCount: Module out of bounds: layer = %d\n",layer);
107     return 0;
108   }
109   return fClusterList[layer].GetNrEntries();
110 }
111 //___________________________________________________________________________________
112 Float_t AliITSUpgradeClusterFinder::GetClusterMeanCol(Int_t layer, UInt_t index) {
113   if (layer>=6) {
114     printf("ERROR: UpgradeClusterFinder::GetClusterMeanCol: Module out of bounds: layer = %d\n",layer);
115     return 0;
116   }
117   return fClusterList[layer].GetColIndex(index);
118 }
119 //___________________________________________________________________________________
120 Float_t AliITSUpgradeClusterFinder::GetClusterMeanRow(Int_t layer, UInt_t index) {
121   if (layer>=6) {
122     printf("ERROR: UpgradeClusterFinder::GetClusterMeanRow: Module out of bounds: layer = %d\n",layer);
123     return 0;
124   }
125   return fClusterList[layer].GetRowIndex(index);
126 }
127 //___________________________________________________________________________________
128 UInt_t AliITSUpgradeClusterFinder::GetClusterSize(Int_t layer, UInt_t index) {
129   if (layer>=6) {
130     printf("ERROR: UpgradeClusterFinder::GetClusterSize: Module out of bounds: layer = %d\n",layer);
131     return 0;
132   }
133   return fClusterList[layer].GetSizeIndex(index);
134 }
135 //___________________________________________________________________________________
136 UInt_t AliITSUpgradeClusterFinder::GetClusterWidthZ(Int_t layer, UInt_t index) {
137   if (layer>=6) {
138     printf("ERROR: UpgradeClusterFinder::GetClusterWidthZ: Module out of bounds: layer = %d\n",layer);
139     return 0;
140   }
141   return fClusterList[layer].GetWidthZIndex(index);
142 }
143 //___________________________________________________________________________________
144 UInt_t AliITSUpgradeClusterFinder::GetClusterWidthPhi(Int_t layer, UInt_t index) {
145   if (layer>=6) {
146     printf("ERROR: UpgradeClusterFinder::GetClusterWidthPhi: Module out of bounds: layer = %d\n",layer);
147     return 0;
148   }
149   return fClusterList[layer].GetWidthPhiIndex(index);
150 }
151 //___________________________________________________________________________________
152 UInt_t AliITSUpgradeClusterFinder::GetClusterType(Int_t layer, UInt_t index) {
153   if (layer>=6) {
154     printf("ERROR: UpgradeClusterFinder::GetClusterType: Module out of bounds: layer = %d\n",layer);
155     return 0;
156   }
157   return fClusterList[layer].GetTypeIndex(index);
158 }
159 //___________________________________________________________________________________
160 void AliITSUpgradeClusterFinder::PrintAllClusters() {
161   for (Int_t layer=0; layer<6; layer++) {
162     PrintClusters(layer);
163   }
164 }
165 //___________________________________________________________________________________
166 void AliITSUpgradeClusterFinder::PrintClusters(Int_t layer) {
167   if (layer>=6) {
168     printf("ERROR: UpgradeClusterFinder::PrintClusters: Out of bounds: layer = %d\n",layer);
169     return;
170   }
171   if(fClusterList[layer].GetNrEntries()==0) {
172     printf("no cluster list entries. Exiting... \n");
173     return;
174   }
175   for (UInt_t c=0; c<fClusterList[layer].GetNrEntries(); c++) {
176   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);  
177 }  
178 }
179 //___________________________________________________________________________________
180 void AliITSUpgradeClusterFinder::NewEvent() {
181   for (Int_t i=0; i<6; i++) {
182     fClusterList[i].Clear();
183   }
184   NewModule();
185   fOldModule = -1;
186 }
187 //___________________________________________________________________________________
188 void AliITSUpgradeClusterFinder::NewModule() {
189   fNhitsLeft=0;
190   memset(fHits,0,999*999*sizeof(Bool_t));
191 }
192 //___________________________________________________________________________________
193 Int_t AliITSUpgradeClusterFinder::DoModuleClustering(Int_t Layer, UShort_t charge) {
194   UInt_t maxHits=fNhitsLeft;
195   for (UInt_t hit=0; hit<maxHits; hit++) {
196     if (fClusterTypeFlag) fFindClusterType = kTRUE;
197
198     fClusterWidthMinCol = fHitCol[hit];
199     fClusterWidthMinRow = fHitRow[hit];
200     fClusterWidthMaxCol = fHitCol[hit];
201     fClusterWidthMaxRow = fHitRow[hit];
202     
203  
204     fClusterTypeOrigCol = fHitCol[hit];
205     fClusterTypeOrigRow = fHitRow[hit];
206     memset(fClusterTypeArea,0,kMAXCLUSTERTYPESIDEZ*kMAXCLUSTERTYPESIDEY*sizeof(Bool_t));
207     fColSum=0;
208     fRowSum=0;
209     UInt_t size = FindClusterRecu(fClusterTypeOrigCol,fClusterTypeOrigRow,charge);
210
211     if(size==1){
212       fCharge=GetPixelCharge(fColSum,fRowSum);
213       AddLabelIndex(fColSum,fRowSum);
214     }
215     if (size>0) {
216      if(size==2) printf("DoModuleClustering, size 2, labels :  %i  %i  %i \n",fLabels[0],fLabels[1],fLabels[2]);
217       fClusterList[Layer].Insert((Float_t)fColSum/size, (Float_t)fRowSum/size, size, GetClusterWidthZ(), GetClusterWidthPhi(), GetClusterType(size), fCharge,fLabels);
218       fCharge=0;
219       for(Int_t i=0; i<10; i++) fLabels[i]=-5;
220     }
221     if (fNhitsLeft==0) break;
222   }
223   return 0;
224 }
225 //___________________________________________________________________________________
226 UInt_t AliITSUpgradeClusterFinder::FindClusterRecu(Int_t col, Int_t row, UShort_t charge) {
227   if (col<0 || !fHits[col][row]) {return 0;}
228   fHits[col][row]=kFALSE;
229   fColSum+=col;
230   fRowSum+=row;
231   fNhitsLeft--;
232   UInt_t size=1;
233
234   if (col>fClusterWidthMaxCol) fClusterWidthMaxCol = col;
235   if (col<fClusterWidthMinCol) fClusterWidthMinCol = col;
236   if (row>fClusterWidthMaxRow) fClusterWidthMaxRow = row;
237   if (row<fClusterWidthMaxRow) fClusterWidthMinRow = row;
238
239   if (fFindClusterType) {
240     Short_t diffz = fClusterTypeOrigCol - col;
241     Short_t diffy = row - fClusterTypeOrigRow;
242     if (diffz>=kMAXCLUSTERTYPESIDEZ || diffy>=kMAXCLUSTERTYPESIDEY) {
243       fFindClusterType=kFALSE;
244     }
245     else {
246       if (diffz==-1) {
247         ShiftClusterTypeArea(kSHIFTRIGHT);
248         diffz=0;
249       }
250       if (diffy==-1) {
251         ShiftClusterTypeArea(kSHIFTDOWN);
252         diffy=0;
253       }
254       fClusterTypeArea[diffz][diffy] = kTRUE;
255     }
256   }
257   // straight:
258   size+=FindClusterRecu(col  ,row-1,charge);
259   fCharge+=GetPixelCharge(col,row-1);
260   AddLabelIndex(col,row-1);
261
262   size+=FindClusterRecu(col-1,row  ,charge);
263   fCharge+=GetPixelCharge(col-1,row);
264   AddLabelIndex(col-1,row);
265
266   size+=FindClusterRecu(col  ,row+1,charge);
267   fCharge+=GetPixelCharge(col,row+1);
268   AddLabelIndex(col,row+1);
269
270   size+=FindClusterRecu(col+1,row ,charge );
271   fCharge+=GetPixelCharge(col+1,row);
272   AddLabelIndex(col+1,row);
273
274   // diagonal:
275   size+=FindClusterRecu(col-1,row-1,charge);
276   fCharge+=GetPixelCharge(col-1,row-1);
277   AddLabelIndex(col-1,row-1);
278
279   size+=FindClusterRecu(col-1,row+1,charge);
280   fCharge+=GetPixelCharge(col-1,row+1);
281   AddLabelIndex(col-1,row+1);
282
283   size+=FindClusterRecu(col+1,row-1,charge);
284   fCharge+=GetPixelCharge(col+1,row-1);
285   AddLabelIndex(col+1,row-1);
286
287   size+=FindClusterRecu(col+1,row+1,charge);
288   fCharge+=GetPixelCharge(col+1,row+1);
289   AddLabelIndex(col+1,row+1);
290
291   return size;
292 }
293 //___________________________________________________________________________________
294 void AliITSUpgradeClusterFinder::ShiftClusterTypeArea(UInt_t direction) {
295   if (direction == kSHIFTRIGHT) {
296     fClusterTypeOrigCol++;
297     Bool_t tmpClusterTypeArea[kMAXCLUSTERTYPESIDEZ][kMAXCLUSTERTYPESIDEY];
298     memset(tmpClusterTypeArea,0,kMAXCLUSTERTYPESIDEZ*kMAXCLUSTERTYPESIDEY*sizeof(Bool_t));
299     for (UInt_t z=0; z<kMAXCLUSTERTYPESIDEZ; z++) {
300       for (UInt_t y=0; y<kMAXCLUSTERTYPESIDEY; y++) {
301         if (fClusterTypeArea[z][y]) {
302           if (z==kMAXCLUSTERTYPESIDEZ-1) {
303             fFindClusterType=kFALSE;
304             return;
305           }
306           else {
307             tmpClusterTypeArea[z+1][y] = kTRUE;
308           }
309         }
310       }
311     }
312     memcpy(fClusterTypeArea,tmpClusterTypeArea,kMAXCLUSTERTYPESIDEZ*kMAXCLUSTERTYPESIDEY*sizeof(Bool_t));
313   }
314   else if (direction == kSHIFTDOWN) {
315     fClusterTypeOrigRow--;
316     Bool_t tmpClusterTypeArea[kMAXCLUSTERTYPESIDEZ][kMAXCLUSTERTYPESIDEY];
317     memset(tmpClusterTypeArea,0,kMAXCLUSTERTYPESIDEZ*kMAXCLUSTERTYPESIDEY*sizeof(Bool_t));
318     for (UInt_t z=0; z<kMAXCLUSTERTYPESIDEZ; z++) {
319       for (UInt_t y=0; y<kMAXCLUSTERTYPESIDEY; y++) {
320         if (fClusterTypeArea[z][y]) {
321           if (y==kMAXCLUSTERTYPESIDEY-1) {
322             fFindClusterType=kFALSE;
323             return;
324           }
325           else {
326             tmpClusterTypeArea[z][y+1] = kTRUE;
327           }
328         }
329       }
330     }
331     memcpy(fClusterTypeArea,tmpClusterTypeArea,kMAXCLUSTERTYPESIDEZ*kMAXCLUSTERTYPESIDEY*sizeof(Bool_t));
332   }
333 }
334
335
336
337 //___________________________________________________________________________________
338 UShort_t AliITSUpgradeClusterFinder::GetCharge(Int_t layer,UInt_t index ) {
339   return fClusterList[layer].GetCharge(index);
340 }
341 //___________________________________________________________________________________
342 Int_t * AliITSUpgradeClusterFinder::GetLabels(Int_t layer,UInt_t index) {
343        return fClusterList[layer].GetLabels(index);
344 }
345
346 //___________________________________________________________________________________
347 UInt_t AliITSUpgradeClusterFinder::GetClusterWidthZ() {
348   return fClusterWidthMaxCol-fClusterWidthMinCol+1;
349 }
350 //___________________________________________________________________________________
351 UInt_t AliITSUpgradeClusterFinder::GetClusterWidthPhi() {
352   return fClusterWidthMaxRow-fClusterWidthMinRow+1;
353 }
354 //___________________________________________________________________________________
355 UInt_t AliITSUpgradeClusterFinder::GetClusterType(UInt_t size) {
356   // category 'other':
357   if (!fFindClusterType || size>4) 
358     return 0;
359
360   // X
361   if (size==1) 
362     return 1;
363
364   // X
365   // X
366   if (size==2 && 
367       fClusterTypeArea[0][1] && 
368       fClusterTypeArea[0][0] ) 
369     return 2;
370
371   // XX
372   if (size==2 && 
373       fClusterTypeArea[1][0] &&
374       fClusterTypeArea[0][0] ) 
375     return 3;
376
377   // X
378   // X
379   // X
380   if (size==3 && 
381       fClusterTypeArea[0][2] && 
382       fClusterTypeArea[0][1] && 
383       fClusterTypeArea[0][0] ) 
384     return 4;
385   
386   // XX
387   // X      and equivalent...
388   if (size==3 && 
389       (
390        (
391         fClusterTypeArea[0][0] && 
392         fClusterTypeArea[0][1] && 
393         (fClusterTypeArea[1][1] || 
394          fClusterTypeArea[1][0])
395         )
396        ||
397        (
398         fClusterTypeArea[1][0] && 
399         fClusterTypeArea[1][1] && 
400         (fClusterTypeArea[0][1] || 
401          fClusterTypeArea[0][0])
402         )
403        )
404       ) 
405     return 5;
406   
407   // XXX
408   if (size==3 && 
409       fClusterTypeArea[2][0] && 
410       fClusterTypeArea[1][0] && 
411       fClusterTypeArea[0][0] )
412     return 6;
413   
414   // X
415   // X
416   // X
417   // X
418   if (size==4 &&
419       fClusterTypeArea[0][3] && 
420       fClusterTypeArea[0][2] &&
421       fClusterTypeArea[0][1] && 
422       fClusterTypeArea[0][0] ) 
423     return 7;
424   
425   // XX
426   // XX
427   if (size==4 && 
428       fClusterTypeArea[1][1] && 
429       fClusterTypeArea[1][0] && 
430       fClusterTypeArea[0][1] &&
431       fClusterTypeArea[0][0] ) 
432     return 8;
433
434   // XX
435   // X
436   // X     and equivalent...
437   if (size==4 &&
438       (
439        (
440         fClusterTypeArea[1][0] && 
441         fClusterTypeArea[1][1] && 
442         fClusterTypeArea[1][2] && 
443         fClusterTypeArea[0][0]
444         )
445        ||
446        (
447         fClusterTypeArea[1][0] && 
448         fClusterTypeArea[1][1] && 
449         fClusterTypeArea[1][2] && 
450         fClusterTypeArea[0][2]
451         )
452        ||
453        (
454         fClusterTypeArea[0][0] && 
455         fClusterTypeArea[0][1] && 
456         fClusterTypeArea[0][2] && 
457         fClusterTypeArea[1][0]
458         )
459        ||
460        (
461         fClusterTypeArea[0][0] && 
462         fClusterTypeArea[0][1] && 
463         fClusterTypeArea[0][2] && 
464         fClusterTypeArea[1][2]
465         )
466        )
467       )
468     return 9;
469
470   // X
471   // XX
472   // X     and equivalent...
473   if (size==4 &&
474       (
475        (
476         fClusterTypeArea[1][0] && 
477         fClusterTypeArea[1][1] && 
478         fClusterTypeArea[1][2] && 
479         fClusterTypeArea[0][1]
480         )
481        ||
482        (
483         fClusterTypeArea[0][0] && 
484         fClusterTypeArea[0][1] && 
485         fClusterTypeArea[0][2] && 
486         fClusterTypeArea[1][1]
487         )
488        )
489       )
490     return 10;
491
492   // X
493   // XXX    and equivalent...
494   if (size==4 &&
495       (
496        (
497         fClusterTypeArea[2][0] && 
498         fClusterTypeArea[2][1] && 
499         fClusterTypeArea[1][1] && 
500         fClusterTypeArea[0][1]
501         )
502        ||
503        (
504         fClusterTypeArea[0][0] && 
505         fClusterTypeArea[2][1] && 
506         fClusterTypeArea[1][1] && 
507         fClusterTypeArea[0][1]
508         )
509        ||
510        (
511         fClusterTypeArea[2][1] && 
512         fClusterTypeArea[2][0] && 
513         fClusterTypeArea[1][0] && 
514         fClusterTypeArea[0][0]
515         )
516        ||
517        (
518         fClusterTypeArea[0][1] && 
519         fClusterTypeArea[2][0] && 
520         fClusterTypeArea[1][0] && 
521         fClusterTypeArea[0][0]
522         )
523        )
524       )
525     return 11;
526
527   //  X
528   // XXX     and equivalent...
529   if (size==4 &&
530       (
531        (
532         fClusterTypeArea[1][0] && 
533         fClusterTypeArea[2][1] && 
534         fClusterTypeArea[1][1] && 
535         fClusterTypeArea[0][1]
536         )
537        ||
538        (
539         fClusterTypeArea[1][1] && 
540         fClusterTypeArea[2][0] && 
541         fClusterTypeArea[1][0] && 
542         fClusterTypeArea[0][0]
543         )
544        )
545       )
546     return 12;
547
548   //  X
549   // X      and equivalent...
550   if (size==2 &&
551       (
552        (
553         fClusterTypeArea[0][0] &&
554         fClusterTypeArea[1][1] 
555         )
556        ||
557        (
558         fClusterTypeArea[1][0] &&
559         fClusterTypeArea[0][1] 
560         )
561        )
562       )
563     return 13;
564
565   //  X
566   //  X
567   // X      and equivalent...
568   if (size==3 &&
569       (
570        (
571         fClusterTypeArea[0][0] &&
572         fClusterTypeArea[0][1] && 
573         fClusterTypeArea[1][2] 
574         )
575        ||
576        (
577         fClusterTypeArea[1][0] &&
578         fClusterTypeArea[1][1] && 
579         fClusterTypeArea[0][2] 
580         )
581        ||
582        (
583         fClusterTypeArea[1][0] &&
584         fClusterTypeArea[0][1] && 
585         fClusterTypeArea[0][2] 
586         )
587        ||
588        (
589         fClusterTypeArea[0][0] &&
590         fClusterTypeArea[1][1] && 
591         fClusterTypeArea[1][2] 
592         )
593        )
594       )
595     return 14;
596
597   //  XX
598   // X      and equivalent...
599   if (size==3 &&
600       (
601        (
602         fClusterTypeArea[0][0] &&
603         fClusterTypeArea[1][0] && 
604         fClusterTypeArea[2][1] 
605         )
606        ||
607        (
608         fClusterTypeArea[0][1] &&
609         fClusterTypeArea[1][1] && 
610         fClusterTypeArea[2][0] 
611         )
612        ||
613        (
614         fClusterTypeArea[0][0] &&
615         fClusterTypeArea[1][1] && 
616         fClusterTypeArea[2][1] 
617         )
618        ||
619        (
620         fClusterTypeArea[0][1] &&
621         fClusterTypeArea[1][0] && 
622         fClusterTypeArea[2][0] 
623         )
624        )
625       )
626     return 15;
627
628
629
630   return 0;
631 }
632
633 //___________________________________________________________________________________________________
634
635 UInt_t AliITSUpgradeClusterFinder::GetPixelCharge(UInt_t col, UInt_t row){
636   Int_t q=0;
637   for(Int_t entry =0; entry < fChargeArray->GetEntries(); entry++){
638     TObjString *s = (TObjString*)fChargeArray->At(entry);
639     TString name = s->GetString();
640     if(!name.Contains(Form("%i %i",col,row)))
641       continue;
642     TObjArray *array = name.Tokenize(" ");
643     TString charge = ((TObjString*)array->At(2))->String();
644     TString rowS, colS;
645     rowS = ((TObjString*)array->At(0))->String();
646     colS = ((TObjString*)array->At(1))->String();
647     delete array;
648     q=charge.Atoi();
649     return q;
650
651   }
652   return q;
653 }
654 //____________________________________________________
655
656 void AliITSUpgradeClusterFinder::AddLabelIndex(UInt_t col, UInt_t row){
657
658  for(Int_t entry =0; entry < fChargeArray->GetEntries(); entry++){
659     TObjString *s = (TObjString*)fChargeArray->At(entry);
660     TString name = s->GetString();
661  if(!name.Contains(Form("%i %i",col,row)))
662   continue;
663  TObjArray *array = name.Tokenize(" ");
664  TString index[3];
665  Int_t label[3];
666   for(Int_t i=0; i<3; i++){
667   index[i]=((TObjString*)array->At(3+i))->String();
668   label[i]=index[i].Atoi();
669
670   }
671  SetLabels(label);
672  delete array;
673  }
674 }
675 //____________________________________________________
676
677 void AliITSUpgradeClusterFinder::SetLabels(Int_t label[3]){
678
679
680  Bool_t isAssigned[3]={kFALSE,kFALSE,kFALSE};
681
682  for(Int_t i=0; i<10; i++){
683   for(Int_t k=0; k<3; k++){
684    //printf(" fLabels[%i]=%i  label[%i]=%i \n",i,fLabels[i],k,label[k]);
685     if(fLabels[i]!=label[k] && label[k]>-1 && fLabels[i]<0) {
686       if(!isAssigned[k]) {
687      //  printf("assignign...\n.");
688     printf("Assignign  fLabels[%i]=%i  label[%i]=%i \n",i,fLabels[i],k,label[k]);
689         fLabels[i+k]=label[k];
690         isAssigned[k]=kTRUE;
691       }
692      }
693    }
694  }
695 }
696