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), fTimeList(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]), fTimeList(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]), fTimeList(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 fTimeList[i] = rp.fTimeList[i];
128 fAbsIdList[i] = rp.fAbsIdList[i];
131 for(Int_t i = 0; i < rp.fMulTrack; i++) fDETracksList[i] = rp.fDETracksList[i];
133 for(Int_t i = 0; i < rp.fMulParent; i++) {
134 fParentsList[i] = rp.fParentsList[i];
135 fDEParentsList[i] = rp.fDEParentsList[i];
139 //____________________________________________________________________________
140 AliEMCALRecPoint::~AliEMCALRecPoint()
144 delete[] fEnergyList ;
148 delete[] fAbsIdList ;
150 delete[] fDETracksList;
152 delete[] fParentsList;
154 delete[] fDEParentsList;
156 delete [] fDigitsList ;
157 delete [] fTracksList ;
160 //____________________________________________________________________________
161 AliEMCALRecPoint& AliEMCALRecPoint::operator= (const AliEMCALRecPoint &rp)
163 // assignment operator
165 if(&rp == this) return *this;
167 fGeomPtr = rp.fGeomPtr;
169 fIndexInList = rp.fIndexInList;
170 fGlobPos = rp.fGlobPos;
171 fLocPos = rp.fLocPos;
172 fMaxDigit = rp.fMaxDigit;
173 fMulDigit = rp.fMulDigit;
174 fMaxTrack = rp.fMaxTrack;
175 fMulTrack = rp.fMaxTrack;
176 for(Int_t i = 0; i<fMaxDigit; i++) fDigitsList[i] = rp.fDigitsList[i];
177 for(Int_t i = 0; i<fMaxTrack; i++) fTracksList[i] = rp.fTracksList[i];
178 fClusterType = rp.fClusterType;
179 fCoreEnergy = rp.fCoreEnergy;
180 fDispersion = rp.fDispersion;
181 for(Int_t i = 0; i<fMaxDigit; i++) {
182 fEnergyList[i] = rp.fEnergyList[i];
183 fTimeList[i] = rp.fTimeList[i];
184 fAbsIdList[i] = rp.fAbsIdList[i];
187 fNExMax = rp.fNExMax;
188 fCoreRadius = rp.fCoreRadius;
189 for(Int_t i = 0; i < fMaxTrack; i++) fDETracksList[i] = rp.fDETracksList[i];
190 fMulParent = rp.fMulParent;
191 fMaxParent = rp.fMaxParent;
192 for(Int_t i = 0; i < fMaxParent; i++) {
193 fParentsList[i] = rp.fParentsList[i];
194 fDEParentsList[i] = rp.fDEParentsList[i];
196 fSuperModuleNumber = rp.fSuperModuleNumber;
197 fDigitIndMax = rp.fDigitIndMax;
199 fLambda[0] = rp.fLambda[0];
200 fLambda[1] = rp.fLambda[1];
202 fDistToBadTower = rp.fDistToBadTower;
203 fSharedCluster = rp.fSharedCluster;
209 //____________________________________________________________________________
210 void AliEMCALRecPoint::AddDigit(AliEMCALDigit & digit, const Float_t energy, const Bool_t shared)
212 // Adds a digit to the RecPoint
213 // and accumulates the total amplitude and the multiplicity
216 fEnergyList = new Float_t[fMaxDigit];
218 fTimeList = new Float_t[fMaxDigit];
219 if(fAbsIdList == 0) {
220 fAbsIdList = new Int_t [fMaxDigit];
223 if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists
225 Int_t * tempo = new Int_t [fMaxDigit];
226 Float_t * tempoE = new Float_t[fMaxDigit];
227 Float_t * tempoT = new Float_t[fMaxDigit];
228 Int_t * tempoId = new Int_t [fMaxDigit];
231 for ( index = 0 ; index < fMulDigit ; index++ ){
232 tempo [index] = fDigitsList[index] ;
233 tempoE [index] = fEnergyList[index] ;
234 tempoT [index] = fTimeList [index] ;
235 tempoId[index] = fAbsIdList [index] ;
238 delete [] fDigitsList ;
239 delete [] fEnergyList ;
240 delete [] fTimeList ;
241 delete [] fAbsIdList ;
244 fEnergyList = tempoE;
246 fAbsIdList = tempoId;
249 fDigitsList[fMulDigit] = digit.GetIndexInList() ;
250 fEnergyList[fMulDigit] = energy ;
251 fTimeList [fMulDigit] = digit.GetTime();
252 fAbsIdList [fMulDigit] = digit.GetId();
256 if(shared) fSharedCluster = kTRUE;
258 //____________________________________________________________________________
259 Bool_t AliEMCALRecPoint::AreNeighbours(AliEMCALDigit * digit1, AliEMCALDigit * digit2 ) const
261 // Tells if (true) or not (false) two digits are neighbours
262 // A neighbour is defined as being two digits which share a corner
263 // ONLY USED IN CASE OF UNFOLDING
265 Bool_t areNeighbours = kFALSE ;
266 Int_t nSupMod=0, nModule=0, nIphi=0, nIeta=0;
267 Int_t nSupMod1=0, nModule1=0, nIphi1=0, nIeta1=0;
268 Int_t relid1[2] , relid2[2] ; // ieta, iphi
269 Int_t rowdiff=0, coldiff=0;
271 areNeighbours = kFALSE ;
273 fGeomPtr->GetCellIndex(digit1->GetId(), nSupMod,nModule,nIphi,nIeta);
274 fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nModule,nIphi,nIeta, relid1[0],relid1[1]);
276 fGeomPtr->GetCellIndex(digit2->GetId(), nSupMod1,nModule1,nIphi1,nIeta1);
277 fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod1,nModule1,nIphi1,nIeta1, relid2[0],relid2[1]);
279 // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2-1
280 // C Side impair SM, nSupMod%2=1; A side pair SM nSupMod%2=0
282 if(nSupMod1%2) relid1[1]+=AliEMCALGeoParams::fgkEMCALCols;
283 else relid2[1]+=AliEMCALGeoParams::fgkEMCALCols;
286 rowdiff = TMath::Abs( relid1[0] - relid2[0] ) ;
287 coldiff = TMath::Abs( relid1[1] - relid2[1] ) ;
289 if (( coldiff <= 1 ) && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0))
290 areNeighbours = kTRUE ;
292 return areNeighbours;
295 //____________________________________________________________________________
296 Int_t AliEMCALRecPoint::Compare(const TObject * obj) const
298 // Compares two RecPoints according to their position in the EMCAL modules
300 Float_t delta = 1 ; //Width of "Sorting row".
304 AliEMCALRecPoint * clu = (AliEMCALRecPoint *)obj ;
307 GetLocalPosition(locpos1);
309 clu->GetLocalPosition(locpos2);
311 Int_t rowdif = (Int_t)(TMath::Ceil(locpos1.X()/delta)-TMath::Ceil(locpos2.X()/delta)) ;
316 else if(locpos1.Y()>locpos2.Y())
324 //___________________________________________________________________________
325 void AliEMCALRecPoint::Draw(Option_t *option)
327 // Draw this AliEMCALRecPoint with its current attributes
332 //____________________________________________________________________________
333 void AliEMCALRecPoint::EvalAll(Float_t logWeight,TClonesArray * digits, const Bool_t justClusters)
335 // Evaluates cluster parameters
337 // First calculate the index of digit with maximum amplitude and get
338 // the supermodule number where it sits.
340 fDigitIndMax = GetMaximalEnergyIndex();
341 fSuperModuleNumber = fGeomPtr->GetSuperModuleNumber(GetAbsIdMaxDigit());
343 //Evaluate global and local position
344 EvalGlobalPosition(logWeight, digits) ;
345 EvalLocalPosition(logWeight, digits) ;
347 //Evaluate shower parameters
348 EvalElipsAxis(logWeight, digits) ;
349 EvalDispersion(logWeight, digits) ;
351 //EvalCoreEnergy(logWeight, digits);
353 EvalPrimaries(digits) ;
356 //Called last because it sets the global position of the cluster?
357 //Do not call it when recalculating clusters out of standard reconstruction
359 EvalLocal2TrackingCSTransform();
364 //____________________________________________________________________________
365 void AliEMCALRecPoint::EvalDispersion(Float_t logWeight, TClonesArray * digits)
367 // Calculates the dispersion of the shower at the origin of the RecPoint
368 // in cell units - Nov 16,2006
370 Double_t d = 0., wtot = 0., w = 0.;
371 Int_t iDigit=0, nstat=0;
372 AliEMCALDigit * digit=0;
374 // Calculates the dispersion in cell units
375 Double_t etai, phii, etaMean=0.0, phiMean=0.0;
376 int nSupMod=0, nModule=0, nIphi=0, nIeta=0;
378 // Calculate mean values
379 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
380 digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit]) ;
382 if (fAmp>0 && fEnergyList[iDigit]>0) {
383 fGeomPtr->GetCellIndex(digit->GetId(), nSupMod,nModule,nIphi,nIeta);
384 fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nModule,nIphi,nIeta, iphi,ieta);
386 // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2
387 // C Side impair SM, nSupMod%2=1; A side pair SM nSupMod%2=0
388 if(fSharedCluster && nSupMod%2) ieta+=AliEMCALGeoParams::fgkEMCALCols;
392 w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
404 } else AliError(Form("Wrong weight %f\n", wtot));
406 // Calculate dispersion
407 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
408 digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit]) ;
410 if (fAmp>0 && fEnergyList[iDigit]>0) {
411 fGeomPtr->GetCellIndex(digit->GetId(), nSupMod,nModule,nIphi,nIeta);
412 fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nModule,nIphi,nIeta, iphi,ieta);
414 // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2
415 // C Side impair SM, nSupMod%2=1; A side pair SM, nSupMod%2=0
416 if(fSharedCluster && nSupMod%2) ieta+=AliEMCALGeoParams::fgkEMCALCols;
420 w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
424 d += w*((etai-etaMean)*(etai-etaMean)+(phii-phiMean)*(phii-phiMean));
429 if ( wtot > 0 && nstat>1) d /= wtot ;
432 fDispersion = TMath::Sqrt(d) ;
433 //printf("AliEMCALRecPoint::EvalDispersion() : Dispersion %f \n",fDispersion);
436 //____________________________________________________________________________
437 void AliEMCALRecPoint::EvalDistanceToBadChannels(AliCaloCalibPedestal* caloped)
439 //For each EMC rec. point set the distance to the nearest bad channel.
440 //AliInfo(Form("%d bad channel(s) found.\n", caloped->GetDeadTowerCount()));
441 //It is done in cell units and not in global or local position as before (Sept 2010)
443 if(!caloped->GetDeadTowerCount()) return;
445 //Get channels map of the supermodule where the cluster is.
446 TH2D* hMap = caloped->GetDeadMap(fSuperModuleNumber);
449 Float_t minDist = 10000.;
451 Int_t nSupMod, nModule;
454 fDigitIndMax = GetMaximalEnergyIndex();
455 fGeomPtr->GetCellIndex(fAbsIdList[fDigitIndMax], nSupMod,nModule,nIphi,nIeta);
456 fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nModule,nIphi,nIeta, iphi,ieta);
458 //Loop on tower status map
459 for(Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++){
460 for(Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++){
461 //Check if tower is bad.
462 if(hMap->GetBinContent(icol,irow)==AliCaloCalibPedestal::kAlive) continue;
463 //printf("AliEMCALRecPoint::EvalDistanceToBadChannels() - Bad channel in SM %d, col %d, row %d\n",iSM,icol, irow);
465 dRrow=TMath::Abs(irow-iphi);
466 dReta=TMath::Abs(icol-ieta);
467 dist=TMath::Sqrt(dRrow*dRrow+dReta*dReta);
468 if(dist < minDist) minDist = dist;
473 //In case the cluster is shared by 2 SuperModules, need to check the map of the second Super Module
474 if (fSharedCluster) {
478 //The only possible combinations are (0,1), (2,3) ... (10,11)
479 if(fSuperModuleNumber%2) nSupMod2 = fSuperModuleNumber-1;
480 else nSupMod2 = fSuperModuleNumber+1;
481 hMap2 = caloped->GetDeadMap(nSupMod2);
483 //Loop on tower status map of second super module
484 for(Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++){
485 for(Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++){
486 //Check if tower is bad.
487 if(hMap2->GetBinContent(icol,irow)==AliCaloCalibPedestal::kAlive) continue;
488 //printf("AliEMCALRecPoint::EvalDistanceToBadChannels() - Bad channel in SM %d, col %d, row %d\n",iSM,icol, irow);
489 dRrow=TMath::Abs(irow-iphi);
491 if(fSuperModuleNumber%2) {
492 dReta=TMath::Abs(icol-(AliEMCALGeoParams::fgkEMCALCols+ieta));
495 dReta=TMath::Abs(AliEMCALGeoParams::fgkEMCALCols+icol-ieta);
498 dist=TMath::Sqrt(dRrow*dRrow+dReta*dReta);
499 if(dist < minDist) minDist = dist;
504 }// shared cluster in 2 SuperModules
506 fDistToBadTower = minDist;
507 //printf("AliEMCALRecPoint::EvalDistanceToBadChannel() - Distance to Bad is %f cm, shared cluster? %d \n",fDistToBadTower,fSharedCluster);
511 //____________________________________________________________________________
512 void AliEMCALRecPoint::EvalLocalPosition(Float_t logWeight, TClonesArray * digits)
514 // Calculates the center of gravity in the local EMCAL-module coordinates
515 // Info("Print", " logWeight %f : cluster energy %f ", logWeight, fAmp); // for testing
517 AliEMCALDigit * digit=0;
520 Double_t dist = TmaxInCm(Double_t(fAmp));
521 //Int_t idMax = GetAbsIdMaxDigit();// idMax is not used at all in RelPosCellInSModule, why use it?
523 Double_t clXYZ[3]={0.,0.,0.}, clRmsXYZ[3]={0.,0.,0.}, xyzi[3], wtot=0., w=0.;
525 //printf(" dist : %f e : %f \n", dist, fAmp);
526 for(Int_t iDigit=0; iDigit<fMulDigit; iDigit++) {
527 digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
530 AliError("No Digit!!");
534 fGeomPtr->RelPosCellInSModule(digit->GetId(), dist, xyzi[0], xyzi[1], xyzi[2]);
536 //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
537 if(fSharedCluster && fSuperModuleNumber != fGeomPtr->GetSuperModuleNumber(digit->GetId())) xyzi[1]*=-1;
539 //printf("EvalLocalPosition Cell: Id %i, SM %i : dist %f Local x,y,z %f %f %f \n",
540 // digit->GetId(), fGeomPtr->GetSuperModuleNumber(digit->GetId()), dist, xyzi[0], xyzi[1], xyzi[2]);
542 if(logWeight > 0.0) w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ));
543 else w = fEnergyList[iDigit]; // just energy
548 for(i=0; i<3; i++ ) {
549 clXYZ[i] += (w*xyzi[i]);
550 clRmsXYZ[i] += (w*xyzi[i]*xyzi[i]);
554 // cout << " wtot " << wtot << endl;
556 // xRMS = TMath::Sqrt(x2m - xMean*xMean);
557 for(i=0; i<3; i++ ) {
560 clRmsXYZ[i] /= (wtot*wtot);
561 clRmsXYZ[i] = clRmsXYZ[i] - clXYZ[i]*clXYZ[i];
562 if(clRmsXYZ[i] > 0.0) {
563 clRmsXYZ[i] = TMath::Sqrt(clRmsXYZ[i]);
564 } else clRmsXYZ[i] = 0;
565 } else clRmsXYZ[i] = 0;
568 for(i=0; i<3; i++ ) {
569 clXYZ[i] = clRmsXYZ[i] = -1.;
573 // // Cluster of one single digit, smear the position to avoid discrete position
574 // // smear x and z with +- 3 cm to uniform (avoid discrete effects). Tower size is approx 6 cm.
575 // // Rndm generates a number in ]0,1]
576 // if (fMulDigit==1) {
577 // clXYZ[0] += fGeomPtr->GetPhiTileSize()*(0.5 - gRandom->Rndm());
578 // clXYZ[2] += fGeomPtr->GetEtaTileSize()*(0.5 - gRandom->Rndm());
581 //Set position in local vector
582 fLocPos.SetX(clXYZ[0]);
583 fLocPos.SetY(clXYZ[1]);
584 fLocPos.SetZ(clXYZ[2]);
587 printf("EvalLocalPosition Cluster: Local (x,y,z) = (%f,%f,%f) \n", fLocPos.X(), fLocPos.Y(), fLocPos.Z()) ;
592 //____________________________________________________________________________
593 void AliEMCALRecPoint::EvalGlobalPosition(Float_t logWeight, TClonesArray * digits)
595 // Calculates the center of gravity in the global ALICE coordinates
596 // Info("Print", " logWeight %f : cluster energy %f ", logWeight, fAmp); // for testing
598 AliEMCALDigit * digit=0;
601 Double_t dist = TmaxInCm(Double_t(fAmp));
602 //Int_t idMax = GetAbsIdMaxDigit();// idMax is not used at all in RelPosCellInSModule, why use it?
604 Double_t clXYZ[3]={0.,0.,0.}, clRmsXYZ[3]={0.,0.,0.}, lxyzi[3], xyzi[3], wtot=0., w=0.;
606 //printf(" dist : %f e : %f \n", dist, fAmp);
607 for(Int_t iDigit=0; iDigit<fMulDigit; iDigit++) {
608 digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
611 AliError("No Digit!!");
615 //Get the local coordinates of the cell
616 fGeomPtr->RelPosCellInSModule(digit->GetId(), dist, lxyzi[0], lxyzi[1], lxyzi[2]);
618 //Now get the global coordinate
619 fGeomPtr->GetGlobal(lxyzi,xyzi, fGeomPtr->GetSuperModuleNumber(digit->GetId()));
620 //TVector3 pos(xyzi[0], xyzi[1], xyzi[2]);
621 //printf("EvalGlobalPosition Cell: Id %i, SM %i : dist %f Local (x,y,z) = (%f %f %f), eta %f, phi%f \n",
622 // digit->GetId(), fGeomPtr->GetSuperModuleNumber(digit->GetId()),dist, xyzi[0], xyzi[1], xyzi[2],pos.Eta(),pos.Phi()*TMath::RadToDeg());
624 if(logWeight > 0.0) w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ));
625 else w = fEnergyList[iDigit]; // just energy
630 for(i=0; i<3; i++ ) {
631 clXYZ[i] += (w*xyzi[i]);
632 clRmsXYZ[i] += (w*xyzi[i]*xyzi[i]);
636 // cout << " wtot " << wtot << endl;
638 // xRMS = TMath::Sqrt(x2m - xMean*xMean);
639 for(i=0; i<3; i++ ) {
642 clRmsXYZ[i] /= (wtot*wtot);
643 clRmsXYZ[i] = clRmsXYZ[i] - clXYZ[i]*clXYZ[i];
644 if(clRmsXYZ[i] > 0.0) {
645 clRmsXYZ[i] = TMath::Sqrt(clRmsXYZ[i]);
646 } else clRmsXYZ[i] = 0;
647 } else clRmsXYZ[i] = 0;
650 for(i=0; i<3; i++ ) {
651 clXYZ[i] = clRmsXYZ[i] = -1.;
655 // // Cluster of one single digit, smear the position to avoid discrete position
656 // // smear x and z with +- 3 cm to uniform (avoid discrete effects). Tower size is approx 6 cm.
657 // // Rndm generates a number in ]0,1]
658 // if (fMulDigit==1) {
659 // clXYZ[0] += fGeomPtr->GetPhiTileSize()*(0.5 - gRandom->Rndm());
660 // clXYZ[2] += fGeomPtr->GetEtaTileSize()*(0.5 - gRandom->Rndm());
663 //Set position in global vector
664 fGlobPos.SetX(clXYZ[0]);
665 fGlobPos.SetY(clXYZ[1]);
666 fGlobPos.SetZ(clXYZ[2]);
669 printf("EvalGlobalPosition Cluster: (x ,y ,z) = (%f,%f,%f), eta %f,phi %f\n",
670 fGlobPos.X(), fGlobPos.Y(), fGlobPos.Z(),fGlobPos.Eta(),fGlobPos.Phi()*TMath::RadToDeg()) ;
673 //____________________________________________________________________________
674 void AliEMCALRecPoint::EvalLocalPositionFit(Double_t deff, Double_t logWeight,
675 Double_t phiSlope, TClonesArray * digits)
677 // Evaluates local position of clusters in SM
680 AliEMCALDigit *digit=0;
682 Double_t clXYZ[3]={0.,0.,0.}, clRmsXYZ[3]={0.,0.,0.}, xyzi[3], wtot=0., w=0.;
684 Double_t dist = TmaxInCm(Double_t(fAmp));
685 //Int_t idMax = GetAbsIdMaxDigit();// idMax is not used at all in RelPosCellInSModule, why use it?
687 for(Int_t iDigit=0; iDigit<digits->GetEntries(); iDigit++) {
688 digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
691 //fGeomPtr->RelPosCellInSModule(digit->GetId(), idMax, dist, xyzi[0], xyzi[1], xyzi[2]);
692 fGeomPtr->RelPosCellInSModule(digit->GetId(), dist, xyzi[0], xyzi[1], xyzi[2]);
694 if(logWeight > 0.0) w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ));
695 else w = fEnergyList[iDigit]; // just energy
700 for(i=0; i<3; i++ ) {
701 clXYZ[i] += (w*xyzi[i]);
702 clRmsXYZ[i] += (w*xyzi[i]*xyzi[i]);
705 }else AliError("Digit null");
707 // cout << " wtot " << wtot << endl;
709 // xRMS = TMath::Sqrt(x2m - xMean*xMean);
710 for(i=0; i<3; i++ ) {
713 clRmsXYZ[i] /= (wtot*wtot);
714 clRmsXYZ[i] = clRmsXYZ[i] - clXYZ[i]*clXYZ[i];
715 if(clRmsXYZ[i] > 0.0) {
716 clRmsXYZ[i] = TMath::Sqrt(clRmsXYZ[i]);
717 } else clRmsXYZ[i] = 0;
718 } else clRmsXYZ[i] = 0;
721 for(i=0; i<3; i++ ) {
722 clXYZ[i] = clRmsXYZ[i] = -1.;
726 if(phiSlope != 0.0 && logWeight > 0.0 && wtot) {
727 // Correction in phi direction (y - coords here); Aug 16;
728 // May be put to global level or seperate method
729 ycorr = clXYZ[1] * (1. + phiSlope);
730 //printf(" y %f : ycorr %f : slope %f \n", clXYZ[1], ycorr, phiSlope);
734 fLocPos.SetX(clXYZ[0]);
735 fLocPos.SetY(clXYZ[1]);
736 fLocPos.SetZ(clXYZ[2]);
739 // printf("EvalLocalPosition: eta,phi,r = %f,%f,%f", fLocPos.X(), fLocPos.Y(), fLocPos.Z()) ;
742 //_____________________________________________________________________________
743 Bool_t AliEMCALRecPoint::EvalLocalPosition2(TClonesArray * digits, TArrayD &ed)
745 // Evaluated local position of rec.point using digits
746 // and parametrisation of w0 and deff
747 //printf(" <I> AliEMCALRecPoint::EvalLocalPosition2() \n");
748 return AliEMCALRecPoint::EvalLocalPositionFromDigits(digits, ed, fLocPos);
751 //_____________________________________________________________________________
752 Bool_t AliEMCALRecPoint::EvalLocalPositionFromDigits(TClonesArray *digits, TArrayD &ed, TVector3 &locPos)
754 // Used when digits should be recalibrated
755 Double_t deff=0, w0=0, esum=0;
757 // AliEMCALDigit *digit;
759 if(ed.GetSize() && (digits->GetEntries()!=ed.GetSize())) return kFALSE;
761 // Calculate sum energy of digits
763 for(iDigit=0; iDigit<ed.GetSize(); iDigit++) esum += ed[iDigit];
765 GetDeffW0(esum, deff, w0);
767 return EvalLocalPositionFromDigits(esum, deff, w0, digits, ed, locPos);
770 //_____________________________________________________________________________
771 Bool_t AliEMCALRecPoint::EvalLocalPositionFromDigits(const Double_t esum, const Double_t deff, const Double_t w0, TClonesArray *digits, TArrayD &ed, TVector3 &locPos)
773 //Evaluate position of digits in supermodule.
774 AliEMCALDigit *digit=0;
777 Double_t clXYZ[3]={0.,0.,0.}, xyzi[3], wtot=0., w=0.;
778 //Int_t idMax = GetAbsIdMaxDigit();// idMax is not used at all in RelPosCellInSModule, why use it?
780 // Get pointer to EMCAL geometry
781 // (can't use fGeomPtr in static method)
782 AliEMCALGeometry* geo = AliEMCALGeometry::GetInstance();
784 for(Int_t iDigit=0; iDigit<digits->GetEntries(); iDigit++) {
785 digit = dynamic_cast<AliEMCALDigit *>(digits->At(iDigit));
787 //geo->RelPosCellInSModule(digit->GetId(), idMax, deff, xyzi[0], xyzi[1], xyzi[2]);
788 geo->RelPosCellInSModule(digit->GetId(), deff, xyzi[0], xyzi[1], xyzi[2]);
790 if(w0 > 0.0) w = TMath::Max( 0., w0 + TMath::Log(ed[iDigit] / esum));
791 else w = ed[iDigit]; // just energy
796 for(i=0; i<3; i++ ) {
797 clXYZ[i] += (w*xyzi[i]);
800 }else AliError("Digit null");
802 // cout << " wtot " << wtot << endl;
804 for(i=0; i<3; i++ ) {
807 locPos.SetX(clXYZ[0]);
808 locPos.SetY(clXYZ[1]);
809 locPos.SetZ(clXYZ[2]);
817 //_____________________________________________________________________________
818 void AliEMCALRecPoint::GetDeffW0(const Double_t esum , Double_t &deff, Double_t &w0)
822 // Applied for simulation data with threshold 3 adc
823 // Calculate efective distance (deff) and weigh parameter (w0)
824 // for coordinate calculation; 0.5 GeV < esum <100 GeV.
825 // Look to: http://rhic.physics.wayne.edu/~pavlinov/ALICE/SHISHKEBAB/RES/CALIB/GEOMCORR/deffandW0VaEgamma_2.gif
828 const Double_t kdp0=9.25147, kdp1=1.16700; // Hard coded now
829 const Double_t kwp0=4.83713, kwp1=-2.77970e-01, kwp2 = 4.41116;
831 // No extrapolation here
832 e = esum<0.5?0.5:esum;
835 deff = kdp0 + kdp1*TMath::Log(e);
836 w0 = kwp0 / (1. + TMath::Exp(kwp1*(e+kwp2)));
837 //printf("<I> AliEMCALRecPoint::GetDeffW0 esum %5.2f : deff %5.2f : w0 %5.2f \n", esum, deff, w0);
840 //______________________________________________________________________________
841 void AliEMCALRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits)
843 // This function calculates energy in the core,
844 // i.e. within a radius rad = fCoreEnergy around the center. Beyond this radius
845 // in accordance with shower profile the energy deposition
846 // should be less than 2%
847 // Unfinished - Nov 15,2006
848 // Distance is calculate in (phi,eta) units
850 AliEMCALDigit * digit = 0 ;
854 if (!fLocPos.Mag()) {
855 EvalLocalPosition(logWeight, digits);
858 Double_t phiPoint = fLocPos.Phi(), etaPoint = fLocPos.Eta();
859 Double_t eta, phi, distance;
860 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
861 digit = (AliEMCALDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
864 fGeomPtr->EtaPhiFromIndex(digit->GetId(),eta, phi) ;
865 phi = phi * TMath::DegToRad();
867 distance = TMath::Sqrt((eta-etaPoint)*(eta-etaPoint)+(phi-phiPoint)*(phi-phiPoint));
868 if(distance < fCoreRadius)
869 fCoreEnergy += fEnergyList[iDigit] ;
873 //____________________________________________________________________________
874 void AliEMCALRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
876 // Calculates the axis of the shower ellipsoid in eta and phi
879 TString gn(fGeomPtr->GetName());
888 AliEMCALDigit * digit = 0;
890 Double_t etai =0, phii=0, w=0;
891 int nSupMod=0, nModule=0, nIphi=0, nIeta=0;
893 for(Int_t iDigit=0; iDigit<fMulDigit; iDigit++) {
894 digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit]) ;
896 // Nov 15,2006 - use cell numbers as coordinates
897 // Copied for shish-kebab geometry, ieta,iphi is cast as double as eta,phi
898 // We can use the eta,phi(or coordinates) of cell
899 nSupMod = nModule = nIphi = nIeta = iphi = ieta = 0;
901 fGeomPtr->GetCellIndex(digit->GetId(), nSupMod,nModule,nIphi,nIeta);
902 fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nModule,nIphi,nIeta, iphi,ieta);
904 // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2
905 // C Side impair SM, nSupMod%2=1; A side pair SM, nSupMod%2=0
906 if(fSharedCluster && nSupMod%2) ieta+=AliEMCALGeoParams::fgkEMCALCols;
911 w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
912 // fAmp summed amplitude of digits, i.e. energy of recpoint
913 // Gives smaller value of lambda than log weight
914 // w = fEnergyList[iDigit] / fAmp; // Nov 16, 2006 - try just energy
916 dxx += w * etai * etai ;
918 dzz += w * phii * phii ;
921 dxz += w * etai * phii ;
936 fLambda[0] = 0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
938 fLambda[0] = TMath::Sqrt(fLambda[0]) ;
942 fLambda[1] = 0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
944 if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
945 fLambda[1] = TMath::Sqrt(fLambda[1]) ;
953 //printf("AliEMCALRecPoint::EvalElipsAxis() lambdas = %f,%f \n", fLambda[0],fLambda[1]) ;
957 //______________________________________________________________________________
958 void AliEMCALRecPoint::EvalPrimaries(TClonesArray * digits)
960 // Constructs the list of primary particles (tracks) which
961 // have contributed to this RecPoint and calculate deposited energy
964 AliEMCALDigit * digit =0;
965 Int_t * primArray = new Int_t[fMaxTrack] ;
966 memset(primArray,-1,sizeof(Int_t)*fMaxTrack);
967 Float_t * dEPrimArray = new Float_t[fMaxTrack] ;
968 memset(dEPrimArray,-1,sizeof(Int_t)*fMaxTrack);
971 for ( index = 0 ; index < GetDigitsMultiplicity() ; index++ ) { // all digits
972 digit = dynamic_cast<AliEMCALDigit *>(digits->At( fDigitsList[index] )) ;
974 AliError("No Digit!!");
978 Int_t nprimaries = digit->GetNprimary() ;
979 if ( nprimaries == 0 ) continue ;
981 for ( jndex = 0 ; jndex < nprimaries ; jndex++ ) { // all primaries in digit
982 if ( fMulTrack > fMaxTrack ) {
983 fMulTrack = fMaxTrack ;
984 Error("EvalPrimaries", "increase fMaxTrack ") ;
987 Int_t newPrimary = digit->GetPrimary(jndex+1);
988 Float_t dEPrimary = digit->GetDEPrimary(jndex+1);
990 Bool_t already = kFALSE ;
991 for ( kndex = 0 ; kndex < fMulTrack ; kndex++ ) { //check if not already stored
992 if ( newPrimary == primArray[kndex] ){
994 dEPrimArray[kndex] += dEPrimary;
998 if ( !already && (fMulTrack < fMaxTrack)) { // store it
999 primArray[fMulTrack] = newPrimary ;
1000 dEPrimArray[fMulTrack] = dEPrimary ;
1003 } // all primaries in digit
1006 Int_t *sortIdx = new Int_t[fMulTrack];
1007 TMath::Sort(fMulTrack,dEPrimArray,sortIdx);
1008 for(index = 0; index < fMulTrack; index++) {
1009 fTracksList[index] = primArray[sortIdx[index]] ;
1010 fDETracksList[index] = dEPrimArray[sortIdx[index]] ;
1013 delete [] primArray ;
1014 delete [] dEPrimArray ;
1018 //______________________________________________________________________________
1019 void AliEMCALRecPoint::EvalParents(TClonesArray * digits)
1021 // Constructs the list of parent particles (tracks) which have contributed to this RecPoint
1023 AliEMCALDigit * digit=0 ;
1024 Int_t * parentArray = new Int_t[fMaxTrack] ;
1025 memset(parentArray,-1,sizeof(Int_t)*fMaxTrack);
1026 Float_t * dEParentArray = new Float_t[fMaxTrack] ;
1027 memset(dEParentArray,-1,sizeof(Int_t)*fMaxTrack);
1030 for ( index = 0 ; index < GetDigitsMultiplicity() ; index++ ) { // all digits
1031 if (fDigitsList[index] >= digits->GetEntries() || fDigitsList[index] < 0)
1032 AliError(Form("Trying to get invalid digit %d (idx in WriteRecPoint %d)",fDigitsList[index],index));
1033 digit = dynamic_cast<AliEMCALDigit *>(digits->At( fDigitsList[index] )) ;
1035 AliError("No Digit!!");
1039 Int_t nparents = digit->GetNiparent() ;
1040 if ( nparents == 0 ) continue ;
1043 for ( jndex = 0 ; jndex < nparents ; jndex++ ) { // all primaries in digit
1044 if ( fMulParent > fMaxParent ) {
1046 Error("EvalParents", "increase fMaxParent") ;
1049 Int_t newParent = digit->GetIparent(jndex+1) ;
1050 Float_t newdEParent = digit->GetDEParent(jndex+1) ;
1052 Bool_t already = kFALSE ;
1053 for ( kndex = 0 ; kndex < fMulParent ; kndex++ ) { //check if not already stored
1054 if ( newParent == parentArray[kndex] ){
1055 dEParentArray[kndex] += newdEParent;
1060 if ( !already && (fMulParent < fMaxParent)) { // store it
1061 parentArray[fMulParent] = newParent ;
1062 dEParentArray[fMulParent] = newdEParent ;
1065 } // all parents in digit
1069 Int_t *sortIdx = new Int_t[fMulParent];
1070 TMath::Sort(fMulParent,dEParentArray,sortIdx);
1071 for(index = 0; index < fMulParent; index++) {
1072 fParentsList[index] = parentArray[sortIdx[index]] ;
1073 fDEParentsList[index] = dEParentArray[sortIdx[index]] ;
1078 delete [] parentArray;
1079 delete [] dEParentArray;
1082 //____________________________________________________________________________
1083 void AliEMCALRecPoint::GetLocalPosition(TVector3 & lpos) const
1085 // returns the position of the cluster in the local reference system
1086 // of the sub-detector
1091 //____________________________________________________________________________
1092 void AliEMCALRecPoint::GetGlobalPosition(TVector3 & gpos) const
1094 // returns the position of the cluster in the global reference system of ALICE
1095 // These are now the Cartesian X, Y and Z
1096 // cout<<" geom "<<geom<<endl;
1097 // fGeomPtr->GetGlobal(fLocPos, gpos, fSuperModuleNumber);
1102 //____________________________________________________________________________
1103 //void AliEMCALRecPoint::GetGlobalPosition(TVector3 & gpos, TMatrixF & gmat) const
1105 // // returns the position of the cluster in the global reference system of ALICE
1106 // // These are now the Cartesian X, Y and Z
1107 // // cout<<" geom "<<geom<<endl;
1109 // //To be implemented
1110 // fGeomPtr->GetGlobalEMCAL(this, gpos, gmat);
1114 //_____________________________________________________________________________
1115 void AliEMCALRecPoint::EvalLocal2TrackingCSTransform()
1117 //Evaluates local to "tracking" c.s. transformation (B.P.).
1118 //All evaluations should be completed before calling for this
1120 //See ALICE PPR Chapter 5 p.18 for "tracking" c.s. definition,
1121 //or just ask Jouri Belikov. :)
1123 SetVolumeId(AliGeomManager::LayerToVolUID(AliGeomManager::kEMCAL,GetSuperModuleNumber()));
1125 const TGeoHMatrix* tr2loc = GetTracking2LocalMatrix();
1126 if(!tr2loc) AliFatal(Form("No Tracking2LocalMatrix found."));
1128 Double_t lxyz[3] = {fLocPos.X(),fLocPos.Y(),fLocPos.Z()};
1129 Double_t txyz[3] = {0,0,0};
1131 tr2loc->MasterToLocal(lxyz,txyz);
1132 SetX(txyz[0]); SetY(txyz[1]); SetZ(txyz[2]);
1134 if(AliLog::GetGlobalDebugLevel()>0) {
1135 TVector3 gpos; //TMatrixF gmat;
1136 //GetGlobalPosition(gpos,gmat); //Not doing anythin special, replace by next line.
1137 fGeomPtr->GetGlobal(fLocPos, gpos, GetSuperModuleNumber());
1141 AliInfo(Form("lCS-->(%.3f,%.3f,%.3f), tCS-->(%.3f,%.3f,%.3f), gCS-->(%.3f,%.3f,%.3f), gCScalc-\
1142 ->(%.3f,%.3f,%.3f), supermodule %d",
1143 fLocPos.X(),fLocPos.Y(),fLocPos.Z(),
1144 GetX(),GetY(),GetZ(),
1145 gpos.X(),gpos.Y(),gpos.Z(),
1146 gxyz[0],gxyz[1],gxyz[2],GetSuperModuleNumber()));
1151 //____________________________________________________________________________
1152 Float_t AliEMCALRecPoint::GetMaximalEnergy(void) const
1154 // Finds the maximum energy in the cluster
1156 Float_t menergy = 0. ;
1159 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
1161 if(fEnergyList[iDigit] > menergy)
1162 menergy = fEnergyList[iDigit] ;
1167 //____________________________________________________________________________
1168 Int_t AliEMCALRecPoint::GetMaximalEnergyIndex(void) const
1170 // Finds the maximum energy in the cluster
1172 Float_t menergy = 0. ;
1176 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
1178 if(fEnergyList[iDigit] > menergy){
1179 menergy = fEnergyList[iDigit] ;
1182 }//loop on cluster digits
1188 //____________________________________________________________________________
1189 Int_t AliEMCALRecPoint::GetMultiplicityAtLevel(Float_t H) const
1191 // Calculates the multiplicity of digits with energy larger than H*energy
1195 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
1197 if(fEnergyList[iDigit] > H * fAmp)
1203 //____________________________________________________________________________
1204 Int_t AliEMCALRecPoint::GetNumberOfLocalMax(AliEMCALDigit ** maxAt, Float_t * maxAtEnergy,
1205 Float_t locMaxCut,TClonesArray * digits) const
1207 // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
1208 // energy difference between two local maxima
1210 AliEMCALDigit * digit = 0;
1211 AliEMCALDigit * digitN = 0;
1216 for(iDigit = 0; iDigit < fMulDigit; iDigit++)
1217 maxAt[iDigit] = (AliEMCALDigit*) digits->At(fDigitsList[iDigit]) ;
1219 for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {
1221 digit = maxAt[iDigit] ;
1223 for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {
1224 digitN = (AliEMCALDigit *) digits->At(fDigitsList[iDigitN]) ;
1226 if ( AreNeighbours(digit, digitN) ) {
1227 if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {
1228 maxAt[iDigitN] = 0 ;
1229 // but may be digit too is not local max ?
1230 if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut)
1235 // but may be digitN too is not local max ?
1236 if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut)
1237 maxAt[iDigitN] = 0 ;
1239 } // if Areneighbours
1245 for(iDigit = 0; iDigit < fMulDigit; iDigit++) {
1247 maxAt[iDigitN] = maxAt[iDigit] ;
1248 maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
1255 //____________________________________________________________________________
1256 Int_t AliEMCALRecPoint::GetPrimaryIndex() const
1258 // Get the primary track index in TreeK which deposits the most energy
1259 // in Digits which forms RecPoint.
1262 return fTracksList[0];
1266 //____________________________________________________________________________
1267 void AliEMCALRecPoint::EvalTime(TClonesArray * digits){
1268 // time is set to the time of the digit with the maximum energy
1272 for(Int_t idig=0; idig < fMulDigit; idig++){
1273 if(fEnergyList[idig] > maxE){
1274 maxE = fEnergyList[idig] ;
1278 fTime = ((AliEMCALDigit*) digits->At(fDigitsList[maxAt]))->GetTime() ;
1282 //______________________________________________________________________________
1283 void AliEMCALRecPoint::Paint(Option_t *)
1285 // Paint this ALiRecPoint as a TMarker with its current attributes
1287 TVector3 pos(0.,0.,0.) ;
1288 GetLocalPosition(pos) ;
1289 Coord_t x = pos.X() ;
1290 Coord_t y = pos.Z() ;
1291 Color_t markercolor = 1 ;
1292 Size_t markersize = 1.;
1293 Style_t markerstyle = 5 ;
1295 if (!gPad->IsBatch()) {
1296 gVirtualX->SetMarkerColor(markercolor) ;
1297 gVirtualX->SetMarkerSize (markersize) ;
1298 gVirtualX->SetMarkerStyle(markerstyle) ;
1300 gPad->SetAttMarkerPS(markercolor,markerstyle,markersize) ;
1301 gPad->PaintPolyMarker(1,&x,&y,"") ;
1304 //_____________________________________________________________________
1305 Double_t AliEMCALRecPoint::TmaxInCm(const Double_t e , const Int_t key)
1308 // key = 0(gamma, default)
1310 const Double_t ca = 4.82; // shower max parameter - first guess; ca=TMath::Log(1000./8.07)
1311 Double_t tmax = 0.; // position of electromagnetic shower max in cm
1313 Double_t x0 = 1.31; // radiation lenght (cm)
1314 //If old geometry in use
1315 if(!((fGeomPtr->GetEMCGeometry()->GetGeoName()).Contains("V1"))) x0 = 1.28;
1318 tmax = TMath::Log(e) + ca;
1319 if (key==0) tmax += 0.5;
1321 tmax *= x0; // convert to cm
1326 //______________________________________________________________________________
1327 Float_t AliEMCALRecPoint::EtaToTheta(Float_t arg) const
1329 //Converts Theta (Radians) to Eta(Radians)
1330 return (2.*TMath::ATan(TMath::Exp(-arg)));
1333 //______________________________________________________________________________
1334 Float_t AliEMCALRecPoint::ThetaToEta(Float_t arg) const
1336 //Converts Eta (Radians) to Theta(Radians)
1337 return (-1 * TMath::Log(TMath::Tan(0.5 * arg)));
1340 //____________________________________________________________________________
1341 void AliEMCALRecPoint::Print(Option_t *opt) const
1343 // Print the list of digits belonging to the cluster
1344 if(strlen(opt)==0) return;
1346 message = "AliEMCALRecPoint:\n" ;
1347 message += " digits # = " ;
1348 AliInfo(message.Data()) ;
1351 for(iDigit=0; iDigit<fMulDigit; iDigit++)
1352 printf(" %d ", fDigitsList[iDigit] ) ;
1355 AliInfo(" Energies = ") ;
1356 for(iDigit=0; iDigit<fMulDigit; iDigit++)
1357 printf(" %f ", fEnergyList[iDigit] ) ;
1360 AliInfo("\n Abs Ids = ") ;
1361 for(iDigit=0; iDigit<fMulDigit; iDigit++)
1362 printf(" %i ", fAbsIdList[iDigit] ) ;
1365 AliInfo(" Primaries ") ;
1366 for(iDigit = 0;iDigit < fMulTrack; iDigit++)
1367 printf(" %d ", fTracksList[iDigit]) ;
1369 printf("\n Local x %6.2f y %7.2f z %7.1f \n", fLocPos[0], fLocPos[1], fLocPos[2]);
1371 message = " ClusterType = %d" ;
1372 message += " Multiplicity = %d" ;
1373 message += " Cluster Energy = %f" ;
1374 message += " Core energy = %f" ;
1375 message += " Core radius = %f" ;
1376 message += " Number of primaries %d" ;
1377 message += " Stored at position %d" ;
1378 AliInfo(Form(message.Data(), fClusterType, fMulDigit, fAmp, fCoreEnergy, fCoreRadius, fMulTrack, GetIndexInList()) ) ;
1381 //___________________________________________________________
1382 Double_t AliEMCALRecPoint::GetPointEnergy() const
1384 //Returns energy ....
1386 for(int ic=0; ic<GetMultiplicity(); ic++) e += double(fEnergyList[ic]);