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 **************************************************************************/
16 //_________________________________________________________________________
17 // Reconstructed Points for the EMCAL
18 // A RecPoint is a cluster of digits
21 //*-- Author: Yves Schutz (SUBATECH)
22 //*-- Author: Dmitri Peressounko (RRC KI & SUBATECH)
23 //*-- Author: Heather Gray (LBL) merged AliEMCALRecPoint and AliEMCALTowerRecPoint 02/04
25 // --- ROOT system ---
28 #include "TPaveText.h"
29 #include "TClonesArray.h"
31 #include "TGeoMatrix.h"
32 #include "TGeoManager.h"
33 #include "TGeoPhysicalNode.h"
36 // --- Standard library ---
37 #include <Riostream.h>
39 // --- AliRoot header files ---
40 //#include "AliGenerator.h"
44 #include "AliGeomManager.h"
45 #include "AliEMCALGeometry.h"
46 #include "AliEMCALHit.h"
47 #include "AliEMCALDigit.h"
48 #include "AliEMCALRecPoint.h"
49 #include "AliCaloCalibPedestal.h"
50 #include "AliEMCALGeoParams.h"
52 ClassImp(AliEMCALRecPoint)
54 //____________________________________________________________________________
55 AliEMCALRecPoint::AliEMCALRecPoint()
56 : AliCluster(), fGeomPtr(0),
57 fAmp(0), fIndexInList(-1), //to be set when the point is already stored
58 fGlobPos(0,0,0),fLocPos(0,0,0),
59 fMaxDigit(100), fMulDigit(0), fMaxTrack(200),
60 fMulTrack(0), fDigitsList(0), fTracksList(0),
61 fClusterType(-1), fCoreEnergy(0), fDispersion(0),
62 fEnergyList(0), fAbsIdList(0),
63 fTime(0.), fNExMax(0), fCoreRadius(10), //HG check this
64 fDETracksList(0), fMulParent(0), fMaxParent(0),
65 fParentsList(0), fDEParentsList(0), fSuperModuleNumber(0),
66 fDigitIndMax(-1), fDistToBadTower(-1), fSharedCluster(kFALSE)
69 fGeomPtr = AliEMCALGeometry::GetInstance();
76 //____________________________________________________________________________
77 AliEMCALRecPoint::AliEMCALRecPoint(const char *)
78 : AliCluster(), fGeomPtr(0),
79 fAmp(0), fIndexInList(-1), //to be set when the point is already stored
80 fGlobPos(0,0,0), fLocPos(0,0,0),
81 fMaxDigit(100), fMulDigit(0), fMaxTrack(1000), fMulTrack(0),
82 fDigitsList(new Int_t[fMaxDigit]), fTracksList(new Int_t[fMaxTrack]),
83 fClusterType(-1), fCoreEnergy(0), fDispersion(0),
84 fEnergyList(new Float_t[fMaxDigit]),
85 fAbsIdList(new Int_t[fMaxDigit]), fTime(-1.), fNExMax(0), fCoreRadius(10),
86 fDETracksList(new Float_t[fMaxTrack]), fMulParent(0), fMaxParent(1000),
87 fParentsList(new Int_t[fMaxParent]), fDEParentsList(new Float_t[fMaxParent]),
88 fSuperModuleNumber(0), fDigitIndMax(-1), fDistToBadTower(-1),fSharedCluster(kFALSE)
91 for (Int_t i = 0; i < fMaxTrack; i++)
93 for (Int_t i = 0; i < fMaxParent; i++) {
95 fDEParentsList[i] = 0;
98 fGeomPtr = AliEMCALGeometry::GetInstance();
103 //____________________________________________________________________________
104 AliEMCALRecPoint::AliEMCALRecPoint(const AliEMCALRecPoint & rp)
105 : AliCluster(rp), fGeomPtr(rp.fGeomPtr),
106 fAmp(rp.fAmp), fIndexInList(rp.fIndexInList),
107 fGlobPos(rp.fGlobPos),fLocPos(rp.fLocPos),
108 fMaxDigit(rp.fMaxDigit), fMulDigit(rp.fMulDigit),
109 fMaxTrack(rp.fMaxTrack), fMulTrack(rp.fMaxTrack),
110 fDigitsList(new Int_t[rp.fMaxDigit]), fTracksList(new Int_t[rp.fMaxTrack]),
111 fClusterType(rp.fClusterType), fCoreEnergy(rp.fCoreEnergy),
112 fDispersion(rp.fDispersion),
113 fEnergyList(new Float_t[rp.fMaxDigit]),
114 fAbsIdList(new Int_t[rp.fMaxDigit]), fTime(rp.fTime), fNExMax(rp.fNExMax),fCoreRadius(rp.fCoreRadius),
115 fDETracksList(new Float_t[rp.fMaxTrack]), fMulParent(rp.fMulParent),
116 fMaxParent(rp.fMaxParent), fParentsList(new Int_t[rp.fMaxParent]),
117 fDEParentsList(new Float_t[rp.fMaxParent]),
118 fSuperModuleNumber(rp.fSuperModuleNumber), fDigitIndMax(rp.fDigitIndMax),
119 fDistToBadTower(rp.fDistToBadTower), fSharedCluster(rp.fSharedCluster)
122 fLambda[0] = rp.fLambda[0];
123 fLambda[1] = rp.fLambda[1];
125 for(Int_t i = 0; i < rp.fMulDigit; i++) {
126 fEnergyList[i] = rp.fEnergyList[i];
127 fAbsIdList[i] = rp.fAbsIdList[i];
130 for(Int_t i = 0; i < rp.fMulTrack; i++) fDETracksList[i] = rp.fDETracksList[i];
132 for(Int_t i = 0; i < rp.fMulParent; i++) {
133 fParentsList[i] = rp.fParentsList[i];
134 fDEParentsList[i] = rp.fDEParentsList[i];
138 //____________________________________________________________________________
139 AliEMCALRecPoint::~AliEMCALRecPoint()
143 delete[] fEnergyList ;
145 delete[] fAbsIdList ;
147 delete[] fDETracksList;
149 delete[] fParentsList;
151 delete[] fDEParentsList;
153 delete [] fDigitsList ;
154 delete [] fTracksList ;
157 //____________________________________________________________________________
158 AliEMCALRecPoint& AliEMCALRecPoint::operator= (const AliEMCALRecPoint &rp)
160 // assignment operator
162 if(&rp == this) return *this;
164 fGeomPtr = rp.fGeomPtr;
166 fIndexInList = rp.fIndexInList;
167 fGlobPos = rp.fGlobPos;
168 fLocPos = rp.fLocPos;
169 fMaxDigit = rp.fMaxDigit;
170 fMulDigit = rp.fMulDigit;
171 fMaxTrack = rp.fMaxTrack;
172 fMulTrack = rp.fMulTrack;
174 if(fDigitsList) delete [] fDigitsList;
175 fDigitsList = new Int_t[rp.fMaxDigit];
176 if(fTracksList) delete [] fTracksList;
177 fTracksList = new Int_t[rp.fMaxTrack];
178 for(Int_t i = 0; i<fMaxDigit; i++) fDigitsList[i] = rp.fDigitsList[i];
179 for(Int_t i = 0; i<fMaxTrack; i++) fTracksList[i] = rp.fTracksList[i];
181 fClusterType = rp.fClusterType;
182 fCoreEnergy = rp.fCoreEnergy;
183 fDispersion = rp.fDispersion;
186 if(fEnergyList) delete [] fEnergyList;
187 fEnergyList = new Float_t[rp.fMaxDigit];
188 if(fAbsIdList) delete [] fAbsIdList;
189 fAbsIdList = new Int_t[rp.fMaxDigit];
190 for(Int_t i = 0; i<fMaxDigit; i++) {
191 fEnergyList[i] = rp.fEnergyList[i];
192 fAbsIdList[i] = rp.fAbsIdList[i];
196 fNExMax = rp.fNExMax;
197 fCoreRadius = rp.fCoreRadius;
199 if(fDETracksList) delete [] fDETracksList;
200 fDETracksList = new Float_t[rp.fMaxTrack];
201 for(Int_t i = 0; i < fMaxTrack; i++) fDETracksList[i] = rp.fDETracksList[i];
203 fMulParent = rp.fMulParent;
204 fMaxParent = rp.fMaxParent;
206 if(fParentsList) delete [] fParentsList;
207 fParentsList = new Int_t[rp.fMaxParent];
208 if(fDEParentsList) delete [] fDEParentsList;
209 fDEParentsList = new Float_t[rp.fMaxParent];
210 for(Int_t i = 0; i < fMaxParent; i++) {
211 fParentsList[i] = rp.fParentsList[i];
212 fDEParentsList[i] = rp.fDEParentsList[i];
215 fSuperModuleNumber = rp.fSuperModuleNumber;
216 fDigitIndMax = rp.fDigitIndMax;
218 fLambda[0] = rp.fLambda[0];
219 fLambda[1] = rp.fLambda[1];
221 fDistToBadTower = rp.fDistToBadTower;
222 fSharedCluster = rp.fSharedCluster;
228 //____________________________________________________________________________
229 void AliEMCALRecPoint::AddDigit(AliEMCALDigit & digit, const Float_t energy, const Bool_t shared)
231 // Adds a digit to the RecPoint
232 // and accumulates the total amplitude and the multiplicity
235 fEnergyList = new Float_t[fMaxDigit];
237 if(fAbsIdList == 0) {
238 fAbsIdList = new Int_t [fMaxDigit];
241 if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists
243 Int_t * tempo = new Int_t [fMaxDigit];
244 Float_t * tempoE = new Float_t[fMaxDigit];
245 Int_t * tempoId = new Int_t [fMaxDigit];
248 for ( index = 0 ; index < fMulDigit ; index++ ){
249 tempo [index] = fDigitsList[index] ;
250 tempoE [index] = fEnergyList[index] ;
251 tempoId[index] = fAbsIdList [index] ;
254 delete [] fDigitsList ;
255 delete [] fEnergyList ;
256 delete [] fAbsIdList ;
259 fEnergyList = tempoE;
260 fAbsIdList = tempoId;
263 fDigitsList[fMulDigit] = digit.GetIndexInList() ;
264 fEnergyList[fMulDigit] = energy ;
265 fAbsIdList [fMulDigit] = digit.GetId();
269 if(shared) fSharedCluster = kTRUE;
271 //____________________________________________________________________________
272 Bool_t AliEMCALRecPoint::AreNeighbours(AliEMCALDigit * digit1, AliEMCALDigit * digit2 ) const
274 // Tells if (true) or not (false) two digits are neighbours
275 // A neighbour is defined as being two digits which share a corner
276 // ONLY USED IN CASE OF UNFOLDING
278 Bool_t areNeighbours = kFALSE ;
279 Int_t nSupMod=0, nModule=0, nIphi=0, nIeta=0;
280 Int_t nSupMod1=0, nModule1=0, nIphi1=0, nIeta1=0;
281 Int_t relid1[2] , relid2[2] ; // ieta, iphi
282 Int_t rowdiff=0, coldiff=0;
284 areNeighbours = kFALSE ;
286 fGeomPtr->GetCellIndex(digit1->GetId(), nSupMod,nModule,nIphi,nIeta);
287 fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nModule,nIphi,nIeta, relid1[0],relid1[1]);
289 fGeomPtr->GetCellIndex(digit2->GetId(), nSupMod1,nModule1,nIphi1,nIeta1);
290 fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod1,nModule1,nIphi1,nIeta1, relid2[0],relid2[1]);
292 // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2-1
293 // C Side impair SM, nSupMod%2=1; A side pair SM nSupMod%2=0
295 //printf("Shared cluster in 2 SMs!\n");
297 // if(nSupMod1%2) relid1[1]+=AliEMCALGeoParams::fgkEMCALCols;//bad
298 // else relid2[1]+=AliEMCALGeoParams::fgkEMCALCols;//bad
299 if(nSupMod1%2) relid2[1]+=AliEMCALGeoParams::fgkEMCALCols;
300 else relid1[1]+=AliEMCALGeoParams::fgkEMCALCols;
303 rowdiff = TMath::Abs( relid1[0] - relid2[0] ) ;
304 coldiff = TMath::Abs( relid1[1] - relid2[1] ) ;
306 //if (( coldiff <= 1 ) && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0))
307 if ((coldiff + rowdiff == 1 ))
308 areNeighbours = kTRUE ;
310 return areNeighbours;
313 //____________________________________________________________________________
314 Int_t AliEMCALRecPoint::Compare(const TObject * obj) const
316 // Compares two RecPoints according to their position in the EMCAL modules
318 Float_t delta = 1 ; //Width of "Sorting row".
322 AliEMCALRecPoint * clu = (AliEMCALRecPoint *)obj ;
325 GetLocalPosition(locpos1);
327 clu->GetLocalPosition(locpos2);
329 Int_t rowdif = (Int_t)(TMath::Ceil(locpos1.X()/delta)-TMath::Ceil(locpos2.X()/delta)) ;
334 else if(locpos1.Y()>locpos2.Y())
342 //___________________________________________________________________________
343 void AliEMCALRecPoint::Draw(Option_t *option)
345 // Draw this AliEMCALRecPoint with its current attributes
350 //____________________________________________________________________________
351 void AliEMCALRecPoint::EvalAll(Float_t logWeight,TClonesArray * digits, const Bool_t justClusters)
353 // Evaluates cluster parameters
355 // First calculate the index of digit with maximum amplitude and get
356 // the supermodule number where it sits.
358 fDigitIndMax = GetMaximalEnergyIndex();
359 fSuperModuleNumber = fGeomPtr->GetSuperModuleNumber(GetAbsIdMaxDigit());
361 //Evaluate global and local position
362 EvalGlobalPosition(logWeight, digits) ;
363 EvalLocalPosition(logWeight, digits) ;
365 //Evaluate shower parameters
366 EvalElipsAxis(logWeight, digits) ;
367 EvalDispersion(logWeight, digits) ;
369 //EvalCoreEnergy(logWeight, digits);
371 EvalPrimaries(digits) ;
374 //Called last because it sets the global position of the cluster?
375 //Do not call it when recalculating clusters out of standard reconstruction
377 EvalLocal2TrackingCSTransform();
382 //____________________________________________________________________________
383 void AliEMCALRecPoint::EvalDispersion(Float_t logWeight, TClonesArray * digits)
385 // Calculates the dispersion of the shower at the origin of the RecPoint
386 // in cell units - Nov 16,2006
388 Double_t d = 0., wtot = 0., w = 0.;
389 Int_t iDigit=0, nstat=0;
390 AliEMCALDigit * digit=0;
392 // Calculates the dispersion in cell units
393 Double_t etai, phii, etaMean=0.0, phiMean=0.0;
394 int nSupMod=0, nModule=0, nIphi=0, nIeta=0;
396 // Calculate mean values
397 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
398 digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit]) ;
400 if (fAmp>0 && fEnergyList[iDigit]>0) {
401 fGeomPtr->GetCellIndex(digit->GetId(), nSupMod,nModule,nIphi,nIeta);
402 fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nModule,nIphi,nIeta, iphi,ieta);
404 // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2
405 // C Side impair SM, nSupMod%2=1; A side pair SM nSupMod%2=0
406 if(fSharedCluster && nSupMod%2) ieta+=AliEMCALGeoParams::fgkEMCALCols;
410 w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
422 } else AliError(Form("Wrong weight %f\n", wtot));
424 // Calculate dispersion
425 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
426 digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit]) ;
428 if (fAmp>0 && fEnergyList[iDigit]>0) {
429 fGeomPtr->GetCellIndex(digit->GetId(), nSupMod,nModule,nIphi,nIeta);
430 fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nModule,nIphi,nIeta, iphi,ieta);
432 // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2
433 // C Side impair SM, nSupMod%2=1; A side pair SM, nSupMod%2=0
434 if(fSharedCluster && nSupMod%2) ieta+=AliEMCALGeoParams::fgkEMCALCols;
438 w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
442 d += w*((etai-etaMean)*(etai-etaMean)+(phii-phiMean)*(phii-phiMean));
447 if ( wtot > 0 && nstat>1) d /= wtot ;
450 fDispersion = TMath::Sqrt(d) ;
451 //printf("AliEMCALRecPoint::EvalDispersion() : Dispersion %f \n",fDispersion);
454 //____________________________________________________________________________
455 void AliEMCALRecPoint::EvalDistanceToBadChannels(AliCaloCalibPedestal* caloped)
457 //For each EMC rec. point set the distance to the nearest bad channel.
458 //AliInfo(Form("%d bad channel(s) found.\n", caloped->GetDeadTowerCount()));
459 //It is done in cell units and not in global or local position as before (Sept 2010)
461 if(!caloped->GetDeadTowerCount()) return;
463 //Get channels map of the supermodule where the cluster is.
464 TH2D* hMap = caloped->GetDeadMap(fSuperModuleNumber);
467 Float_t minDist = 10000.;
469 Int_t nSupMod, nModule;
472 fDigitIndMax = GetMaximalEnergyIndex();
473 fGeomPtr->GetCellIndex(fAbsIdList[fDigitIndMax], nSupMod,nModule,nIphi,nIeta);
474 fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nModule,nIphi,nIeta, iphi,ieta);
476 //Loop on tower status map
477 for(Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++){
478 for(Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++){
479 //Check if tower is bad.
480 if(hMap->GetBinContent(icol,irow)==AliCaloCalibPedestal::kAlive) continue;
481 //printf("AliEMCALRecPoint::EvalDistanceToBadChannels() - Bad channel in SM %d, col %d, row %d\n",iSM,icol, irow);
483 dRrow=TMath::Abs(irow-iphi);
484 dReta=TMath::Abs(icol-ieta);
485 dist=TMath::Sqrt(dRrow*dRrow+dReta*dReta);
486 if(dist < minDist) minDist = dist;
491 //In case the cluster is shared by 2 SuperModules, need to check the map of the second Super Module
492 if (fSharedCluster) {
496 //The only possible combinations are (0,1), (2,3) ... (10,11)
497 if(fSuperModuleNumber%2) nSupMod2 = fSuperModuleNumber-1;
498 else nSupMod2 = fSuperModuleNumber+1;
499 hMap2 = caloped->GetDeadMap(nSupMod2);
501 //Loop on tower status map of second super module
502 for(Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++){
503 for(Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++){
504 //Check if tower is bad.
505 if(hMap2->GetBinContent(icol,irow)==AliCaloCalibPedestal::kAlive) continue;
506 //printf("AliEMCALRecPoint::EvalDistanceToBadChannels() - Bad channel in SM %d, col %d, row %d\n",iSM,icol, irow);
507 dRrow=TMath::Abs(irow-iphi);
509 if(fSuperModuleNumber%2) {
510 dReta=TMath::Abs(icol-(AliEMCALGeoParams::fgkEMCALCols+ieta));
513 dReta=TMath::Abs(AliEMCALGeoParams::fgkEMCALCols+icol-ieta);
516 dist=TMath::Sqrt(dRrow*dRrow+dReta*dReta);
517 if(dist < minDist) minDist = dist;
522 }// shared cluster in 2 SuperModules
524 fDistToBadTower = minDist;
525 //printf("AliEMCALRecPoint::EvalDistanceToBadChannel() - Distance to Bad is %f cm, shared cluster? %d \n",fDistToBadTower,fSharedCluster);
529 //____________________________________________________________________________
530 void AliEMCALRecPoint::EvalLocalPosition(Float_t logWeight, TClonesArray * digits)
532 // Calculates the center of gravity in the local EMCAL-module coordinates
533 // Info("Print", " logWeight %f : cluster energy %f ", logWeight, fAmp); // for testing
535 AliEMCALDigit * digit=0;
538 Double_t dist = TmaxInCm(Double_t(fAmp));
539 //Int_t idMax = GetAbsIdMaxDigit();// idMax is not used at all in RelPosCellInSModule, why use it?
541 Double_t clXYZ[3]={0.,0.,0.}, clRmsXYZ[3]={0.,0.,0.}, xyzi[3], wtot=0., w=0.;
543 //printf(" dist : %f e : %f \n", dist, fAmp);
544 for(Int_t iDigit=0; iDigit<fMulDigit; iDigit++) {
545 digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
548 AliError("No Digit!!");
552 fGeomPtr->RelPosCellInSModule(digit->GetId(), dist, xyzi[0], xyzi[1], xyzi[2]);
554 //Temporal patch, due to mapping problem, need to swap "y" in one of the 2 SM, although no effect in position calculation. GCB 05/2010
555 if(fSharedCluster && fSuperModuleNumber != fGeomPtr->GetSuperModuleNumber(digit->GetId())) xyzi[1]*=-1;
557 //printf("EvalLocalPosition Cell: Id %i, SM %i : dist %f Local x,y,z %f %f %f \n",
558 // digit->GetId(), fGeomPtr->GetSuperModuleNumber(digit->GetId()), dist, xyzi[0], xyzi[1], xyzi[2]);
560 if(logWeight > 0.0) w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ));
561 else w = fEnergyList[iDigit]; // just energy
566 for(i=0; i<3; i++ ) {
567 clXYZ[i] += (w*xyzi[i]);
568 clRmsXYZ[i] += (w*xyzi[i]*xyzi[i]);
572 // cout << " wtot " << wtot << endl;
574 // xRMS = TMath::Sqrt(x2m - xMean*xMean);
575 for(i=0; i<3; i++ ) {
578 clRmsXYZ[i] /= (wtot*wtot);
579 clRmsXYZ[i] = clRmsXYZ[i] - clXYZ[i]*clXYZ[i];
580 if(clRmsXYZ[i] > 0.0) {
581 clRmsXYZ[i] = TMath::Sqrt(clRmsXYZ[i]);
582 } else clRmsXYZ[i] = 0;
583 } else clRmsXYZ[i] = 0;
586 for(i=0; i<3; i++ ) {
587 clXYZ[i] = clRmsXYZ[i] = -1.;
591 // // Cluster of one single digit, smear the position to avoid discrete position
592 // // smear x and z with +- 3 cm to uniform (avoid discrete effects). Tower size is approx 6 cm.
593 // // Rndm generates a number in ]0,1]
594 // if (fMulDigit==1) {
595 // clXYZ[0] += fGeomPtr->GetPhiTileSize()*(0.5 - gRandom->Rndm());
596 // clXYZ[2] += fGeomPtr->GetEtaTileSize()*(0.5 - gRandom->Rndm());
599 //Set position in local vector
600 fLocPos.SetX(clXYZ[0]);
601 fLocPos.SetY(clXYZ[1]);
602 fLocPos.SetZ(clXYZ[2]);
605 printf("EvalLocalPosition Cluster: Local (x,y,z) = (%f,%f,%f) \n", fLocPos.X(), fLocPos.Y(), fLocPos.Z()) ;
610 //____________________________________________________________________________
611 void AliEMCALRecPoint::EvalGlobalPosition(Float_t logWeight, TClonesArray * digits)
613 // Calculates the center of gravity in the global ALICE coordinates
614 // Info("Print", " logWeight %f : cluster energy %f ", logWeight, fAmp); // for testing
616 AliEMCALDigit * digit=0;
619 Double_t dist = TmaxInCm(Double_t(fAmp));
620 //Int_t idMax = GetAbsIdMaxDigit();// idMax is not used at all in RelPosCellInSModule, why use it?
622 Double_t clXYZ[3]={0.,0.,0.}, clRmsXYZ[3]={0.,0.,0.}, lxyzi[3], xyzi[3], wtot=0., w=0.;
624 //printf(" dist : %f e : %f \n", dist, fAmp);
625 for(Int_t iDigit=0; iDigit<fMulDigit; iDigit++) {
626 digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
629 AliError("No Digit!!");
633 //Get the local coordinates of the cell
634 fGeomPtr->RelPosCellInSModule(digit->GetId(), dist, lxyzi[0], lxyzi[1], lxyzi[2]);
636 //Now get the global coordinate
637 fGeomPtr->GetGlobal(lxyzi,xyzi, fGeomPtr->GetSuperModuleNumber(digit->GetId()));
638 //TVector3 pos(xyzi[0], xyzi[1], xyzi[2]);
639 //printf("EvalGlobalPosition Cell: Id %i, SM %i : dist %f Local (x,y,z) = (%f %f %f), eta %f, phi%f \n",
640 // digit->GetId(), fGeomPtr->GetSuperModuleNumber(digit->GetId()),dist, xyzi[0], xyzi[1], xyzi[2],pos.Eta(),pos.Phi()*TMath::RadToDeg());
642 if(logWeight > 0.0) w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ));
643 else w = fEnergyList[iDigit]; // just energy
648 for(i=0; i<3; i++ ) {
649 clXYZ[i] += (w*xyzi[i]);
650 clRmsXYZ[i] += (w*xyzi[i]*xyzi[i]);
654 // cout << " wtot " << wtot << endl;
656 // xRMS = TMath::Sqrt(x2m - xMean*xMean);
657 for(i=0; i<3; i++ ) {
660 clRmsXYZ[i] /= (wtot*wtot);
661 clRmsXYZ[i] = clRmsXYZ[i] - clXYZ[i]*clXYZ[i];
662 if(clRmsXYZ[i] > 0.0) {
663 clRmsXYZ[i] = TMath::Sqrt(clRmsXYZ[i]);
664 } else clRmsXYZ[i] = 0;
665 } else clRmsXYZ[i] = 0;
668 for(i=0; i<3; i++ ) {
669 clXYZ[i] = clRmsXYZ[i] = -1.;
673 // // Cluster of one single digit, smear the position to avoid discrete position
674 // // smear x and z with +- 3 cm to uniform (avoid discrete effects). Tower size is approx 6 cm.
675 // // Rndm generates a number in ]0,1]
676 // if (fMulDigit==1) {
677 // clXYZ[0] += fGeomPtr->GetPhiTileSize()*(0.5 - gRandom->Rndm());
678 // clXYZ[2] += fGeomPtr->GetEtaTileSize()*(0.5 - gRandom->Rndm());
681 //Set position in global vector
682 fGlobPos.SetX(clXYZ[0]);
683 fGlobPos.SetY(clXYZ[1]);
684 fGlobPos.SetZ(clXYZ[2]);
687 printf("EvalGlobalPosition Cluster: (x ,y ,z) = (%f,%f,%f), eta %f,phi %f\n",
688 fGlobPos.X(), fGlobPos.Y(), fGlobPos.Z(),fGlobPos.Eta(),fGlobPos.Phi()*TMath::RadToDeg()) ;
691 //____________________________________________________________________________
692 void AliEMCALRecPoint::EvalLocalPositionFit(Double_t deff, Double_t logWeight,
693 Double_t phiSlope, TClonesArray * digits)
695 // Evaluates local position of clusters in SM
698 AliEMCALDigit *digit=0;
700 Double_t clXYZ[3]={0.,0.,0.}, clRmsXYZ[3]={0.,0.,0.}, xyzi[3], wtot=0., w=0.;
702 Double_t dist = TmaxInCm(Double_t(fAmp));
703 //Int_t idMax = GetAbsIdMaxDigit();// idMax is not used at all in RelPosCellInSModule, why use it?
705 for(Int_t iDigit=0; iDigit<digits->GetEntries(); iDigit++) {
706 digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
709 //fGeomPtr->RelPosCellInSModule(digit->GetId(), idMax, dist, xyzi[0], xyzi[1], xyzi[2]);
710 fGeomPtr->RelPosCellInSModule(digit->GetId(), dist, xyzi[0], xyzi[1], xyzi[2]);
712 if(logWeight > 0.0) w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ));
713 else w = fEnergyList[iDigit]; // just energy
718 for(i=0; i<3; i++ ) {
719 clXYZ[i] += (w*xyzi[i]);
720 clRmsXYZ[i] += (w*xyzi[i]*xyzi[i]);
723 }else AliError("Digit null");
725 // cout << " wtot " << wtot << endl;
727 // xRMS = TMath::Sqrt(x2m - xMean*xMean);
728 for(i=0; i<3; i++ ) {
731 clRmsXYZ[i] /= (wtot*wtot);
732 clRmsXYZ[i] = clRmsXYZ[i] - clXYZ[i]*clXYZ[i];
733 if(clRmsXYZ[i] > 0.0) {
734 clRmsXYZ[i] = TMath::Sqrt(clRmsXYZ[i]);
735 } else clRmsXYZ[i] = 0;
736 } else clRmsXYZ[i] = 0;
739 for(i=0; i<3; i++ ) {
740 clXYZ[i] = clRmsXYZ[i] = -1.;
744 if(phiSlope != 0.0 && logWeight > 0.0 && wtot) {
745 // Correction in phi direction (y - coords here); Aug 16;
746 // May be put to global level or seperate method
747 ycorr = clXYZ[1] * (1. + phiSlope);
748 //printf(" y %f : ycorr %f : slope %f \n", clXYZ[1], ycorr, phiSlope);
752 fLocPos.SetX(clXYZ[0]);
753 fLocPos.SetY(clXYZ[1]);
754 fLocPos.SetZ(clXYZ[2]);
757 // printf("EvalLocalPosition: eta,phi,r = %f,%f,%f", fLocPos.X(), fLocPos.Y(), fLocPos.Z()) ;
760 //_____________________________________________________________________________
761 Bool_t AliEMCALRecPoint::EvalLocalPosition2(TClonesArray * digits, TArrayD &ed)
763 // Evaluated local position of rec.point using digits
764 // and parametrisation of w0 and deff
765 //printf(" <I> AliEMCALRecPoint::EvalLocalPosition2() \n");
766 return AliEMCALRecPoint::EvalLocalPositionFromDigits(digits, ed, fLocPos);
769 //_____________________________________________________________________________
770 Bool_t AliEMCALRecPoint::EvalLocalPositionFromDigits(TClonesArray *digits, TArrayD &ed, TVector3 &locPos)
772 // Used when digits should be recalibrated
773 Double_t deff=0, w0=0, esum=0;
775 // AliEMCALDigit *digit;
777 if(ed.GetSize() && (digits->GetEntries()!=ed.GetSize())) return kFALSE;
779 // Calculate sum energy of digits
781 for(iDigit=0; iDigit<ed.GetSize(); iDigit++) esum += ed[iDigit];
783 GetDeffW0(esum, deff, w0);
785 return EvalLocalPositionFromDigits(esum, deff, w0, digits, ed, locPos);
788 //_____________________________________________________________________________
789 Bool_t AliEMCALRecPoint::EvalLocalPositionFromDigits(const Double_t esum, const Double_t deff, const Double_t w0, TClonesArray *digits, TArrayD &ed, TVector3 &locPos)
791 //Evaluate position of digits in supermodule.
792 AliEMCALDigit *digit=0;
795 Double_t clXYZ[3]={0.,0.,0.}, xyzi[3], wtot=0., w=0.;
796 //Int_t idMax = GetAbsIdMaxDigit();// idMax is not used at all in RelPosCellInSModule, why use it?
798 // Get pointer to EMCAL geometry
799 // (can't use fGeomPtr in static method)
800 AliEMCALGeometry* geo = AliEMCALGeometry::GetInstance();
802 for(Int_t iDigit=0; iDigit<digits->GetEntries(); iDigit++) {
803 digit = dynamic_cast<AliEMCALDigit *>(digits->At(iDigit));
805 //geo->RelPosCellInSModule(digit->GetId(), idMax, deff, xyzi[0], xyzi[1], xyzi[2]);
806 geo->RelPosCellInSModule(digit->GetId(), deff, xyzi[0], xyzi[1], xyzi[2]);
808 if(w0 > 0.0) w = TMath::Max( 0., w0 + TMath::Log(ed[iDigit] / esum));
809 else w = ed[iDigit]; // just energy
814 for(i=0; i<3; i++ ) {
815 clXYZ[i] += (w*xyzi[i]);
818 }else AliError("Digit null");
820 // cout << " wtot " << wtot << endl;
822 for(i=0; i<3; i++ ) {
825 locPos.SetX(clXYZ[0]);
826 locPos.SetY(clXYZ[1]);
827 locPos.SetZ(clXYZ[2]);
835 //_____________________________________________________________________________
836 void AliEMCALRecPoint::GetDeffW0(const Double_t esum , Double_t &deff, Double_t &w0)
840 // Applied for simulation data with threshold 3 adc
841 // Calculate efective distance (deff) and weigh parameter (w0)
842 // for coordinate calculation; 0.5 GeV < esum <100 GeV.
843 // Look to: http://rhic.physics.wayne.edu/~pavlinov/ALICE/SHISHKEBAB/RES/CALIB/GEOMCORR/deffandW0VaEgamma_2.gif
846 const Double_t kdp0=9.25147, kdp1=1.16700; // Hard coded now
847 const Double_t kwp0=4.83713, kwp1=-2.77970e-01, kwp2 = 4.41116;
849 // No extrapolation here
850 e = esum<0.5?0.5:esum;
853 deff = kdp0 + kdp1*TMath::Log(e);
854 w0 = kwp0 / (1. + TMath::Exp(kwp1*(e+kwp2)));
855 //printf("<I> AliEMCALRecPoint::GetDeffW0 esum %5.2f : deff %5.2f : w0 %5.2f \n", esum, deff, w0);
858 //______________________________________________________________________________
859 void AliEMCALRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits)
861 // This function calculates energy in the core,
862 // i.e. within a radius rad = fCoreEnergy around the center. Beyond this radius
863 // in accordance with shower profile the energy deposition
864 // should be less than 2%
865 // Unfinished - Nov 15,2006
866 // Distance is calculate in (phi,eta) units
868 AliEMCALDigit * digit = 0 ;
872 if (!fLocPos.Mag()) {
873 EvalLocalPosition(logWeight, digits);
876 Double_t phiPoint = fLocPos.Phi(), etaPoint = fLocPos.Eta();
877 Double_t eta, phi, distance;
878 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
879 digit = (AliEMCALDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
882 fGeomPtr->EtaPhiFromIndex(digit->GetId(),eta, phi) ;
883 phi = phi * TMath::DegToRad();
885 distance = TMath::Sqrt((eta-etaPoint)*(eta-etaPoint)+(phi-phiPoint)*(phi-phiPoint));
886 if(distance < fCoreRadius)
887 fCoreEnergy += fEnergyList[iDigit] ;
891 //____________________________________________________________________________
892 void AliEMCALRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
894 // Calculates the axis of the shower ellipsoid in eta and phi
897 TString gn(fGeomPtr->GetName());
906 AliEMCALDigit * digit = 0;
908 Double_t etai =0, phii=0, w=0;
909 int nSupMod=0, nModule=0, nIphi=0, nIeta=0;
911 for(Int_t iDigit=0; iDigit<fMulDigit; iDigit++) {
912 digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit]) ;
914 // Nov 15,2006 - use cell numbers as coordinates
915 // Copied for shish-kebab geometry, ieta,iphi is cast as double as eta,phi
916 // We can use the eta,phi(or coordinates) of cell
917 nSupMod = nModule = nIphi = nIeta = iphi = ieta = 0;
919 fGeomPtr->GetCellIndex(digit->GetId(), nSupMod,nModule,nIphi,nIeta);
920 fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nModule,nIphi,nIeta, iphi,ieta);
922 // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2
923 // C Side impair SM, nSupMod%2=1; A side pair SM, nSupMod%2=0
924 if(fSharedCluster && nSupMod%2) ieta+=AliEMCALGeoParams::fgkEMCALCols;
929 w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
930 // fAmp summed amplitude of digits, i.e. energy of recpoint
931 // Gives smaller value of lambda than log weight
932 // w = fEnergyList[iDigit] / fAmp; // Nov 16, 2006 - try just energy
934 dxx += w * etai * etai ;
936 dzz += w * phii * phii ;
939 dxz += w * etai * phii ;
954 fLambda[0] = 0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
956 fLambda[0] = TMath::Sqrt(fLambda[0]) ;
960 fLambda[1] = 0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
962 if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
963 fLambda[1] = TMath::Sqrt(fLambda[1]) ;
971 //printf("AliEMCALRecPoint::EvalElipsAxis() lambdas = %f,%f \n", fLambda[0],fLambda[1]) ;
975 //______________________________________________________________________________
976 void AliEMCALRecPoint::EvalPrimaries(TClonesArray * digits)
978 // Constructs the list of primary particles (tracks) which
979 // have contributed to this RecPoint and calculate deposited energy
982 AliEMCALDigit * digit =0;
983 Int_t * primArray = new Int_t[fMaxTrack] ;
984 memset(primArray,-1,sizeof(Int_t)*fMaxTrack);
985 Float_t * dEPrimArray = new Float_t[fMaxTrack] ;
986 memset(dEPrimArray,-1,sizeof(Int_t)*fMaxTrack);
989 for ( index = 0 ; index < GetDigitsMultiplicity() ; index++ ) { // all digits
990 digit = dynamic_cast<AliEMCALDigit *>(digits->At( fDigitsList[index] )) ;
992 AliError("No Digit!!");
996 Int_t nprimaries = digit->GetNprimary() ;
997 if ( nprimaries == 0 ) continue ;
999 for ( jndex = 0 ; jndex < nprimaries ; jndex++ ) { // all primaries in digit
1000 if ( fMulTrack > fMaxTrack ) {
1001 fMulTrack = fMaxTrack ;
1002 Error("EvalPrimaries", "increase fMaxTrack ") ;
1005 Int_t newPrimary = digit->GetPrimary(jndex+1);
1006 Float_t dEPrimary = digit->GetDEPrimary(jndex+1);
1008 Bool_t already = kFALSE ;
1009 for ( kndex = 0 ; kndex < fMulTrack ; kndex++ ) { //check if not already stored
1010 if ( newPrimary == primArray[kndex] ){
1012 dEPrimArray[kndex] += dEPrimary;
1016 if ( !already && (fMulTrack < fMaxTrack)) { // store it
1017 primArray[fMulTrack] = newPrimary ;
1018 dEPrimArray[fMulTrack] = dEPrimary ;
1021 } // all primaries in digit
1024 Int_t *sortIdx = new Int_t[fMulTrack];
1025 TMath::Sort(fMulTrack,dEPrimArray,sortIdx);
1026 for(index = 0; index < fMulTrack; index++) {
1027 fTracksList[index] = primArray[sortIdx[index]] ;
1028 fDETracksList[index] = dEPrimArray[sortIdx[index]] ;
1031 delete [] primArray ;
1032 delete [] dEPrimArray ;
1036 //______________________________________________________________________________
1037 void AliEMCALRecPoint::EvalParents(TClonesArray * digits)
1039 // Constructs the list of parent particles (tracks) which have contributed to this RecPoint
1041 AliEMCALDigit * digit=0 ;
1042 Int_t * parentArray = new Int_t[fMaxTrack] ;
1043 memset(parentArray,-1,sizeof(Int_t)*fMaxTrack);
1044 Float_t * dEParentArray = new Float_t[fMaxTrack] ;
1045 memset(dEParentArray,-1,sizeof(Int_t)*fMaxTrack);
1048 for ( index = 0 ; index < GetDigitsMultiplicity() ; index++ ) { // all digits
1049 if (fDigitsList[index] >= digits->GetEntries() || fDigitsList[index] < 0)
1050 AliError(Form("Trying to get invalid digit %d (idx in WriteRecPoint %d)",fDigitsList[index],index));
1051 digit = dynamic_cast<AliEMCALDigit *>(digits->At( fDigitsList[index] )) ;
1053 AliError("No Digit!!");
1057 Int_t nparents = digit->GetNiparent() ;
1058 if ( nparents == 0 ) continue ;
1061 for ( jndex = 0 ; jndex < nparents ; jndex++ ) { // all primaries in digit
1062 if ( fMulParent > fMaxParent ) {
1064 Error("EvalParents", "increase fMaxParent") ;
1067 Int_t newParent = digit->GetIparent(jndex+1) ;
1068 Float_t newdEParent = digit->GetDEParent(jndex+1) ;
1070 Bool_t already = kFALSE ;
1071 for ( kndex = 0 ; kndex < fMulParent ; kndex++ ) { //check if not already stored
1072 if ( newParent == parentArray[kndex] ){
1073 dEParentArray[kndex] += newdEParent;
1078 if ( !already && (fMulParent < fMaxParent)) { // store it
1079 parentArray[fMulParent] = newParent ;
1080 dEParentArray[fMulParent] = newdEParent ;
1083 } // all parents in digit
1087 Int_t *sortIdx = new Int_t[fMulParent];
1088 TMath::Sort(fMulParent,dEParentArray,sortIdx);
1089 for(index = 0; index < fMulParent; index++) {
1090 fParentsList[index] = parentArray[sortIdx[index]] ;
1091 fDEParentsList[index] = dEParentArray[sortIdx[index]] ;
1096 delete [] parentArray;
1097 delete [] dEParentArray;
1100 //____________________________________________________________________________
1101 void AliEMCALRecPoint::GetLocalPosition(TVector3 & lpos) const
1103 // returns the position of the cluster in the local reference system
1104 // of the sub-detector
1109 //____________________________________________________________________________
1110 void AliEMCALRecPoint::GetGlobalPosition(TVector3 & gpos) const
1112 // returns the position of the cluster in the global reference system of ALICE
1113 // These are now the Cartesian X, Y and Z
1114 // cout<<" geom "<<geom<<endl;
1115 // fGeomPtr->GetGlobal(fLocPos, gpos, fSuperModuleNumber);
1120 //____________________________________________________________________________
1121 //void AliEMCALRecPoint::GetGlobalPosition(TVector3 & gpos, TMatrixF & gmat) const
1123 // // returns the position of the cluster in the global reference system of ALICE
1124 // // These are now the Cartesian X, Y and Z
1125 // // cout<<" geom "<<geom<<endl;
1127 // //To be implemented
1128 // fGeomPtr->GetGlobalEMCAL(this, gpos, gmat);
1132 //_____________________________________________________________________________
1133 void AliEMCALRecPoint::EvalLocal2TrackingCSTransform()
1135 //Evaluates local to "tracking" c.s. transformation (B.P.).
1136 //All evaluations should be completed before calling for this
1138 //See ALICE PPR Chapter 5 p.18 for "tracking" c.s. definition,
1139 //or just ask Jouri Belikov. :)
1141 SetVolumeId(AliGeomManager::LayerToVolUID(AliGeomManager::kEMCAL,GetSuperModuleNumber()));
1143 const TGeoHMatrix* tr2loc = GetTracking2LocalMatrix();
1144 if(!tr2loc) AliFatal(Form("No Tracking2LocalMatrix found."));
1146 Double_t lxyz[3] = {fLocPos.X(),fLocPos.Y(),fLocPos.Z()};
1147 Double_t txyz[3] = {0,0,0};
1149 tr2loc->MasterToLocal(lxyz,txyz);
1150 SetX(txyz[0]); SetY(txyz[1]); SetZ(txyz[2]);
1152 if(AliLog::GetGlobalDebugLevel()>0) {
1153 TVector3 gpos; //TMatrixF gmat;
1154 //GetGlobalPosition(gpos,gmat); //Not doing anythin special, replace by next line.
1155 fGeomPtr->GetGlobal(fLocPos, gpos, GetSuperModuleNumber());
1159 AliInfo(Form("lCS-->(%.3f,%.3f,%.3f), tCS-->(%.3f,%.3f,%.3f), gCS-->(%.3f,%.3f,%.3f), gCScalc-\
1160 ->(%.3f,%.3f,%.3f), supermodule %d",
1161 fLocPos.X(),fLocPos.Y(),fLocPos.Z(),
1162 GetX(),GetY(),GetZ(),
1163 gpos.X(),gpos.Y(),gpos.Z(),
1164 gxyz[0],gxyz[1],gxyz[2],GetSuperModuleNumber()));
1169 //____________________________________________________________________________
1170 Float_t AliEMCALRecPoint::GetMaximalEnergy(void) const
1172 // Finds the maximum energy in the cluster
1174 Float_t menergy = 0. ;
1177 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
1179 if(fEnergyList[iDigit] > menergy)
1180 menergy = fEnergyList[iDigit] ;
1185 //____________________________________________________________________________
1186 Int_t AliEMCALRecPoint::GetMaximalEnergyIndex(void) const
1188 // Finds the maximum energy in the cluster
1190 Float_t menergy = 0. ;
1194 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
1196 if(fEnergyList[iDigit] > menergy){
1197 menergy = fEnergyList[iDigit] ;
1200 }//loop on cluster digits
1206 //____________________________________________________________________________
1207 Int_t AliEMCALRecPoint::GetMultiplicityAtLevel(Float_t H) const
1209 // Calculates the multiplicity of digits with energy larger than H*energy
1213 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
1215 if(fEnergyList[iDigit] > H * fAmp)
1221 //____________________________________________________________________________
1222 Int_t AliEMCALRecPoint::GetNumberOfLocalMax(AliEMCALDigit ** maxAt, Float_t * maxAtEnergy,
1223 Float_t locMaxCut,TClonesArray * digits) const
1225 // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
1226 // energy difference between two local maxima
1228 AliEMCALDigit * digit = 0;
1229 AliEMCALDigit * digitN = 0;
1234 for(iDigit = 0; iDigit < fMulDigit; iDigit++)
1235 maxAt[iDigit] = (AliEMCALDigit*) digits->At(fDigitsList[iDigit]) ;
1237 for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {
1239 digit = maxAt[iDigit] ;
1241 for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {
1242 if(iDigitN == iDigit) continue;//the same digit
1243 digitN = (AliEMCALDigit *) digits->At(fDigitsList[iDigitN]) ;
1245 if ( AreNeighbours(digit, digitN) ) {
1246 if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {
1247 maxAt[iDigitN] = 0 ;
1248 // but may be digit too is not local max ?
1249 if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut)
1253 // but may be digitN too is not local max ?
1254 if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut)
1255 maxAt[iDigitN] = 0 ;
1257 } // if Areneighbours
1263 for(iDigit = 0; iDigit < fMulDigit; iDigit++) {
1265 maxAt[iDigitN] = maxAt[iDigit] ;
1266 maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
1274 //____________________________________________________________________________
1275 Int_t AliEMCALRecPoint::GetPrimaryIndex() const
1277 // Get the primary track index in TreeK which deposits the most energy
1278 // in Digits which forms RecPoint.
1281 return fTracksList[0];
1285 //____________________________________________________________________________
1286 void AliEMCALRecPoint::EvalTime(TClonesArray * digits){
1287 // time is set to the time of the digit with the maximum energy
1291 for(Int_t idig=0; idig < fMulDigit; idig++){
1292 if(fEnergyList[idig] > maxE){
1293 maxE = fEnergyList[idig] ;
1297 fTime = ((AliEMCALDigit*) digits->At(fDigitsList[maxAt]))->GetTime() ;
1301 //______________________________________________________________________________
1302 void AliEMCALRecPoint::Paint(Option_t *)
1304 // Paint this ALiRecPoint as a TMarker with its current attributes
1306 TVector3 pos(0.,0.,0.) ;
1307 GetLocalPosition(pos) ;
1308 Coord_t x = pos.X() ;
1309 Coord_t y = pos.Z() ;
1310 Color_t markercolor = 1 ;
1311 Size_t markersize = 1.;
1312 Style_t markerstyle = 5 ;
1314 if (!gPad->IsBatch()) {
1315 gVirtualX->SetMarkerColor(markercolor) ;
1316 gVirtualX->SetMarkerSize (markersize) ;
1317 gVirtualX->SetMarkerStyle(markerstyle) ;
1319 gPad->SetAttMarkerPS(markercolor,markerstyle,markersize) ;
1320 gPad->PaintPolyMarker(1,&x,&y,"") ;
1323 //_____________________________________________________________________
1324 Double_t AliEMCALRecPoint::TmaxInCm(const Double_t e , const Int_t key)
1327 // key = 0(gamma, default)
1329 const Double_t ca = 4.82; // shower max parameter - first guess; ca=TMath::Log(1000./8.07)
1330 Double_t tmax = 0.; // position of electromagnetic shower max in cm
1332 Double_t x0 = 1.31; // radiation lenght (cm)
1333 //If old geometry in use
1334 if(!((fGeomPtr->GetEMCGeometry()->GetGeoName()).Contains("V1"))) x0 = 1.28;
1337 tmax = TMath::Log(e) + ca;
1338 if (key==0) tmax += 0.5;
1340 tmax *= x0; // convert to cm
1345 //______________________________________________________________________________
1346 Float_t AliEMCALRecPoint::EtaToTheta(Float_t arg) const
1348 //Converts Theta (Radians) to Eta(Radians)
1349 return (2.*TMath::ATan(TMath::Exp(-arg)));
1352 //______________________________________________________________________________
1353 Float_t AliEMCALRecPoint::ThetaToEta(Float_t arg) const
1355 //Converts Eta (Radians) to Theta(Radians)
1356 return (-1 * TMath::Log(TMath::Tan(0.5 * arg)));
1359 //____________________________________________________________________________
1360 void AliEMCALRecPoint::Print(Option_t *opt) const
1362 // Print the list of digits belonging to the cluster
1363 if(strlen(opt)==0) return;
1365 message = "AliEMCALRecPoint:\n" ;
1366 message += " digits # = " ;
1367 AliInfo(message.Data()) ;
1370 for(iDigit=0; iDigit<fMulDigit; iDigit++)
1371 printf(" %d ", fDigitsList[iDigit] ) ;
1374 AliInfo(" Energies = ") ;
1375 for(iDigit=0; iDigit<fMulDigit; iDigit++)
1376 printf(" %f ", fEnergyList[iDigit] ) ;
1379 AliInfo("\n Abs Ids = ") ;
1380 for(iDigit=0; iDigit<fMulDigit; iDigit++)
1381 printf(" %i ", fAbsIdList[iDigit] ) ;
1384 AliInfo(" Primaries ") ;
1385 for(iDigit = 0;iDigit < fMulTrack; iDigit++)
1386 printf(" %d ", fTracksList[iDigit]) ;
1388 printf("\n Local x %6.2f y %7.2f z %7.1f \n", fLocPos[0], fLocPos[1], fLocPos[2]);
1390 message = " ClusterType = %d" ;
1391 message += " Multiplicity = %d" ;
1392 message += " Cluster Energy = %f" ;
1393 message += " Core energy = %f" ;
1394 message += " Core radius = %f" ;
1395 message += " Number of primaries %d" ;
1396 message += " Stored at position %d" ;
1397 AliInfo(Form(message.Data(), fClusterType, fMulDigit, fAmp, fCoreEnergy, fCoreRadius, fMulTrack, GetIndexInList()) ) ;
1400 //___________________________________________________________
1401 Double_t AliEMCALRecPoint::GetPointEnergy() const
1403 //Returns energy ....
1405 for(int ic=0; ic<GetMultiplicity(); ic++) e += double(fEnergyList[ic]);