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 "AliITSRecPoint.h"
28 #include "AliITSDigitUpgrade.h"
29 #include "AliITSRawStreamSPD.h"
32 #include <TObjString.h>
33 #include <TClonesArray.h>
37 AliITSUpgradeClusterFinder::AliITSUpgradeClusterFinder() :
40 fClusterTypeFlag(kTRUE),
41 fFindClusterType(kFALSE),
42 fClusterTypeOrigCol(0),
43 fClusterTypeOrigRow(0),
44 fColSum(0),fRowSum(0),
46 fClusterWidthMaxCol(0),
47 fClusterWidthMinCol(0),
48 fClusterWidthMaxRow(0),
49 fClusterWidthMinRow(0),
54 // Default constructor
57 fChargeArray = new TObjArray();
58 fRecPoints = new TClonesArray("AliITSRecPoint",3000);
62 for(int il=0; il<10;il++) fLabels[il]=-5;
65 //___________________________________________________________________________________
66 AliITSUpgradeClusterFinder::~AliITSUpgradeClusterFinder() {
67 if(fChargeArray) delete fChargeArray;
68 if(fRecPoints) delete fRecPoints;
70 //___________________________________________________________________________________
71 void AliITSUpgradeClusterFinder::StartEvent() {
74 //___________________________________________________________________________________
75 Int_t AliITSUpgradeClusterFinder::ProcessHit(Int_t layer , UInt_t col, UInt_t row, UShort_t charge, Int_t label[3]) {
77 // Adds one pixel to the cluster
81 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]);
85 // do we need to perform clustering on previous module?
86 if (fOldModule!=-1 && (Int_t)layer!=fOldModule) {
87 fChargeArray->AddLast(new TObjString(Form("%i %i %i %i %i %i",col,row,charge,label[0],label[1],label[2])));
88 DoModuleClustering(fOldModule,charge);
92 fChargeArray->AddLast(new TObjString(Form("%i %i %i %i %i %i",col,row,charge,label[0],label[1],label[2])));
95 fHitCol[fNhitsLeft]=col;
96 fHitRow[fNhitsLeft]=row;
97 fHits[col][row]=kTRUE;
98 fTmpLabel[0]=label[0];
99 fTmpLabel[1]=label[1];
100 fTmpLabel[2]=label[2];
104 //___________________________________________________________________________________
105 void AliITSUpgradeClusterFinder::FinishEvent() {
107 DoModuleClustering(fOldModule,fCharge);
110 //___________________________________________________________________________________
111 UInt_t AliITSUpgradeClusterFinder::GetClusterCount(Int_t layer) const {
113 // number of clusters in the layer
116 printf("ERROR: UpgradeClusterFinder::GetClusterCount: Module out of bounds: layer = %d\n",layer);
119 return fClusterList[layer].GetNrEntries();
121 //___________________________________________________________________________________
122 Float_t AliITSUpgradeClusterFinder::GetClusterMeanCol(Int_t layer, UInt_t index) {
124 // cluster position in terms of colums : mean column ID
127 printf("ERROR: UpgradeClusterFinder::GetClusterMeanCol: Module out of bounds: layer = %d\n",layer);
130 return fClusterList[layer].GetColIndex(index);
132 //___________________________________________________________________________________
133 Float_t AliITSUpgradeClusterFinder::GetClusterMeanRow(Int_t layer, UInt_t index) {
135 // cluster position in terms of rows : mean row ID
138 printf("ERROR: UpgradeClusterFinder::GetClusterMeanRow: Module out of bounds: layer = %d\n",layer);
141 return fClusterList[layer].GetRowIndex(index);
143 //___________________________________________________________________________________
144 UInt_t AliITSUpgradeClusterFinder::GetClusterSize(Int_t layer, UInt_t index) {
146 // total number of pixels of the cluster
149 printf("ERROR: UpgradeClusterFinder::GetClusterSize: Module out of bounds: layer = %d\n",layer);
152 return fClusterList[layer].GetSizeIndex(index);
154 //___________________________________________________________________________________
155 UInt_t AliITSUpgradeClusterFinder::GetClusterWidthZ(Int_t layer, UInt_t index) {
157 // # pixels of the cluster in Z direction
161 printf("ERROR: UpgradeClusterFinder::GetClusterWidthZ: Module out of bounds: layer = %d\n",layer);
164 return fClusterList[layer].GetWidthZIndex(index);
166 //___________________________________________________________________________________
167 UInt_t AliITSUpgradeClusterFinder::GetClusterWidthPhi(Int_t layer, UInt_t index) {
169 // # pixels of the cluster in phi direction (XY plane)
173 printf("ERROR: UpgradeClusterFinder::GetClusterWidthPhi: Module out of bounds: layer = %d\n",layer);
176 return fClusterList[layer].GetWidthPhiIndex(index);
178 //___________________________________________________________________________________
179 UInt_t AliITSUpgradeClusterFinder::GetClusterType(Int_t layer, UInt_t index) {
185 printf("ERROR: UpgradeClusterFinder::GetClusterType: Module out of bounds: layer = %d\n",layer);
188 return fClusterList[layer].GetTypeIndex(index);
190 //___________________________________________________________________________________
191 void AliITSUpgradeClusterFinder::PrintAllClusters() {
193 // printout of the cluster information
196 for (Int_t layer=0; layer<6; layer++) {
197 PrintClusters(layer);
200 //___________________________________________________________________________________
201 void AliITSUpgradeClusterFinder::PrintClusters(Int_t layer) {
203 // printout of the cluster information
207 printf("ERROR: UpgradeClusterFinder::PrintClusters: Out of bounds: layer = %d\n",layer);
210 if(fClusterList[layer].GetNrEntries()==0) {
211 printf("no cluster list entries. Exiting... \n");
214 for (UInt_t c=0; c<fClusterList[layer].GetNrEntries(); c++) {
215 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);
218 //___________________________________________________________________________________
219 void AliITSUpgradeClusterFinder::NewEvent() {
221 // Cleaning up and preparation for the clustering procedure
224 for (Int_t i=0; i<6; i++) {
225 fClusterList[i].Clear();
230 //___________________________________________________________________________________
231 void AliITSUpgradeClusterFinder::NewModule() {
237 memset(fHits,0,999*999*sizeof(Bool_t));
239 //___________________________________________________________________________________
240 Int_t AliITSUpgradeClusterFinder::DoModuleClustering(Int_t Layer, UShort_t charge) {
242 // Clustering and cluster-list container filling
245 UInt_t maxHits=fNhitsLeft;
246 for (UInt_t hit=0; hit<maxHits; hit++) {
247 if (fClusterTypeFlag) fFindClusterType = kTRUE;
249 fClusterWidthMinCol = fHitCol[hit];
250 fClusterWidthMinRow = fHitRow[hit];
251 fClusterWidthMaxCol = fHitCol[hit];
252 fClusterWidthMaxRow = fHitRow[hit];
255 fClusterTypeOrigCol = fHitCol[hit];
256 fClusterTypeOrigRow = fHitRow[hit];
257 memset(fClusterTypeArea,0,kMAXCLUSTERTYPESIDEZ*kMAXCLUSTERTYPESIDEY*sizeof(Bool_t));
260 UInt_t size = FindClusterRecu(fClusterTypeOrigCol,fClusterTypeOrigRow,charge);
263 fCharge=GetPixelCharge(fColSum,fRowSum);
264 AddLabelIndex(fColSum,fRowSum);
267 if(size>1) AliInfo(Form("DoModuleClustering, size %i , labels : %i %i %i \n",size,fLabels[0],fLabels[1],fLabels[2]));
268 fClusterList[Layer].Insert((Float_t)fColSum/size, (Float_t)fRowSum/size, size, GetClusterWidthZ(), GetClusterWidthPhi(), GetClusterType(size), fCharge,fLabels);
270 for(Int_t i=0; i<10; i++) fLabels[i]=-5;
272 if (fNhitsLeft==0) break;
276 //___________________________________________________________________________________
277 UInt_t AliITSUpgradeClusterFinder::FindClusterRecu(Int_t col, Int_t row, UShort_t charge) {
279 // Pixel selection for the clusters (adjiacent pixels or pixels at the corners)
282 if (col<0 || !fHits[col][row]) {return 0;}
283 fHits[col][row]=kFALSE;
289 if (col>fClusterWidthMaxCol) fClusterWidthMaxCol = col;
290 if (col<fClusterWidthMinCol) fClusterWidthMinCol = col;
291 if (row>fClusterWidthMaxRow) fClusterWidthMaxRow = row;
292 if (row<fClusterWidthMaxRow) fClusterWidthMinRow = row;
294 if (fFindClusterType) {
295 Short_t diffz = fClusterTypeOrigCol - col;
296 Short_t diffy = row - fClusterTypeOrigRow;
297 if (diffz>=kMAXCLUSTERTYPESIDEZ || diffy>=kMAXCLUSTERTYPESIDEY) {
298 fFindClusterType=kFALSE;
302 ShiftClusterTypeArea(kSHIFTRIGHT);
306 ShiftClusterTypeArea(kSHIFTDOWN);
309 fClusterTypeArea[diffz][diffy] = kTRUE;
313 size+=FindClusterRecu(col ,row-1,charge);
314 fCharge+=GetPixelCharge(col,row-1);
315 AddLabelIndex(col,row-1);
317 size+=FindClusterRecu(col-1,row ,charge);
318 fCharge+=GetPixelCharge(col-1,row);
319 AddLabelIndex(col-1,row);
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);
330 size+=FindClusterRecu(col-1,row-1,charge);
331 fCharge+=GetPixelCharge(col-1,row-1);
332 AddLabelIndex(col-1,row-1);
334 size+=FindClusterRecu(col-1,row+1,charge);
335 fCharge+=GetPixelCharge(col-1,row+1);
336 AddLabelIndex(col-1,row+1);
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);
348 //___________________________________________________________________________________
349 void AliITSUpgradeClusterFinder::ShiftClusterTypeArea(UInt_t direction) {
354 if (direction == kSHIFTRIGHT) {
355 fClusterTypeOrigCol++;
356 Bool_t tmpClusterTypeArea[kMAXCLUSTERTYPESIDEZ][kMAXCLUSTERTYPESIDEY];
357 memset(tmpClusterTypeArea,0,kMAXCLUSTERTYPESIDEZ*kMAXCLUSTERTYPESIDEY*sizeof(Bool_t));
358 for (UInt_t z=0; z<kMAXCLUSTERTYPESIDEZ; z++) {
359 for (UInt_t y=0; y<kMAXCLUSTERTYPESIDEY; y++) {
360 if (fClusterTypeArea[z][y]) {
361 if (z==kMAXCLUSTERTYPESIDEZ-1) {
362 fFindClusterType=kFALSE;
366 tmpClusterTypeArea[z+1][y] = kTRUE;
371 memcpy(fClusterTypeArea,tmpClusterTypeArea,kMAXCLUSTERTYPESIDEZ*kMAXCLUSTERTYPESIDEY*sizeof(Bool_t));
373 else if (direction == kSHIFTDOWN) {
374 fClusterTypeOrigRow--;
375 Bool_t tmpClusterTypeArea[kMAXCLUSTERTYPESIDEZ][kMAXCLUSTERTYPESIDEY];
376 memset(tmpClusterTypeArea,0,kMAXCLUSTERTYPESIDEZ*kMAXCLUSTERTYPESIDEY*sizeof(Bool_t));
377 for (UInt_t z=0; z<kMAXCLUSTERTYPESIDEZ; z++) {
378 for (UInt_t y=0; y<kMAXCLUSTERTYPESIDEY; y++) {
379 if (fClusterTypeArea[z][y]) {
380 if (y==kMAXCLUSTERTYPESIDEY-1) {
381 fFindClusterType=kFALSE;
385 tmpClusterTypeArea[z][y+1] = kTRUE;
390 memcpy(fClusterTypeArea,tmpClusterTypeArea,kMAXCLUSTERTYPESIDEZ*kMAXCLUSTERTYPESIDEY*sizeof(Bool_t));
393 //___________________________________________________________________________________
394 UShort_t AliITSUpgradeClusterFinder::GetCharge(Int_t layer,UInt_t index ) {
395 return fClusterList[layer].GetCharge(index);
397 //___________________________________________________________________________________
398 Int_t * AliITSUpgradeClusterFinder::GetLabels(Int_t layer,UInt_t index) {
399 return fClusterList[layer].GetLabels(index);
402 //___________________________________________________________________________________
403 UInt_t AliITSUpgradeClusterFinder::GetClusterWidthZ() {
404 return fClusterWidthMaxCol-fClusterWidthMinCol+1;
406 //___________________________________________________________________________________
407 UInt_t AliITSUpgradeClusterFinder::GetClusterWidthPhi() {
408 return fClusterWidthMaxRow-fClusterWidthMinRow+1;
410 //___________________________________________________________________________________
411 UInt_t AliITSUpgradeClusterFinder::GetClusterType(UInt_t size) {
415 if (!fFindClusterType || size>4)
425 fClusterTypeArea[0][1] &&
426 fClusterTypeArea[0][0] )
431 fClusterTypeArea[1][0] &&
432 fClusterTypeArea[0][0] )
439 fClusterTypeArea[0][2] &&
440 fClusterTypeArea[0][1] &&
441 fClusterTypeArea[0][0] )
445 // X and equivalent...
449 fClusterTypeArea[0][0] &&
450 fClusterTypeArea[0][1] &&
451 (fClusterTypeArea[1][1] ||
452 fClusterTypeArea[1][0])
456 fClusterTypeArea[1][0] &&
457 fClusterTypeArea[1][1] &&
458 (fClusterTypeArea[0][1] ||
459 fClusterTypeArea[0][0])
467 fClusterTypeArea[2][0] &&
468 fClusterTypeArea[1][0] &&
469 fClusterTypeArea[0][0] )
477 fClusterTypeArea[0][3] &&
478 fClusterTypeArea[0][2] &&
479 fClusterTypeArea[0][1] &&
480 fClusterTypeArea[0][0] )
486 fClusterTypeArea[1][1] &&
487 fClusterTypeArea[1][0] &&
488 fClusterTypeArea[0][1] &&
489 fClusterTypeArea[0][0] )
494 // X and equivalent...
498 fClusterTypeArea[1][0] &&
499 fClusterTypeArea[1][1] &&
500 fClusterTypeArea[1][2] &&
501 fClusterTypeArea[0][0]
505 fClusterTypeArea[1][0] &&
506 fClusterTypeArea[1][1] &&
507 fClusterTypeArea[1][2] &&
508 fClusterTypeArea[0][2]
512 fClusterTypeArea[0][0] &&
513 fClusterTypeArea[0][1] &&
514 fClusterTypeArea[0][2] &&
515 fClusterTypeArea[1][0]
519 fClusterTypeArea[0][0] &&
520 fClusterTypeArea[0][1] &&
521 fClusterTypeArea[0][2] &&
522 fClusterTypeArea[1][2]
530 // X and equivalent...
534 fClusterTypeArea[1][0] &&
535 fClusterTypeArea[1][1] &&
536 fClusterTypeArea[1][2] &&
537 fClusterTypeArea[0][1]
541 fClusterTypeArea[0][0] &&
542 fClusterTypeArea[0][1] &&
543 fClusterTypeArea[0][2] &&
544 fClusterTypeArea[1][1]
551 // XXX and equivalent...
555 fClusterTypeArea[2][0] &&
556 fClusterTypeArea[2][1] &&
557 fClusterTypeArea[1][1] &&
558 fClusterTypeArea[0][1]
562 fClusterTypeArea[0][0] &&
563 fClusterTypeArea[2][1] &&
564 fClusterTypeArea[1][1] &&
565 fClusterTypeArea[0][1]
569 fClusterTypeArea[2][1] &&
570 fClusterTypeArea[2][0] &&
571 fClusterTypeArea[1][0] &&
572 fClusterTypeArea[0][0]
576 fClusterTypeArea[0][1] &&
577 fClusterTypeArea[2][0] &&
578 fClusterTypeArea[1][0] &&
579 fClusterTypeArea[0][0]
586 // XXX and equivalent...
590 fClusterTypeArea[1][0] &&
591 fClusterTypeArea[2][1] &&
592 fClusterTypeArea[1][1] &&
593 fClusterTypeArea[0][1]
597 fClusterTypeArea[1][1] &&
598 fClusterTypeArea[2][0] &&
599 fClusterTypeArea[1][0] &&
600 fClusterTypeArea[0][0]
607 // X and equivalent...
611 fClusterTypeArea[0][0] &&
612 fClusterTypeArea[1][1]
616 fClusterTypeArea[1][0] &&
617 fClusterTypeArea[0][1]
625 // X and equivalent...
629 fClusterTypeArea[0][0] &&
630 fClusterTypeArea[0][1] &&
631 fClusterTypeArea[1][2]
635 fClusterTypeArea[1][0] &&
636 fClusterTypeArea[1][1] &&
637 fClusterTypeArea[0][2]
641 fClusterTypeArea[1][0] &&
642 fClusterTypeArea[0][1] &&
643 fClusterTypeArea[0][2]
647 fClusterTypeArea[0][0] &&
648 fClusterTypeArea[1][1] &&
649 fClusterTypeArea[1][2]
656 // X and equivalent...
660 fClusterTypeArea[0][0] &&
661 fClusterTypeArea[1][0] &&
662 fClusterTypeArea[2][1]
666 fClusterTypeArea[0][1] &&
667 fClusterTypeArea[1][1] &&
668 fClusterTypeArea[2][0]
672 fClusterTypeArea[0][0] &&
673 fClusterTypeArea[1][1] &&
674 fClusterTypeArea[2][1]
678 fClusterTypeArea[0][1] &&
679 fClusterTypeArea[1][0] &&
680 fClusterTypeArea[2][0]
691 //___________________________________________________________________________________________________
692 UInt_t AliITSUpgradeClusterFinder::GetPixelCharge(UInt_t col, UInt_t row){
697 for(Int_t entry =0; entry < fChargeArray->GetEntries(); entry++){
698 TObjString *s = (TObjString*)fChargeArray->At(entry);
699 TString name = s->GetString();
700 if(!name.Contains(Form("%i %i",col,row)))
702 TObjArray *array = name.Tokenize(" ");
703 TString charge = ((TObjString*)array->At(2))->String();
705 rowS = ((TObjString*)array->At(0))->String();
706 colS = ((TObjString*)array->At(1))->String();
714 //____________________________________________________
716 void AliITSUpgradeClusterFinder::AddLabelIndex(UInt_t col, UInt_t row){
718 // Adding cluster labels
721 for(Int_t entry =0; entry < fChargeArray->GetEntries(); entry++){
722 TObjString *s = (TObjString*)fChargeArray->At(entry);
723 TString name = s->GetString();
724 if(!name.Contains(Form("%i %i",col,row)))
726 TObjArray *array = name.Tokenize(" ");
729 for(Int_t i=0; i<3; i++){
730 index[i]=((TObjString*)array->At(3+i))->String();
731 label[i]=index[i].Atoi();
738 //____________________________________________________
740 void AliITSUpgradeClusterFinder::SetLabels(Int_t label[3]){
743 //Set the MC lables to the cluster
745 Bool_t isAssigned[3]={kFALSE,kFALSE,kFALSE};
747 for(Int_t i=0; i<10; i++){
748 for(Int_t k=0; k<3; k++){
749 if(fLabels[i]!=label[k] && label[k]>-1 && fLabels[i]<0) {
751 // printf("assignign...\n.");
752 // printf("Assignign fLabels[%i]=%i label[%i]=%i \n",i,fLabels[i],k,label[k]);
753 fLabels[i+k]=label[k];
761 //____________________________________________________
762 void AliITSUpgradeClusterFinder::MakeRecPointBranch(TTree *treeR){
764 // Creating the branch (see AliITSUpgradeReconstructor::Reconstruct)
767 if(!fRecPoints)fRecPoints = new TClonesArray("AliITSRecPoint",1000);
769 TBranch *branch = treeR->GetBranch("ITSRecPoints");
771 else branch = treeR->Branch("ITSRecPoints",&fRecPoints);
774 //____________________________________________________
775 void AliITSUpgradeClusterFinder::SetRecPointTreeAddress(TTree *treeR){
777 // Addressing the branch (see AliITSUpgradeReconstructor::Reconstruct)
780 if(!fRecPoints) fRecPoints = new TClonesArray("AliITSRecPoint",1000);
783 branch = treeR->GetBranch("ITSRecPoints");
785 branch->SetAddress(&fRecPoints);
786 } else AliError("No RecPoint branch available");
789 //____________________________________________________
790 void AliITSUpgradeClusterFinder::DigitsToRecPoints(const TObjArray *digList) {
792 // the clusterization is performed here
794 AliITSsegmentationUpgrade *segmentation = new AliITSsegmentationUpgrade();
795 AliITSRecPoint recpnt;
797 TClonesArray &lrecp = *fRecPoints;
799 for(Int_t ilayer=0; ilayer < 6 ;ilayer ++){
800 TClonesArray *pArrDig= (TClonesArray*)digList->At(ilayer);
802 AliDebug(1,Form("layer %i : # digits %i",ilayer,pArrDig->GetEntries()));
803 for(Int_t ientr =0; ientr < pArrDig->GetEntries() ; ientr++){
804 AliITSDigitUpgrade *dig = (AliITSDigitUpgrade*)pArrDig->At(ientr);
805 Int_t colz=dig->GetzPixelNumber();
806 Int_t rowx=dig->GetxPixelNumber();
807 Double_t hitcharge= (dig->GetNelectrons());
808 ProcessHit(ilayer,colz, rowx,(Short_t)hitcharge,dig->GetTracks());
812 for(UInt_t nClu = 0; nClu < GetClusterCount(ilayer); nClu++){
813 UShort_t charge = GetCharge(ilayer, nClu);
815 recpnt.SetLayer(ilayer);
816 Int_t *lab=GetLabels(ilayer,nClu);
817 for(Int_t l=0; l<3; l++) {recpnt.SetLabel(lab[l],l);}
823 Double_t xzl2[2]={0.,0.};
824 Double_t xPixC2,zPixC2=0.;
826 xPixC2 = GetClusterMeanRow(ilayer, nClu);
827 zPixC2 = GetClusterMeanCol(ilayer, nClu);
828 xzl2[0] = xPixC2*(segmentation->GetCellSizeX(ilayer))+0.5*(segmentation-> GetCellSizeX(ilayer));
829 xzl2[1] = zPixC2*(segmentation->GetCellSizeZ(ilayer))+0.5*(segmentation->GetCellSizeZ(ilayer))-(segmentation->GetHalfLength(ilayer));
830 check2 = segmentation->DetToGlobal(ilayer,xzl2[0], xzl2[1],xcheck2,ycheck2,zcheck2);
831 recpnt.SetType(GetClusterType(ilayer,nClu ));
832 // recpnt.SetLocalCoord(xzl2[0],xzl2[1]); //temporary solution (no LocalToTrack Matrix)
833 //from global to tracking system coordinate
834 // global coordinate -> local coordinate getting alpha angle of the recpoint
835 Float_t xclg = xcheck2;//upgrade clusters global coordinate ( ITS official: GetX tracking coordinate)
836 Float_t yclg = ycheck2;
837 Float_t zclg = zcheck2;
838 Double_t phiclu1rad, phiclu1deg;
839 phiclu1rad=TMath::ATan2(yclg,xclg);//cluster phi angle (rad)
840 if (phiclu1rad<0) phiclu1rad+=TMath::TwoPi();//from 0 to 360
841 else if (phiclu1rad>=TMath::TwoPi()) phiclu1rad-=TMath::TwoPi();//
843 phiclu1deg=180.*phiclu1rad/TMath::Pi();// in deg
844 Int_t ladder;// virtual segmentation starting from the cluster phi
846 ladder=(Int_t)phiclu1deg/18;// in which ladder the cluster is
847 Double_t alpha= (ladder*18.+9.)*TMath::Pi()/180.;//angle at the center of the ladder (rad)
850 Float_t xclu1Tr = xclg*TMath::Cos(alpha)-yclg*TMath::Sin(alpha);
851 Float_t yclu1 = yclg*TMath::Cos(alpha)+xclg*TMath::Sin(alpha);
852 Float_t xclu1 = TMath::Sqrt(xclu1Tr*xclu1Tr+yclu1*yclu1);
853 Float_t zclu1 = zclg;
854 Double_t phiTrk= (phiclu1rad-alpha);// cluster angle in the rotated system (rad)
856 yclu1=xclu1*phiTrk; // tracking system coordinate: r*phi
861 Double_t xsize, zsize;
862 segmentation->GetSegmentation(ilayer,xsize, zsize);
863 recpnt.SetSigmaY2(xsize/TMath::Sqrt(12)*xsize/TMath::Sqrt(12));
864 recpnt.SetSigmaZ2(zsize/TMath::Sqrt(12)*zsize/TMath::Sqrt(12));
865 new(lrecp[nClusters++]) AliITSRecPoint(recpnt);
866 //Int_t idx = fRecPoints->GetEntries();
867 AliInfo(Form("recpoint : Nelectrons %f (entry %i)",recpnt.GetQ(),fRecPoints->GetEntries()));
868 }//cluster list entries
870 if(segmentation) delete segmentation;