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