1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 ////////////////////////////////////////////////////////////////////////////
20 // Base Class used to find //
21 // the reconstructed points for ITS Upgrade //
23 ////////////////////////////////////////////////////////////////////////////
25 #include "AliITSUpgradeClusterFinder.h"
26 #include "AliITSsegmentationUpgrade.h"
27 #include "AliITSRecPointU.h"
28 #include "AliITSDigitUpgrade.h"
29 #include "AliITSRawStreamSPD.h"
32 #include <TObjString.h>
33 #include <TClonesArray.h>
37 AliITSUpgradeClusterFinder::AliITSUpgradeClusterFinder() :
41 fClusterTypeFlag(kTRUE),
42 fFindClusterType(kFALSE),
43 fClusterTypeOrigCol(0),
44 fClusterTypeOrigRow(0),
45 fColSum(0),fRowSum(0),
47 fClusterWidthMaxCol(0),
48 fClusterWidthMinCol(0),
49 fClusterWidthMaxRow(0),
50 fClusterWidthMinRow(0),
56 // Default constructor
58 AliITSsegmentationUpgrade *s = new AliITSsegmentationUpgrade();
59 fNlayers = s->GetNLayers();
62 fChargeArray = new TObjArray();
63 fChargeArray->SetOwner(kTRUE);
64 fRecPoints = new TClonesArray("AliITSRecPointU",3000);
68 for(int il=0; il<10;il++) fLabels[il]=-5;
69 for(int k=0; k<999999; k++){
75 //___________________________________________________________________________________
76 AliITSUpgradeClusterFinder::~AliITSUpgradeClusterFinder() {
77 if(fChargeArray) delete fChargeArray;
78 if(fRecPoints) delete fRecPoints;
80 //___________________________________________________________________________________
81 void AliITSUpgradeClusterFinder::StartEvent() {
84 //___________________________________________________________________________________
85 Int_t AliITSUpgradeClusterFinder::ProcessHit(Int_t layer , UInt_t col, UInt_t row, UShort_t charge, Int_t label[3]) {
87 // Adds one pixel to the cluster
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]);
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);
102 fChargeArray->AddLast(new TObjString(Form("%i %i %i %i %i %i",col,row,charge,label[0],label[1],label[2])));
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];
114 //___________________________________________________________________________________
115 void AliITSUpgradeClusterFinder::FinishEvent() {
117 DoModuleClustering(fOldModule,fCharge);
120 //___________________________________________________________________________________
121 UInt_t AliITSUpgradeClusterFinder::GetClusterCount(Int_t layer) const {
123 // number of clusters in the layer
125 if (layer>=fNlayers ) {
126 printf("ERROR: UpgradeClusterFinder::GetClusterCount: Module out of bounds: layer = %d\n",layer);
129 return fClusterList[layer].GetNrEntries();
131 //___________________________________________________________________________________
132 Float_t AliITSUpgradeClusterFinder::GetClusterMeanCol(Int_t layer, UInt_t index) {
134 // cluster position in terms of colums : mean column ID
136 if (layer>=fNlayers) {
137 printf("ERROR: UpgradeClusterFinder::GetClusterMeanCol: Module out of bounds: layer = %d\n",layer);
140 return fClusterList[layer].GetColIndex(index);
142 //___________________________________________________________________________________
143 Float_t AliITSUpgradeClusterFinder::GetClusterMeanRow(Int_t layer, UInt_t index) {
145 // cluster position in terms of rows : mean row ID
147 if (layer>=fNlayers) {
148 printf("ERROR: UpgradeClusterFinder::GetClusterMeanRow: Module out of bounds: layer = %d\n",layer);
151 return fClusterList[layer].GetRowIndex(index);
153 //___________________________________________________________________________________
154 UInt_t AliITSUpgradeClusterFinder::GetClusterSize(Int_t layer, UInt_t index) {
156 // total number of pixels of the cluster
158 if (layer>=fNlayers) {
159 printf("ERROR: UpgradeClusterFinder::GetClusterSize: Module out of bounds: layer = %d\n",layer);
162 return fClusterList[layer].GetSizeIndex(index);
164 //___________________________________________________________________________________
165 UInt_t AliITSUpgradeClusterFinder::GetClusterWidthZ(Int_t layer, UInt_t index) {
167 // # pixels of the cluster in Z direction
170 if (layer>=fNlayers) {
171 printf("ERROR: UpgradeClusterFinder::GetClusterWidthZ: Module out of bounds: layer = %d\n",layer);
174 return fClusterList[layer].GetWidthZIndex(index);
176 //___________________________________________________________________________________
177 UInt_t AliITSUpgradeClusterFinder::GetClusterWidthPhi(Int_t layer, UInt_t index) {
179 // # pixels of the cluster in phi direction (XY plane)
182 if (layer>=fNlayers) {
183 printf("ERROR: UpgradeClusterFinder::GetClusterWidthPhi: Module out of bounds: layer = %d\n",layer);
186 return fClusterList[layer].GetWidthPhiIndex(index);
188 //___________________________________________________________________________________
189 UInt_t AliITSUpgradeClusterFinder::GetClusterType(Int_t layer, UInt_t index) {
194 if (layer>=fNlayers) {
195 printf("ERROR: UpgradeClusterFinder::GetClusterType: Module out of bounds: layer = %d\n",layer);
198 return fClusterList[layer].GetTypeIndex(index);
200 //___________________________________________________________________________________
201 void AliITSUpgradeClusterFinder::PrintAllClusters() {
203 // printout of the cluster information
206 for (Int_t layer=0; layer<fNlayers; layer++) {
207 PrintClusters(layer);
210 //___________________________________________________________________________________
211 void AliITSUpgradeClusterFinder::PrintClusters(Int_t layer) {
213 // printout of the cluster information
216 if (layer>=fNlayers) {
217 printf("ERROR: UpgradeClusterFinder::PrintClusters: Out of bounds: layer = %d\n",layer);
220 if(fClusterList[layer].GetNrEntries()==0) {
221 printf("no cluster list entries. Exiting... \n");
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);
228 //___________________________________________________________________________________
229 void AliITSUpgradeClusterFinder::NewEvent() {
231 // Cleaning up and preparation for the clustering procedure
234 for (Int_t i=0; i<fNlayers; i++) {
235 fClusterList[i].Clear();
240 //___________________________________________________________________________________
241 void AliITSUpgradeClusterFinder::NewModule() {
246 memset(fHits,0,999*999*sizeof(Bool_t));
247 fChargeArray->Clear();
249 //___________________________________________________________________________________
250 Int_t AliITSUpgradeClusterFinder::DoModuleClustering(Int_t Layer, UShort_t charge) {
252 // Clustering and cluster-list container filling
254 UInt_t maxHits=fNhitsLeft;
255 for (UInt_t hit=0; hit<maxHits; hit++) {
256 if (fClusterTypeFlag) fFindClusterType = kTRUE;
258 fClusterWidthMinCol = fHitCol[hit];
259 fClusterWidthMinRow = fHitRow[hit];
260 fClusterWidthMaxCol = fHitCol[hit];
261 fClusterWidthMaxRow = fHitRow[hit];
263 fClusterTypeOrigCol = fHitCol[hit];
264 fClusterTypeOrigRow = fHitRow[hit];
265 memset(fClusterTypeArea,0,kMAXCLUSTERTYPESIDEZ*kMAXCLUSTERTYPESIDEY*sizeof(Bool_t));
268 UInt_t size = FindClusterRecu(fClusterTypeOrigCol,fClusterTypeOrigRow,charge);
270 fCharge=GetPixelCharge(fColSum,fRowSum);
271 AddLabelIndex(fColSum,fRowSum);
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);
277 for(Int_t i=0; i<10; i++) fLabels[i]=-5;
279 if (fNhitsLeft==0) break;
283 //___________________________________________________________________________________
284 UInt_t AliITSUpgradeClusterFinder::FindClusterRecu(Int_t col, Int_t row, UShort_t charge) {
286 // Pixel selection for the clusters (adjiacent pixels or pixels at the corners)
288 if (col<0 || !fHits[col][row]) {
291 fHits[col][row]=kFALSE;
297 if (col>fClusterWidthMaxCol) fClusterWidthMaxCol = col;
298 if (col<fClusterWidthMinCol) fClusterWidthMinCol = col;
299 if (row>fClusterWidthMaxRow) fClusterWidthMaxRow = row;
300 if (row<fClusterWidthMaxRow) fClusterWidthMinRow = row;
302 if (fFindClusterType) {
303 Short_t diffz = fClusterTypeOrigCol - col;
304 Short_t diffy = row - fClusterTypeOrigRow;
305 if (diffz>=kMAXCLUSTERTYPESIDEZ || diffy>=kMAXCLUSTERTYPESIDEY) {
306 fFindClusterType=kFALSE;
310 ShiftClusterTypeArea(kSHIFTRIGHT);
314 ShiftClusterTypeArea(kSHIFTDOWN);
317 fClusterTypeArea[diffz][diffy] = kTRUE;
321 size+=FindClusterRecu(col ,row-1,charge);
322 fCharge+=GetPixelCharge(col,row-1);
323 AddLabelIndex(col,row-1);
325 size+=FindClusterRecu(col-1,row ,charge);
326 fCharge+=GetPixelCharge(col-1,row);
327 AddLabelIndex(col-1,row);
329 size+=FindClusterRecu(col ,row+1,charge);
330 fCharge+=GetPixelCharge(col,row+1);
331 AddLabelIndex(col,row+1);
333 size+=FindClusterRecu(col+1,row ,charge );
334 fCharge+=GetPixelCharge(col+1,row);
335 AddLabelIndex(col+1,row);
338 size+=FindClusterRecu(col-1,row-1,charge);
339 fCharge+=GetPixelCharge(col-1,row-1);
340 AddLabelIndex(col-1,row-1);
342 size+=FindClusterRecu(col-1,row+1,charge);
343 fCharge+=GetPixelCharge(col-1,row+1);
344 AddLabelIndex(col-1,row+1);
346 size+=FindClusterRecu(col+1,row-1,charge);
347 fCharge+=GetPixelCharge(col+1,row-1);
348 AddLabelIndex(col+1,row-1);
350 size+=FindClusterRecu(col+1,row+1,charge);
351 fCharge+=GetPixelCharge(col+1,row+1);
352 AddLabelIndex(col+1,row+1);
355 //___________________________________________________________________________________
356 void AliITSUpgradeClusterFinder::ShiftClusterTypeArea(UInt_t direction) {
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;
373 tmpClusterTypeArea[z+1][y] = kTRUE;
378 memcpy(fClusterTypeArea,tmpClusterTypeArea,kMAXCLUSTERTYPESIDEZ*kMAXCLUSTERTYPESIDEY*sizeof(Bool_t));
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;
392 tmpClusterTypeArea[z][y+1] = kTRUE;
397 memcpy(fClusterTypeArea,tmpClusterTypeArea,kMAXCLUSTERTYPESIDEZ*kMAXCLUSTERTYPESIDEY*sizeof(Bool_t));
400 //___________________________________________________________________________________
401 UShort_t AliITSUpgradeClusterFinder::GetCharge(Int_t layer,UInt_t index ) {
402 return fClusterList[layer].GetCharge(index);
404 //___________________________________________________________________________________
405 Int_t * AliITSUpgradeClusterFinder::GetLabels(Int_t layer,UInt_t index) {
406 return fClusterList[layer].GetLabels(index);
409 //___________________________________________________________________________________
410 UInt_t AliITSUpgradeClusterFinder::GetClusterWidthZ() {
411 return fClusterWidthMaxCol-fClusterWidthMinCol+1;
413 //___________________________________________________________________________________
414 UInt_t AliITSUpgradeClusterFinder::GetClusterWidthPhi() {
415 return fClusterWidthMaxRow-fClusterWidthMinRow+1;
417 //___________________________________________________________________________________
418 UInt_t AliITSUpgradeClusterFinder::GetClusterType(UInt_t size) {
422 if (!fFindClusterType || size>4)
432 fClusterTypeArea[0][1] &&
433 fClusterTypeArea[0][0] )
438 fClusterTypeArea[1][0] &&
439 fClusterTypeArea[0][0] )
446 fClusterTypeArea[0][2] &&
447 fClusterTypeArea[0][1] &&
448 fClusterTypeArea[0][0] )
452 // X and equivalent...
456 fClusterTypeArea[0][0] &&
457 fClusterTypeArea[0][1] &&
458 (fClusterTypeArea[1][1] ||
459 fClusterTypeArea[1][0])
463 fClusterTypeArea[1][0] &&
464 fClusterTypeArea[1][1] &&
465 (fClusterTypeArea[0][1] ||
466 fClusterTypeArea[0][0])
474 fClusterTypeArea[2][0] &&
475 fClusterTypeArea[1][0] &&
476 fClusterTypeArea[0][0] )
484 fClusterTypeArea[0][3] &&
485 fClusterTypeArea[0][2] &&
486 fClusterTypeArea[0][1] &&
487 fClusterTypeArea[0][0] )
493 fClusterTypeArea[1][1] &&
494 fClusterTypeArea[1][0] &&
495 fClusterTypeArea[0][1] &&
496 fClusterTypeArea[0][0] )
501 // X and equivalent...
505 fClusterTypeArea[1][0] &&
506 fClusterTypeArea[1][1] &&
507 fClusterTypeArea[1][2] &&
508 fClusterTypeArea[0][0]
512 fClusterTypeArea[1][0] &&
513 fClusterTypeArea[1][1] &&
514 fClusterTypeArea[1][2] &&
515 fClusterTypeArea[0][2]
519 fClusterTypeArea[0][0] &&
520 fClusterTypeArea[0][1] &&
521 fClusterTypeArea[0][2] &&
522 fClusterTypeArea[1][0]
526 fClusterTypeArea[0][0] &&
527 fClusterTypeArea[0][1] &&
528 fClusterTypeArea[0][2] &&
529 fClusterTypeArea[1][2]
537 // X and equivalent...
541 fClusterTypeArea[1][0] &&
542 fClusterTypeArea[1][1] &&
543 fClusterTypeArea[1][2] &&
544 fClusterTypeArea[0][1]
548 fClusterTypeArea[0][0] &&
549 fClusterTypeArea[0][1] &&
550 fClusterTypeArea[0][2] &&
551 fClusterTypeArea[1][1]
558 // XXX and equivalent...
562 fClusterTypeArea[2][0] &&
563 fClusterTypeArea[2][1] &&
564 fClusterTypeArea[1][1] &&
565 fClusterTypeArea[0][1]
569 fClusterTypeArea[0][0] &&
570 fClusterTypeArea[2][1] &&
571 fClusterTypeArea[1][1] &&
572 fClusterTypeArea[0][1]
576 fClusterTypeArea[2][1] &&
577 fClusterTypeArea[2][0] &&
578 fClusterTypeArea[1][0] &&
579 fClusterTypeArea[0][0]
583 fClusterTypeArea[0][1] &&
584 fClusterTypeArea[2][0] &&
585 fClusterTypeArea[1][0] &&
586 fClusterTypeArea[0][0]
593 // XXX and equivalent...
597 fClusterTypeArea[1][0] &&
598 fClusterTypeArea[2][1] &&
599 fClusterTypeArea[1][1] &&
600 fClusterTypeArea[0][1]
604 fClusterTypeArea[1][1] &&
605 fClusterTypeArea[2][0] &&
606 fClusterTypeArea[1][0] &&
607 fClusterTypeArea[0][0]
614 // X and equivalent...
618 fClusterTypeArea[0][0] &&
619 fClusterTypeArea[1][1]
623 fClusterTypeArea[1][0] &&
624 fClusterTypeArea[0][1]
632 // X and equivalent...
636 fClusterTypeArea[0][0] &&
637 fClusterTypeArea[0][1] &&
638 fClusterTypeArea[1][2]
642 fClusterTypeArea[1][0] &&
643 fClusterTypeArea[1][1] &&
644 fClusterTypeArea[0][2]
648 fClusterTypeArea[1][0] &&
649 fClusterTypeArea[0][1] &&
650 fClusterTypeArea[0][2]
654 fClusterTypeArea[0][0] &&
655 fClusterTypeArea[1][1] &&
656 fClusterTypeArea[1][2]
663 // X and equivalent...
667 fClusterTypeArea[0][0] &&
668 fClusterTypeArea[1][0] &&
669 fClusterTypeArea[2][1]
673 fClusterTypeArea[0][1] &&
674 fClusterTypeArea[1][1] &&
675 fClusterTypeArea[2][0]
679 fClusterTypeArea[0][0] &&
680 fClusterTypeArea[1][1] &&
681 fClusterTypeArea[2][1]
685 fClusterTypeArea[0][1] &&
686 fClusterTypeArea[1][0] &&
687 fClusterTypeArea[2][0]
698 //___________________________________________________________________________________________________
699 UInt_t AliITSUpgradeClusterFinder::GetPixelCharge(UInt_t col, UInt_t row){
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)))
709 TObjArray *array = name.Tokenize(" ");
710 TString charge = ((TObjString*)array->At(2))->String();
712 rowS = ((TObjString*)array->At(0))->String();
713 colS = ((TObjString*)array->At(1))->String();
721 //____________________________________________________
723 void AliITSUpgradeClusterFinder::AddLabelIndex(UInt_t col, UInt_t row){
725 // Adding cluster labels
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)))
733 TObjArray *array = name.Tokenize(" ");
736 for(Int_t i=0; i<3; i++){
737 index[i]=((TObjString*)array->At(3+i))->String();
738 label[i]=index[i].Atoi();
745 //____________________________________________________
747 void AliITSUpgradeClusterFinder::SetLabels(Int_t label[3]){
749 //Set the MC lables to the cluster
751 Bool_t isAssigned[3]={kFALSE,kFALSE,kFALSE};
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) {
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];
767 //____________________________________________________
768 void AliITSUpgradeClusterFinder::MakeRecPointBranch(TTree *treeR){
770 // Creating the branch (see AliITSUpgradeReconstructor::Reconstruct)
773 if(!fRecPoints)fRecPoints = new TClonesArray("AliITSRecPointU",1000);
775 TBranch *branch = treeR->GetBranch("ITSRecPoints");
777 else branch = treeR->Branch("ITSRecPoints",&fRecPoints);
780 //____________________________________________________
781 void AliITSUpgradeClusterFinder::SetRecPointTreeAddress(TTree *treeR){
783 // Addressing the branch (see AliITSUpgradeReconstructor::Reconstruct)
786 if(!fRecPoints) fRecPoints = new TClonesArray("AliITSRecPointU",1000);
789 branch = treeR->GetBranch("ITSRecPoints");
791 branch->SetAddress(&fRecPoints);
792 } else AliError("No RecPoint branch available");
795 //____________________________________________________
796 void AliITSUpgradeClusterFinder::DigitsToRecPoints(const TObjArray *digList) {
798 // the clusterization is performed here
800 AliITSsegmentationUpgrade *segmentation = new AliITSsegmentationUpgrade();
801 AliITSRecPointU recpnt;
803 TClonesArray &lrecp = *fRecPoints;
805 for(Int_t ilayer=0; ilayer < fNlayers ;ilayer ++){
806 TClonesArray *pArrDig= (TClonesArray*)digList->At(ilayer);
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());
815 ProcessHit(ilayer,colz, rowx,(Short_t)hitcharge,dig->GetTracks());
819 for(UInt_t nClu = 0; nClu < GetClusterCount(ilayer); nClu++){
820 UShort_t charge = GetCharge(ilayer, nClu);
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);}
830 Double_t xzl2[2]={0.,0.};
831 Double_t xPixC2,zPixC2=0.;
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();//
850 phiclu1deg=180.*phiclu1rad/TMath::Pi();// in deg
851 Int_t ladder;// virtual segmentation starting from the cluster phi
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)
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)
863 yclu1=xclu1*phiTrk; // tracking system coordinate: r*phi
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
877 if(segmentation) delete segmentation;