Line data Source code
1 : /**************************************************************************
2 : * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 : * *
4 : * Author: The ALICE Off-line Project. *
5 : * Contributors are mentioned in the code where appropriate. *
6 : * *
7 : * Permission to use, copy, modify and distribute this software and its *
8 : * documentation strictly for non-commercial purposes is hereby granted *
9 : * without fee, provided that the above copyright notice appears in all *
10 : * copies and that both the copyright notice and this permission notice *
11 : * appear in the supporting documentation. The authors make no claims *
12 : * about the suitability of this software for any purpose. It is *
13 : * provided "as is" without express or implied warranty. *
14 : **************************************************************************/
15 : /* $Id$ */
16 : //_________________________________________________________________________
17 : // Reconstructed Points for the EMCAL
18 : // A RecPoint is a cluster of digits
19 : //
20 : //
21 : //*-- Author: Yves Schutz (SUBATECH)
22 : //*-- Author: Dmitri Peressounko (RRC KI & SUBATECH)
23 : //*-- Author: Heather Gray (LBL) merged AliEMCALRecPoint and AliEMCALTowerRecPoint 02/04
24 :
25 : // --- ROOT system ---
26 : #include "TPad.h"
27 : #include "TGraph.h"
28 : #include "TPaveText.h"
29 : #include "TClonesArray.h"
30 : #include "TMath.h"
31 : #include "TGeoMatrix.h"
32 : #include "TGeoManager.h"
33 : #include "TGeoPhysicalNode.h"
34 : #include "TRandom.h"
35 :
36 : // --- Standard library ---
37 : #include <Riostream.h>
38 :
39 : // --- AliRoot header files ---
40 : //#include "AliGenerator.h"
41 : class AliGenerator;
42 : class AliEMCAL;
43 : #include "AliLog.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"
51 :
52 42 : ClassImp(AliEMCALRecPoint)
53 :
54 : //____________________________________________________________________________
55 : AliEMCALRecPoint::AliEMCALRecPoint()
56 132 : : AliCluster(), fGeomPtr(0),
57 66 : fAmp(0), fIndexInList(-1), //to be set when the point is already stored
58 132 : fGlobPos(0,0,0),fLocPos(0,0,0),
59 66 : fMaxDigit(100), fMulDigit(0), fMaxTrack(200),
60 66 : fMulTrack(0), fDigitsList(0), fTracksList(0),
61 66 : fClusterType(-1), fCoreEnergy(0), fDispersion(0),
62 66 : fEnergyList(0), fAbsIdList(0),
63 66 : fTime(0.), fNExMax(0), fCoreRadius(10), //HG check this
64 66 : fDETracksList(0), fMulParent(0), fMaxParent(0),
65 66 : fParentsList(0), fDEParentsList(0), fSuperModuleNumber(0),
66 66 : fDigitIndMax(-1), fDistToBadTower(-1), fSharedCluster(kFALSE)
67 330 : {
68 : // ctor
69 132 : fGeomPtr = AliEMCALGeometry::GetInstance();
70 :
71 66 : fLambda[0] = 0;
72 66 : fLambda[1] = 0;
73 :
74 132 : }
75 :
76 : //____________________________________________________________________________
77 : AliEMCALRecPoint::AliEMCALRecPoint(const char *)
78 66 : : AliCluster(), fGeomPtr(0),
79 33 : fAmp(0), fIndexInList(-1), //to be set when the point is already stored
80 66 : fGlobPos(0,0,0), fLocPos(0,0,0),
81 33 : fMaxDigit(100), fMulDigit(0), fMaxTrack(1000), fMulTrack(0),
82 99 : fDigitsList(new Int_t[fMaxDigit]), fTracksList(new Int_t[fMaxTrack]),
83 33 : fClusterType(-1), fCoreEnergy(0), fDispersion(0),
84 66 : fEnergyList(new Float_t[fMaxDigit]),
85 66 : fAbsIdList(new Int_t[fMaxDigit]), fTime(-1.), fNExMax(0), fCoreRadius(10),
86 66 : fDETracksList(new Float_t[fMaxTrack]), fMulParent(0), fMaxParent(1000),
87 99 : fParentsList(new Int_t[fMaxParent]), fDEParentsList(new Float_t[fMaxParent]),
88 33 : fSuperModuleNumber(0), fDigitIndMax(-1), fDistToBadTower(-1),fSharedCluster(kFALSE)
89 165 : {
90 : // ctor
91 66066 : for (Int_t i = 0; i < fMaxTrack; i++)
92 33000 : fDETracksList[i] = 0;
93 66066 : for (Int_t i = 0; i < fMaxParent; i++) {
94 33000 : fParentsList[i] = -1;
95 33000 : fDEParentsList[i] = 0;
96 : }
97 :
98 66 : fGeomPtr = AliEMCALGeometry::GetInstance();
99 33 : fLambda[0] = 0;
100 33 : fLambda[1] = 0;
101 66 : }
102 :
103 : //____________________________________________________________________________
104 : AliEMCALRecPoint::AliEMCALRecPoint(const AliEMCALRecPoint & rp)
105 0 : : AliCluster(rp), fGeomPtr(rp.fGeomPtr),
106 0 : fAmp(rp.fAmp), fIndexInList(rp.fIndexInList),
107 0 : fGlobPos(rp.fGlobPos),fLocPos(rp.fLocPos),
108 0 : fMaxDigit(rp.fMaxDigit), fMulDigit(rp.fMulDigit),
109 0 : fMaxTrack(rp.fMaxTrack), fMulTrack(rp.fMaxTrack),
110 0 : fDigitsList(new Int_t[rp.fMaxDigit]), fTracksList(new Int_t[rp.fMaxTrack]),
111 0 : fClusterType(rp.fClusterType), fCoreEnergy(rp.fCoreEnergy),
112 0 : fDispersion(rp.fDispersion),
113 0 : fEnergyList(new Float_t[rp.fMaxDigit]),
114 0 : fAbsIdList(new Int_t[rp.fMaxDigit]), fTime(rp.fTime), fNExMax(rp.fNExMax),fCoreRadius(rp.fCoreRadius),
115 0 : fDETracksList(new Float_t[rp.fMaxTrack]), fMulParent(rp.fMulParent),
116 0 : fMaxParent(rp.fMaxParent), fParentsList(new Int_t[rp.fMaxParent]),
117 0 : fDEParentsList(new Float_t[rp.fMaxParent]),
118 0 : fSuperModuleNumber(rp.fSuperModuleNumber), fDigitIndMax(rp.fDigitIndMax),
119 0 : fDistToBadTower(rp.fDistToBadTower), fSharedCluster(rp.fSharedCluster)
120 0 : {
121 : //copy ctor
122 0 : fLambda[0] = rp.fLambda[0];
123 0 : fLambda[1] = rp.fLambda[1];
124 :
125 0 : for(Int_t i = 0; i < rp.fMulDigit; i++) {
126 0 : fEnergyList[i] = rp.fEnergyList[i];
127 0 : fAbsIdList[i] = rp.fAbsIdList[i];
128 : }
129 :
130 0 : for(Int_t i = 0; i < rp.fMulTrack; i++) fDETracksList[i] = rp.fDETracksList[i];
131 :
132 0 : for(Int_t i = 0; i < rp.fMulParent; i++) {
133 0 : fParentsList[i] = rp.fParentsList[i];
134 0 : fDEParentsList[i] = rp.fDEParentsList[i];
135 : }
136 :
137 0 : }
138 : //____________________________________________________________________________
139 : AliEMCALRecPoint::~AliEMCALRecPoint()
140 396 : {
141 : // dtor
142 66 : if ( fEnergyList )
143 132 : delete[] fEnergyList ;
144 66 : if ( fAbsIdList )
145 132 : delete[] fAbsIdList ;
146 66 : if ( fDETracksList)
147 100 : delete[] fDETracksList;
148 66 : if ( fParentsList)
149 100 : delete[] fParentsList;
150 66 : if ( fDEParentsList)
151 100 : delete[] fDEParentsList;
152 :
153 132 : delete [] fDigitsList ;
154 116 : delete [] fTracksList ;
155 198 : }
156 :
157 : //____________________________________________________________________________
158 : AliEMCALRecPoint& AliEMCALRecPoint::operator= (const AliEMCALRecPoint &rp)
159 : {
160 : // assignment operator
161 :
162 0 : if(&rp == this) return *this;
163 :
164 0 : fGeomPtr = rp.fGeomPtr;
165 0 : fAmp = rp.fAmp;
166 0 : fIndexInList = rp.fIndexInList;
167 0 : fGlobPos = rp.fGlobPos;
168 0 : fLocPos = rp.fLocPos;
169 0 : fMaxDigit = rp.fMaxDigit;
170 0 : fMulDigit = rp.fMulDigit;
171 0 : fMaxTrack = rp.fMaxTrack;
172 0 : fMulTrack = rp.fMulTrack;
173 :
174 0 : if(fDigitsList) delete [] fDigitsList;
175 0 : fDigitsList = new Int_t[rp.fMaxDigit];
176 0 : if(fTracksList) delete [] fTracksList;
177 0 : fTracksList = new Int_t[rp.fMaxTrack];
178 0 : for(Int_t i = 0; i<fMaxDigit; i++) fDigitsList[i] = rp.fDigitsList[i];
179 0 : for(Int_t i = 0; i<fMaxTrack; i++) fTracksList[i] = rp.fTracksList[i];
180 :
181 0 : fClusterType = rp.fClusterType;
182 0 : fCoreEnergy = rp.fCoreEnergy;
183 0 : fDispersion = rp.fDispersion;
184 :
185 :
186 0 : if(fEnergyList) delete [] fEnergyList;
187 0 : fEnergyList = new Float_t[rp.fMaxDigit];
188 0 : if(fAbsIdList) delete [] fAbsIdList;
189 0 : fAbsIdList = new Int_t[rp.fMaxDigit];
190 0 : for(Int_t i = 0; i<fMaxDigit; i++) {
191 0 : fEnergyList[i] = rp.fEnergyList[i];
192 0 : fAbsIdList[i] = rp.fAbsIdList[i];
193 : }
194 :
195 0 : fTime = rp.fTime;
196 0 : fNExMax = rp.fNExMax;
197 0 : fCoreRadius = rp.fCoreRadius;
198 :
199 0 : if(fDETracksList) delete [] fDETracksList;
200 0 : fDETracksList = new Float_t[rp.fMaxTrack];
201 0 : for(Int_t i = 0; i < fMaxTrack; i++) fDETracksList[i] = rp.fDETracksList[i];
202 :
203 0 : fMulParent = rp.fMulParent;
204 0 : fMaxParent = rp.fMaxParent;
205 :
206 0 : if(fParentsList) delete [] fParentsList;
207 0 : fParentsList = new Int_t[rp.fMaxParent];
208 0 : if(fDEParentsList) delete [] fDEParentsList;
209 0 : fDEParentsList = new Float_t[rp.fMaxParent];
210 0 : for(Int_t i = 0; i < fMaxParent; i++) {
211 0 : fParentsList[i] = rp.fParentsList[i];
212 0 : fDEParentsList[i] = rp.fDEParentsList[i];
213 : }
214 :
215 0 : fSuperModuleNumber = rp.fSuperModuleNumber;
216 0 : fDigitIndMax = rp.fDigitIndMax;
217 :
218 0 : fLambda[0] = rp.fLambda[0];
219 0 : fLambda[1] = rp.fLambda[1];
220 :
221 0 : fDistToBadTower = rp.fDistToBadTower;
222 0 : fSharedCluster = rp.fSharedCluster;
223 :
224 0 : return *this;
225 :
226 0 : }
227 :
228 : //____________________________________________________________________________
229 : void AliEMCALRecPoint::AddDigit(AliEMCALDigit & digit, const Float_t energy, const Bool_t shared)
230 : {
231 : // Adds a digit to the RecPoint
232 : // and accumulates the total amplitude and the multiplicity
233 :
234 97 : if(fEnergyList == 0)
235 0 : fEnergyList = new Float_t[fMaxDigit];
236 :
237 97 : if(fAbsIdList == 0) {
238 0 : fAbsIdList = new Int_t [fMaxDigit];
239 0 : }
240 :
241 97 : if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists
242 0 : fMaxDigit*=2 ;
243 0 : Int_t * tempo = new Int_t [fMaxDigit];
244 0 : Float_t * tempoE = new Float_t[fMaxDigit];
245 0 : Int_t * tempoId = new Int_t [fMaxDigit];
246 :
247 : Int_t index ;
248 0 : for ( index = 0 ; index < fMulDigit ; index++ ){
249 0 : tempo [index] = fDigitsList[index] ;
250 0 : tempoE [index] = fEnergyList[index] ;
251 0 : tempoId[index] = fAbsIdList [index] ;
252 : }
253 :
254 0 : delete [] fDigitsList ;
255 0 : delete [] fEnergyList ;
256 0 : delete [] fAbsIdList ;
257 :
258 0 : fDigitsList = tempo;
259 0 : fEnergyList = tempoE;
260 0 : fAbsIdList = tempoId;
261 0 : } // if
262 :
263 97 : fDigitsList[fMulDigit] = digit.GetIndexInList() ;
264 97 : fEnergyList[fMulDigit] = energy ;
265 97 : fAbsIdList [fMulDigit] = digit.GetId();
266 97 : fMulDigit++ ;
267 97 : fAmp += energy ;
268 :
269 107 : if(shared) fSharedCluster = kTRUE;
270 97 : }
271 : //____________________________________________________________________________
272 : /// \return true or false if two digits are neighbours
273 : ///
274 : /// A neighbour is defined as being two digits which share a side.
275 : /// If sharing a corner, cells are not neighbours (not the case in early versions).
276 : /// Only used in case of V1 clusterizer, with and without unfolding.
277 : ///
278 : /// \param digit1 check closeness of digit1 with digit2
279 : /// \param digit2 check closeness of digit2 with digit1
280 : ///
281 : Bool_t AliEMCALRecPoint::AreNeighbours(AliEMCALDigit * digit1, AliEMCALDigit * digit2 ) const
282 : {
283 : Bool_t areNeighbours = kFALSE ;
284 588 : Int_t nSupMod = 0, nModule = 0, nIphi = 0, nIeta = 0;
285 294 : Int_t nSupMod1 = 0, nModule1 = 0, nIphi1 = 0, nIeta1 = 0;
286 294 : Int_t relid1[2] , relid2[2] ; // ieta, iphi
287 : Int_t rowdiff = 0, coldiff = 0;
288 :
289 : areNeighbours = kFALSE ;
290 :
291 294 : fGeomPtr->GetCellIndex(digit1->GetId(), nSupMod,nModule,nIphi,nIeta);
292 294 : fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nModule,nIphi,nIeta, relid1[0],relid1[1]);
293 :
294 294 : fGeomPtr->GetCellIndex(digit2->GetId(), nSupMod1,nModule1,nIphi1,nIeta1);
295 294 : fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod1,nModule1,nIphi1,nIeta1, relid2[0],relid2[1]);
296 :
297 : // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2-1
298 : // C Side impair SM, nSupMod%2=1; A side pair SM nSupMod%2=0
299 568 : if ( fSharedCluster && nSupMod1 != nSupMod )
300 : {
301 96 : if(nSupMod1%2) relid2[1]+=AliEMCALGeoParams::fgkEMCALCols;
302 44 : else relid1[1]+=AliEMCALGeoParams::fgkEMCALCols;
303 : }
304 :
305 294 : rowdiff = TMath::Abs( relid1[0] - relid2[0] ) ;
306 294 : coldiff = TMath::Abs( relid1[1] - relid2[1] ) ;
307 :
308 : // With this condition, cells touching one corner are considered neighbours
309 : // which was considered as the behavior expected initially, but changed later.
310 : //if (( coldiff <= 1 ) && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0))
311 :
312 294 : if ((coldiff + rowdiff == 1 ))
313 89 : areNeighbours = kTRUE ;
314 :
315 588 : return areNeighbours;
316 294 : }
317 :
318 : //____________________________________________________________________________
319 : Int_t AliEMCALRecPoint::Compare(const TObject * obj) const
320 : {
321 : // Compares two RecPoints according to their position in the EMCAL modules
322 :
323 : Float_t delta = 1 ; //Width of "Sorting row".
324 :
325 : Int_t rv = 2 ;
326 :
327 136 : AliEMCALRecPoint * clu = (AliEMCALRecPoint *)obj ;
328 :
329 68 : TVector3 locpos1;
330 68 : GetLocalPosition(locpos1);
331 68 : TVector3 locpos2;
332 68 : clu->GetLocalPosition(locpos2);
333 :
334 68 : Int_t rowdif = (Int_t)(TMath::Ceil(locpos1.X()/delta)-TMath::Ceil(locpos2.X()/delta)) ;
335 68 : if (rowdif> 0)
336 35 : rv = 1 ;
337 33 : else if(rowdif < 0)
338 20 : rv = -1 ;
339 13 : else if(locpos1.Y()>locpos2.Y())
340 7 : rv = -1 ;
341 : else
342 : rv = 1 ;
343 :
344 : return rv ;
345 68 : }
346 :
347 : //___________________________________________________________________________
348 : void AliEMCALRecPoint::Draw(Option_t *option)
349 : {
350 : // Draw this AliEMCALRecPoint with its current attributes
351 :
352 0 : AppendPad(option);
353 0 : }
354 :
355 : //____________________________________________________________________________
356 : void AliEMCALRecPoint::EvalAll(Float_t logWeight,TClonesArray * digits, const Bool_t justClusters)
357 : {
358 : // Evaluates cluster parameters
359 :
360 : // First calculate the index of digit with maximum amplitude and get
361 : // the supermodule number where it sits.
362 :
363 33 : fDigitIndMax = GetMaximalEnergyIndex();
364 33 : fSuperModuleNumber = fGeomPtr->GetSuperModuleNumber(GetAbsIdMaxDigit());
365 :
366 : //Evaluate global and local position
367 33 : EvalGlobalPosition(logWeight, digits) ;
368 33 : EvalLocalPosition(logWeight, digits) ;
369 :
370 : //Evaluate shower parameters
371 33 : EvalElipsAxis(logWeight, digits) ;
372 33 : EvalDispersion(logWeight, digits) ;
373 :
374 : //EvalCoreEnergy(logWeight, digits);
375 33 : EvalTime(digits) ;
376 33 : EvalPrimaries(digits) ;
377 33 : EvalParents(digits);
378 :
379 : //Called last because it sets the global position of the cluster?
380 : //Do not call it when recalculating clusters out of standard reconstruction
381 33 : if(!justClusters){
382 33 : EvalLocal2TrackingCSTransform();
383 33 : }
384 :
385 33 : }
386 :
387 : //____________________________________________________________________________
388 : void AliEMCALRecPoint::EvalDispersion(Float_t logWeight, TClonesArray * digits)
389 : {
390 : // Calculates the dispersion of the shower at the origin of the RecPoint
391 : // in cell units - Nov 16,2006
392 :
393 : Double_t d = 0., wtot = 0., w = 0.;
394 : Int_t iDigit=0, nstat=0;
395 : AliEMCALDigit * digit=0;
396 :
397 : // Calculates the dispersion in cell units
398 : Double_t etai, phii, etaMean=0.0, phiMean=0.0;
399 66 : int nSupMod=0, nModule=0, nIphi=0, nIeta=0;
400 33 : int iphi=0, ieta=0;
401 : // Calculate mean values
402 260 : for(iDigit=0; iDigit < fMulDigit; iDigit++) {
403 97 : digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit]) ;
404 :
405 194 : if (fAmp>0 && fEnergyList[iDigit]>0) {
406 97 : fGeomPtr->GetCellIndex(digit->GetId(), nSupMod,nModule,nIphi,nIeta);
407 97 : fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nModule,nIphi,nIeta, iphi,ieta);
408 :
409 : // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2
410 : // C Side impair SM, nSupMod%2=1; A side pair SM nSupMod%2=0
411 207 : if(fSharedCluster && nSupMod%2) ieta+=AliEMCALGeoParams::fgkEMCALCols;
412 :
413 97 : etai=(Double_t)ieta;
414 97 : phii=(Double_t)iphi;
415 97 : w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
416 :
417 97 : if(w>0.0) {
418 88 : phiMean += phii*w;
419 88 : etaMean += etai*w;
420 88 : wtot += w;
421 88 : }
422 : }
423 : }
424 33 : if (wtot>0) {
425 33 : phiMean /= wtot ;
426 33 : etaMean /= wtot ;
427 33 : } else AliError(Form("Wrong weight %f\n", wtot));
428 :
429 : // Calculate dispersion
430 260 : for(iDigit=0; iDigit < fMulDigit; iDigit++) {
431 97 : digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit]) ;
432 :
433 194 : if (fAmp>0 && fEnergyList[iDigit]>0) {
434 97 : fGeomPtr->GetCellIndex(digit->GetId(), nSupMod,nModule,nIphi,nIeta);
435 97 : fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nModule,nIphi,nIeta, iphi,ieta);
436 :
437 : // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2
438 : // C Side impair SM, nSupMod%2=1; A side pair SM, nSupMod%2=0
439 207 : if(fSharedCluster && nSupMod%2) ieta+=AliEMCALGeoParams::fgkEMCALCols;
440 :
441 97 : etai=(Double_t)ieta;
442 97 : phii=(Double_t)iphi;
443 97 : w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
444 :
445 97 : if(w>0.0) {
446 88 : nstat++;
447 88 : d += w*((etai-etaMean)*(etai-etaMean)+(phii-phiMean)*(phii-phiMean));
448 88 : }
449 : }
450 : }
451 :
452 48 : if ( wtot > 0 && nstat>1) d /= wtot ;
453 : else d = 0. ;
454 :
455 33 : fDispersion = TMath::Sqrt(d) ;
456 : //printf("AliEMCALRecPoint::EvalDispersion() : Dispersion %f \n",fDispersion);
457 33 : }
458 :
459 : //____________________________________________________________________________
460 : void AliEMCALRecPoint::EvalDistanceToBadChannels(AliCaloCalibPedestal* caloped)
461 : {
462 : //For each EMC rec. point set the distance to the nearest bad channel.
463 : //AliInfo(Form("%d bad channel(s) found.\n", caloped->GetDeadTowerCount()));
464 : //It is done in cell units and not in global or local position as before (Sept 2010)
465 :
466 66 : if(!caloped->GetDeadTowerCount()) return;
467 :
468 : //Get channels map of the supermodule where the cluster is.
469 0 : TH2D* hMap = caloped->GetDeadMap(fSuperModuleNumber);
470 :
471 : Int_t dRrow, dReta;
472 : Float_t minDist = 10000.;
473 : Float_t dist = 0.;
474 0 : Int_t nSupMod, nModule;
475 0 : Int_t nIphi, nIeta;
476 0 : Int_t iphi, ieta;
477 0 : fDigitIndMax = GetMaximalEnergyIndex();
478 0 : fGeomPtr->GetCellIndex(fAbsIdList[fDigitIndMax], nSupMod,nModule,nIphi,nIeta);
479 0 : fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nModule,nIphi,nIeta, iphi,ieta);
480 :
481 : //Loop on tower status map
482 0 : for(Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++){
483 0 : for(Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++){
484 : //Check if tower is bad.
485 0 : if(hMap->GetBinContent(icol,irow)==AliCaloCalibPedestal::kAlive) continue;
486 : //printf("AliEMCALRecPoint::EvalDistanceToBadChannels() - Bad channel in SM %d, col %d, row %d\n",iSM,icol, irow);
487 :
488 0 : dRrow=TMath::Abs(irow-iphi);
489 0 : dReta=TMath::Abs(icol-ieta);
490 0 : dist=TMath::Sqrt(dRrow*dRrow+dReta*dReta);
491 0 : if(dist < minDist) minDist = dist;
492 :
493 : }
494 : }
495 :
496 : //In case the cluster is shared by 2 SuperModules, need to check the map of the second Super Module
497 0 : if (fSharedCluster) {
498 : TH2D* hMap2 = 0;
499 : Int_t nSupMod2 = -1;
500 :
501 : //The only possible combinations are (0,1), (2,3) ... (10,11)
502 0 : if(fSuperModuleNumber%2) nSupMod2 = fSuperModuleNumber-1;
503 0 : else nSupMod2 = fSuperModuleNumber+1;
504 0 : hMap2 = caloped->GetDeadMap(nSupMod2);
505 :
506 : //Loop on tower status map of second super module
507 0 : for(Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++){
508 0 : for(Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++){
509 : //Check if tower is bad.
510 0 : if(hMap2->GetBinContent(icol,irow)==AliCaloCalibPedestal::kAlive) continue;
511 : //printf("AliEMCALRecPoint::EvalDistanceToBadChannels() - Bad channel in SM %d, col %d, row %d\n",iSM,icol, irow);
512 0 : dRrow=TMath::Abs(irow-iphi);
513 :
514 0 : if(fSuperModuleNumber%2) {
515 0 : dReta=TMath::Abs(icol-(AliEMCALGeoParams::fgkEMCALCols+ieta));
516 0 : }
517 : else {
518 0 : dReta=TMath::Abs(AliEMCALGeoParams::fgkEMCALCols+icol-ieta);
519 : }
520 :
521 0 : dist=TMath::Sqrt(dRrow*dRrow+dReta*dReta);
522 0 : if(dist < minDist) minDist = dist;
523 :
524 : }
525 : }
526 :
527 0 : }// shared cluster in 2 SuperModules
528 :
529 0 : fDistToBadTower = minDist;
530 : //printf("AliEMCALRecPoint::EvalDistanceToBadChannel() - Distance to Bad is %f cm, shared cluster? %d \n",fDistToBadTower,fSharedCluster);
531 33 : }
532 :
533 :
534 : //____________________________________________________________________________
535 : void AliEMCALRecPoint::EvalLocalPosition(Float_t logWeight, TClonesArray * digits)
536 : {
537 : // Calculates the center of gravity in the local EMCAL-module coordinates
538 : // Info("Print", " logWeight %f : cluster energy %f ", logWeight, fAmp); // for testing
539 :
540 : AliEMCALDigit * digit=0;
541 : Int_t i=0, nstat=0;
542 :
543 66 : Double_t dist = TmaxInCm(Double_t(fAmp));
544 : //Int_t idMax = GetAbsIdMaxDigit();// idMax is not used at all in RelPosCellInSModule, why use it?
545 :
546 33 : Double_t clXYZ[3]={0.,0.,0.}, clRmsXYZ[3]={0.,0.,0.}, xyzi[3], wtot=0., w=0.;
547 :
548 : //printf(" dist : %f e : %f \n", dist, fAmp);
549 260 : for(Int_t iDigit=0; iDigit<fMulDigit; iDigit++) {
550 291 : digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
551 :
552 97 : if(!digit) {
553 0 : AliError("No Digit!!");
554 0 : continue;
555 : }
556 :
557 97 : fGeomPtr->RelPosCellInSModule(digit->GetId(), dist, xyzi[0], xyzi[1], xyzi[2]);
558 :
559 : //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
560 170 : if(fSharedCluster && fSuperModuleNumber != fGeomPtr->GetSuperModuleNumber(digit->GetId())) xyzi[1]*=-1;
561 :
562 : //printf("EvalLocalPosition Cell: Id %i, SM %i : dist %f Local x,y,z %f %f %f \n",
563 : // digit->GetId(), fGeomPtr->GetSuperModuleNumber(digit->GetId()), dist, xyzi[0], xyzi[1], xyzi[2]);
564 :
565 194 : if(logWeight > 0.0) w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ));
566 0 : else w = fEnergyList[iDigit]; // just energy
567 :
568 97 : if(w>0.0) {
569 88 : wtot += w ;
570 88 : nstat++;
571 704 : for(i=0; i<3; i++ ) {
572 264 : clXYZ[i] += (w*xyzi[i]);
573 264 : clRmsXYZ[i] += (w*xyzi[i]*xyzi[i]);
574 : }
575 : }
576 : }
577 : // cout << " wtot " << wtot << endl;
578 33 : if ( wtot > 0 ) {
579 : // xRMS = TMath::Sqrt(x2m - xMean*xMean);
580 231 : for(i=0; i<3; i++ ) {
581 99 : clXYZ[i] /= wtot;
582 99 : if(nstat>1) {
583 45 : clRmsXYZ[i] /= (wtot*wtot);
584 45 : clRmsXYZ[i] = clRmsXYZ[i] - clXYZ[i]*clXYZ[i];
585 45 : if(clRmsXYZ[i] > 0.0) {
586 0 : clRmsXYZ[i] = TMath::Sqrt(clRmsXYZ[i]);
587 45 : } else clRmsXYZ[i] = 0;
588 54 : } else clRmsXYZ[i] = 0;
589 : }
590 : } else {
591 0 : for(i=0; i<3; i++ ) {
592 0 : clXYZ[i] = clRmsXYZ[i] = -1.;
593 : }
594 : }
595 :
596 : // // Cluster of one single digit, smear the position to avoid discrete position
597 : // // smear x and z with +- 3 cm to uniform (avoid discrete effects). Tower size is approx 6 cm.
598 : // // Rndm generates a number in ]0,1]
599 : // if (fMulDigit==1) {
600 : // clXYZ[0] += fGeomPtr->GetPhiTileSize()*(0.5 - gRandom->Rndm());
601 : // clXYZ[2] += fGeomPtr->GetEtaTileSize()*(0.5 - gRandom->Rndm());
602 : // }
603 :
604 : //Set position in local vector
605 33 : fLocPos.SetX(clXYZ[0]);
606 33 : fLocPos.SetY(clXYZ[1]);
607 33 : fLocPos.SetZ(clXYZ[2]);
608 :
609 33 : if (gDebug==2)
610 0 : printf("EvalLocalPosition Cluster: Local (x,y,z) = (%f,%f,%f) \n", fLocPos.X(), fLocPos.Y(), fLocPos.Z()) ;
611 :
612 33 : }
613 :
614 :
615 : //____________________________________________________________________________
616 : void AliEMCALRecPoint::EvalGlobalPosition(Float_t logWeight, TClonesArray * digits)
617 : {
618 : // Calculates the center of gravity in the global ALICE coordinates
619 : // Info("Print", " logWeight %f : cluster energy %f ", logWeight, fAmp); // for testing
620 :
621 : AliEMCALDigit * digit=0;
622 : Int_t i=0, nstat=0;
623 :
624 66 : Double_t dist = TmaxInCm(Double_t(fAmp));
625 : //Int_t idMax = GetAbsIdMaxDigit();// idMax is not used at all in RelPosCellInSModule, why use it?
626 :
627 33 : Double_t clXYZ[3]={0.,0.,0.}, clRmsXYZ[3]={0.,0.,0.}, lxyzi[3], xyzi[3], wtot=0., w=0.;
628 :
629 : //printf(" dist : %f e : %f \n", dist, fAmp);
630 260 : for(Int_t iDigit=0; iDigit<fMulDigit; iDigit++) {
631 291 : digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
632 :
633 97 : if(!digit) {
634 0 : AliError("No Digit!!");
635 0 : continue;
636 : }
637 :
638 : //Get the local coordinates of the cell
639 97 : fGeomPtr->RelPosCellInSModule(digit->GetId(), dist, lxyzi[0], lxyzi[1], lxyzi[2]);
640 :
641 : //Now get the global coordinate
642 97 : fGeomPtr->GetGlobal(lxyzi,xyzi, fGeomPtr->GetSuperModuleNumber(digit->GetId()));
643 : //TVector3 pos(xyzi[0], xyzi[1], xyzi[2]);
644 : //printf("EvalGlobalPosition Cell: Id %i, SM %i : dist %f Local (x,y,z) = (%f %f %f), eta %f, phi%f \n",
645 : // digit->GetId(), fGeomPtr->GetSuperModuleNumber(digit->GetId()),dist, xyzi[0], xyzi[1], xyzi[2],pos.Eta(),pos.Phi()*TMath::RadToDeg());
646 :
647 194 : if(logWeight > 0.0) w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ));
648 0 : else w = fEnergyList[iDigit]; // just energy
649 :
650 97 : if(w>0.0) {
651 88 : wtot += w ;
652 88 : nstat++;
653 704 : for(i=0; i<3; i++ ) {
654 264 : clXYZ[i] += (w*xyzi[i]);
655 264 : clRmsXYZ[i] += (w*xyzi[i]*xyzi[i]);
656 : }
657 : }
658 : }
659 : // cout << " wtot " << wtot << endl;
660 33 : if ( wtot > 0 ) {
661 : // xRMS = TMath::Sqrt(x2m - xMean*xMean);
662 231 : for(i=0; i<3; i++ ) {
663 99 : clXYZ[i] /= wtot;
664 99 : if(nstat>1) {
665 45 : clRmsXYZ[i] /= (wtot*wtot);
666 45 : clRmsXYZ[i] = clRmsXYZ[i] - clXYZ[i]*clXYZ[i];
667 45 : if(clRmsXYZ[i] > 0.0) {
668 2 : clRmsXYZ[i] = TMath::Sqrt(clRmsXYZ[i]);
669 45 : } else clRmsXYZ[i] = 0;
670 54 : } else clRmsXYZ[i] = 0;
671 : }
672 : } else {
673 0 : for(i=0; i<3; i++ ) {
674 0 : clXYZ[i] = clRmsXYZ[i] = -1.;
675 : }
676 : }
677 :
678 : // // Cluster of one single digit, smear the position to avoid discrete position
679 : // // smear x and z with +- 3 cm to uniform (avoid discrete effects). Tower size is approx 6 cm.
680 : // // Rndm generates a number in ]0,1]
681 : // if (fMulDigit==1) {
682 : // clXYZ[0] += fGeomPtr->GetPhiTileSize()*(0.5 - gRandom->Rndm());
683 : // clXYZ[2] += fGeomPtr->GetEtaTileSize()*(0.5 - gRandom->Rndm());
684 : // }
685 :
686 : //Set position in global vector
687 33 : fGlobPos.SetX(clXYZ[0]);
688 33 : fGlobPos.SetY(clXYZ[1]);
689 33 : fGlobPos.SetZ(clXYZ[2]);
690 :
691 33 : if (gDebug==2)
692 0 : printf("EvalGlobalPosition Cluster: (x ,y ,z) = (%f,%f,%f), eta %f,phi %f\n",
693 0 : fGlobPos.X(), fGlobPos.Y(), fGlobPos.Z(),fGlobPos.Eta(),fGlobPos.Phi()*TMath::RadToDeg()) ;
694 33 : }
695 :
696 : //____________________________________________________________________________
697 : void AliEMCALRecPoint::EvalLocalPositionFit(Double_t deff, Double_t logWeight,
698 : Double_t phiSlope, TClonesArray * digits)
699 : {
700 : // Evaluates local position of clusters in SM
701 :
702 : Double_t ycorr=0;
703 : AliEMCALDigit *digit=0;
704 : Int_t i=0, nstat=0;
705 0 : Double_t clXYZ[3]={0.,0.,0.}, clRmsXYZ[3]={0.,0.,0.}, xyzi[3], wtot=0., w=0.;
706 :
707 0 : Double_t dist = TmaxInCm(Double_t(fAmp));
708 : //Int_t idMax = GetAbsIdMaxDigit();// idMax is not used at all in RelPosCellInSModule, why use it?
709 :
710 0 : for(Int_t iDigit=0; iDigit<digits->GetEntries(); iDigit++) {
711 0 : digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
712 0 : if(digit){
713 : dist = deff;
714 : //fGeomPtr->RelPosCellInSModule(digit->GetId(), idMax, dist, xyzi[0], xyzi[1], xyzi[2]);
715 0 : fGeomPtr->RelPosCellInSModule(digit->GetId(), dist, xyzi[0], xyzi[1], xyzi[2]);
716 :
717 0 : if(logWeight > 0.0) w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ));
718 0 : else w = fEnergyList[iDigit]; // just energy
719 :
720 0 : if(w>0.0) {
721 0 : wtot += w ;
722 0 : nstat++;
723 0 : for(i=0; i<3; i++ ) {
724 0 : clXYZ[i] += (w*xyzi[i]);
725 0 : clRmsXYZ[i] += (w*xyzi[i]*xyzi[i]);
726 : }
727 : }
728 0 : }else AliError("Digit null");
729 : }//loop
730 : // cout << " wtot " << wtot << endl;
731 0 : if ( wtot > 0 ) {
732 : // xRMS = TMath::Sqrt(x2m - xMean*xMean);
733 0 : for(i=0; i<3; i++ ) {
734 0 : clXYZ[i] /= wtot;
735 0 : if(nstat>1) {
736 0 : clRmsXYZ[i] /= (wtot*wtot);
737 0 : clRmsXYZ[i] = clRmsXYZ[i] - clXYZ[i]*clXYZ[i];
738 0 : if(clRmsXYZ[i] > 0.0) {
739 0 : clRmsXYZ[i] = TMath::Sqrt(clRmsXYZ[i]);
740 0 : } else clRmsXYZ[i] = 0;
741 0 : } else clRmsXYZ[i] = 0;
742 : }
743 : } else {
744 0 : for(i=0; i<3; i++ ) {
745 0 : clXYZ[i] = clRmsXYZ[i] = -1.;
746 : }
747 : }
748 : // clRmsXYZ[i] ??
749 0 : if(phiSlope != 0.0 && logWeight > 0.0 && wtot) {
750 : // Correction in phi direction (y - coords here); Aug 16;
751 : // May be put to global level or seperate method
752 0 : ycorr = clXYZ[1] * (1. + phiSlope);
753 : //printf(" y %f : ycorr %f : slope %f \n", clXYZ[1], ycorr, phiSlope);
754 0 : clXYZ[1] = ycorr;
755 0 : }
756 :
757 0 : fLocPos.SetX(clXYZ[0]);
758 0 : fLocPos.SetY(clXYZ[1]);
759 0 : fLocPos.SetZ(clXYZ[2]);
760 :
761 : // if (gDebug==2)
762 : // printf("EvalLocalPosition: eta,phi,r = %f,%f,%f", fLocPos.X(), fLocPos.Y(), fLocPos.Z()) ;
763 0 : }
764 :
765 : //_____________________________________________________________________________
766 : Bool_t AliEMCALRecPoint::EvalLocalPosition2(TClonesArray * digits, TArrayD &ed)
767 : {
768 : // Evaluated local position of rec.point using digits
769 : // and parametrisation of w0 and deff
770 : //printf(" <I> AliEMCALRecPoint::EvalLocalPosition2() \n");
771 0 : return AliEMCALRecPoint::EvalLocalPositionFromDigits(digits, ed, fLocPos);
772 : }
773 :
774 : //_____________________________________________________________________________
775 : Bool_t AliEMCALRecPoint::EvalLocalPositionFromDigits(TClonesArray *digits, TArrayD &ed, TVector3 &locPos)
776 : {
777 : // Used when digits should be recalibrated
778 0 : Double_t deff=0, w0=0, esum=0;
779 : Int_t iDigit=0;
780 : // AliEMCALDigit *digit;
781 :
782 0 : if(ed.GetSize() && (digits->GetEntries()!=ed.GetSize())) return kFALSE;
783 :
784 : // Calculate sum energy of digits
785 : esum = 0.0;
786 0 : for(iDigit=0; iDigit<ed.GetSize(); iDigit++) esum += ed[iDigit];
787 :
788 0 : GetDeffW0(esum, deff, w0);
789 :
790 0 : return EvalLocalPositionFromDigits(esum, deff, w0, digits, ed, locPos);
791 0 : }
792 :
793 : //_____________________________________________________________________________
794 : Bool_t AliEMCALRecPoint::EvalLocalPositionFromDigits(const Double_t esum, const Double_t deff, const Double_t w0, TClonesArray *digits, TArrayD &ed, TVector3 &locPos)
795 : {
796 : //Evaluate position of digits in supermodule.
797 : AliEMCALDigit *digit=0;
798 :
799 : Int_t i=0, nstat=0;
800 0 : Double_t clXYZ[3]={0.,0.,0.}, xyzi[3], wtot=0., w=0.;
801 : //Int_t idMax = GetAbsIdMaxDigit();// idMax is not used at all in RelPosCellInSModule, why use it?
802 :
803 : // Get pointer to EMCAL geometry
804 : // (can't use fGeomPtr in static method)
805 0 : AliEMCALGeometry* geo = AliEMCALGeometry::GetInstance();
806 :
807 0 : for(Int_t iDigit=0; iDigit<digits->GetEntries(); iDigit++) {
808 0 : digit = dynamic_cast<AliEMCALDigit *>(digits->At(iDigit));
809 0 : if(digit){
810 : //geo->RelPosCellInSModule(digit->GetId(), idMax, deff, xyzi[0], xyzi[1], xyzi[2]);
811 0 : geo->RelPosCellInSModule(digit->GetId(), deff, xyzi[0], xyzi[1], xyzi[2]);
812 :
813 0 : if(w0 > 0.0) w = TMath::Max( 0., w0 + TMath::Log(ed[iDigit] / esum));
814 0 : else w = ed[iDigit]; // just energy
815 :
816 0 : if(w>0.0) {
817 0 : wtot += w ;
818 0 : nstat++;
819 0 : for(i=0; i<3; i++ ) {
820 0 : clXYZ[i] += (w*xyzi[i]);
821 : }
822 : }
823 0 : }else AliError("Digit null");
824 : }//loop
825 : // cout << " wtot " << wtot << endl;
826 0 : if (wtot > 0) {
827 0 : for(i=0; i<3; i++ ) {
828 0 : clXYZ[i] /= wtot;
829 : }
830 0 : locPos.SetX(clXYZ[0]);
831 0 : locPos.SetY(clXYZ[1]);
832 0 : locPos.SetZ(clXYZ[2]);
833 0 : return kTRUE;
834 : } else {
835 0 : return kFALSE;
836 : }
837 :
838 0 : }
839 :
840 : //_____________________________________________________________________________
841 : void AliEMCALRecPoint::GetDeffW0(const Double_t esum , Double_t &deff, Double_t &w0)
842 : {
843 : //
844 : // Aug 31, 2001
845 : // Applied for simulation data with threshold 3 adc
846 : // Calculate efective distance (deff) and weigh parameter (w0)
847 : // for coordinate calculation; 0.5 GeV < esum <100 GeV.
848 : // Look to: http://rhic.physics.wayne.edu/~pavlinov/ALICE/SHISHKEBAB/RES/CALIB/GEOMCORR/deffandW0VaEgamma_2.gif
849 : //
850 : Double_t e=0.0;
851 : const Double_t kdp0=9.25147, kdp1=1.16700; // Hard coded now
852 : const Double_t kwp0=4.83713, kwp1=-2.77970e-01, kwp2 = 4.41116;
853 :
854 : // No extrapolation here
855 0 : e = esum<0.5?0.5:esum;
856 0 : e = e>100.?100.:e;
857 :
858 0 : deff = kdp0 + kdp1*TMath::Log(e);
859 0 : w0 = kwp0 / (1. + TMath::Exp(kwp1*(e+kwp2)));
860 : //printf("<I> AliEMCALRecPoint::GetDeffW0 esum %5.2f : deff %5.2f : w0 %5.2f \n", esum, deff, w0);
861 0 : }
862 :
863 : //______________________________________________________________________________
864 : void AliEMCALRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits)
865 : {
866 : // This function calculates energy in the core,
867 : // i.e. within a radius rad = fCoreEnergy around the center. Beyond this radius
868 : // in accordance with shower profile the energy deposition
869 : // should be less than 2%
870 : // Unfinished - Nov 15,2006
871 : // Distance is calculate in (phi,eta) units
872 :
873 : AliEMCALDigit * digit = 0 ;
874 :
875 : Int_t iDigit=0;
876 :
877 0 : if (!fLocPos.Mag()) {
878 0 : EvalLocalPosition(logWeight, digits);
879 0 : }
880 :
881 0 : Double_t phiPoint = fLocPos.Phi(), etaPoint = fLocPos.Eta();
882 0 : Double_t eta, phi, distance;
883 0 : for(iDigit=0; iDigit < fMulDigit; iDigit++) {
884 0 : digit = (AliEMCALDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
885 :
886 0 : eta = phi = 0.0;
887 0 : fGeomPtr->EtaPhiFromIndex(digit->GetId(),eta, phi) ;
888 0 : phi = phi * TMath::DegToRad();
889 :
890 0 : distance = TMath::Sqrt((eta-etaPoint)*(eta-etaPoint)+(phi-phiPoint)*(phi-phiPoint));
891 0 : if(distance < fCoreRadius)
892 0 : fCoreEnergy += fEnergyList[iDigit] ;
893 : }
894 :
895 0 : }
896 : //____________________________________________________________________________
897 : void AliEMCALRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
898 : {
899 : // Calculates the axis of the shower ellipsoid in eta and phi
900 : // in cell units
901 :
902 66 : TString gn(fGeomPtr->GetName());
903 :
904 : Double_t wtot = 0.;
905 : Double_t x = 0.;
906 : Double_t z = 0.;
907 : Double_t dxx = 0.;
908 : Double_t dzz = 0.;
909 : Double_t dxz = 0.;
910 :
911 : AliEMCALDigit * digit = 0;
912 :
913 : Double_t etai =0, phii=0, w=0;
914 33 : int nSupMod=0, nModule=0, nIphi=0, nIeta=0;
915 33 : int iphi=0, ieta=0;
916 260 : for(Int_t iDigit=0; iDigit<fMulDigit; iDigit++) {
917 194 : digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit]) ;
918 : etai = phii = 0.;
919 : // Nov 15,2006 - use cell numbers as coordinates
920 : // Copied for shish-kebab geometry, ieta,iphi is cast as double as eta,phi
921 : // We can use the eta,phi(or coordinates) of cell
922 97 : nSupMod = nModule = nIphi = nIeta = iphi = ieta = 0;
923 :
924 97 : fGeomPtr->GetCellIndex(digit->GetId(), nSupMod,nModule,nIphi,nIeta);
925 97 : fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nModule,nIphi,nIeta, iphi,ieta);
926 :
927 : // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2
928 : // C Side impair SM, nSupMod%2=1; A side pair SM, nSupMod%2=0
929 207 : if(fSharedCluster && nSupMod%2) ieta+=AliEMCALGeoParams::fgkEMCALCols;
930 :
931 97 : etai=(Double_t)ieta;
932 97 : phii=(Double_t)iphi;
933 :
934 97 : w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
935 : // fAmp summed amplitude of digits, i.e. energy of recpoint
936 : // Gives smaller value of lambda than log weight
937 : // w = fEnergyList[iDigit] / fAmp; // Nov 16, 2006 - try just energy
938 :
939 97 : dxx += w * etai * etai ;
940 97 : x += w * etai ;
941 97 : dzz += w * phii * phii ;
942 97 : z += w * phii ;
943 :
944 97 : dxz += w * etai * phii ;
945 :
946 97 : wtot += w ;
947 : }
948 :
949 33 : if ( wtot > 0 ) {
950 33 : dxx /= wtot ;
951 33 : x /= wtot ;
952 33 : dxx -= x * x ;
953 33 : dzz /= wtot ;
954 33 : z /= wtot ;
955 33 : dzz -= z * z ;
956 33 : dxz /= wtot ;
957 33 : dxz -= x * z ;
958 :
959 33 : fLambda[0] = 0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
960 33 : if(fLambda[0] > 0)
961 15 : fLambda[0] = TMath::Sqrt(fLambda[0]) ;
962 : else
963 18 : fLambda[0] = 0;
964 :
965 33 : fLambda[1] = 0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
966 :
967 33 : if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
968 10 : fLambda[1] = TMath::Sqrt(fLambda[1]) ;
969 : else
970 23 : fLambda[1]= 0. ;
971 : } else {
972 0 : fLambda[0]= 0. ;
973 0 : fLambda[1]= 0. ;
974 : }
975 :
976 : //printf("AliEMCALRecPoint::EvalElipsAxis() lambdas = %f,%f \n", fLambda[0],fLambda[1]) ;
977 :
978 33 : }
979 :
980 : ///
981 : /// Constructs the list of primary particles which
982 : /// have contributed to this RecPoint and calculate deposited energy.
983 : /// List of primaries ordered from highest to lowest energy deposition.
984 : ///
985 : //______________________________________________________________________________
986 : void AliEMCALRecPoint::EvalPrimaries(TClonesArray * digits)
987 : {
988 : AliEMCALDigit * digit =0;
989 66 : Int_t * primArray = new Int_t[fMaxTrack] ;
990 33 : memset(primArray,-1,sizeof(Int_t)*fMaxTrack);
991 33 : Float_t * dEPrimArray = new Float_t[fMaxTrack] ;
992 33 : memset(dEPrimArray,-1,sizeof(Int_t)*fMaxTrack);
993 :
994 : // All digits in rec point
995 : Int_t index ;
996 260 : for ( index = 0 ; index < GetDigitsMultiplicity() ; index++ )
997 : {
998 291 : digit = dynamic_cast<AliEMCALDigit *>(digits->At( fDigitsList[index] )) ;
999 :
1000 97 : if(!digit)
1001 : {
1002 0 : AliError("No Digit!!");
1003 0 : continue;
1004 : }
1005 :
1006 97 : Int_t nprimaries = digit->GetNprimary() ;
1007 150 : if ( nprimaries == 0 ) continue ;
1008 :
1009 : // All primaries in digit
1010 : Int_t jndex ;
1011 176 : for ( jndex = 0 ; jndex < nprimaries ; jndex++ )
1012 : {
1013 44 : if ( fMulTrack > fMaxTrack )
1014 : {
1015 0 : fMulTrack = fMaxTrack ;
1016 0 : Error("EvalPrimaries", "increase fMaxTrack ") ;
1017 0 : break ;
1018 : }
1019 :
1020 44 : Int_t newPrimary = digit->GetPrimary (jndex+1);
1021 44 : Float_t dEPrimary = digit->GetDEPrimary(jndex+1);
1022 :
1023 : // Check if not already stored
1024 : Int_t kndex ;
1025 : Bool_t already = kFALSE ;
1026 88 : for ( kndex = 0 ; kndex < fMulTrack ; kndex++ )
1027 : {
1028 27 : if ( newPrimary == primArray[kndex] ){
1029 : already = kTRUE ;
1030 27 : dEPrimArray[kndex] += dEPrimary;
1031 27 : break ;
1032 : }
1033 : } // end of check
1034 :
1035 : // Store it
1036 61 : if ( !already && (fMulTrack < fMaxTrack))
1037 : {
1038 17 : primArray [fMulTrack] = newPrimary ;
1039 17 : dEPrimArray[fMulTrack] = dEPrimary ;
1040 17 : fMulTrack++ ;
1041 17 : } // store it
1042 : } // all primaries in digit
1043 44 : } // all digits
1044 :
1045 33 : Int_t *sortIdx = new Int_t[fMulTrack];
1046 33 : TMath::Sort(fMulTrack,dEPrimArray,sortIdx);
1047 :
1048 100 : for(index = 0; index < fMulTrack; index++)
1049 : {
1050 17 : fTracksList [index] = primArray [sortIdx[index]] ;
1051 17 : fDETracksList[index] = dEPrimArray[sortIdx[index]] ;
1052 : }
1053 :
1054 66 : delete [] sortIdx;
1055 66 : delete [] primArray ;
1056 66 : delete [] dEPrimArray ;
1057 33 : }
1058 :
1059 : ///
1060 : /// Constructs the list of parent particles which have contributed to this RecPoint.
1061 : /// Access the digits belonging to this recPoint and get the labels and energy deposition.
1062 : /// List of parents ordered from highest to lowest energy deposition.
1063 : ///
1064 : //______________________________________________________________________________
1065 : void AliEMCALRecPoint::EvalParents(TClonesArray * digits)
1066 : {
1067 : AliEMCALDigit * digit=0 ;
1068 66 : Int_t * parentArray = new Int_t[fMaxTrack] ;
1069 33 : memset(parentArray,-1,sizeof(Int_t)*fMaxTrack);
1070 33 : Float_t * dEParentArray = new Float_t[fMaxTrack] ;
1071 33 : memset(dEParentArray,-1,sizeof(Int_t)*fMaxTrack);
1072 :
1073 : // Loop on all digits in rec point (cluster)
1074 : Int_t index ;
1075 260 : for ( index = 0 ; index < GetDigitsMultiplicity() ; index++ )
1076 : {
1077 194 : if (fDigitsList[index] >= digits->GetEntries() || fDigitsList[index] < 0)
1078 0 : AliError(Form("Trying to get invalid digit %d (idx in WriteRecPoint %d)",fDigitsList[index],index));
1079 :
1080 291 : digit = dynamic_cast<AliEMCALDigit *>(digits->At( fDigitsList[index] )) ;
1081 :
1082 97 : if(!digit)
1083 : {
1084 0 : AliError("No Digit!!");
1085 0 : continue;
1086 : }
1087 :
1088 97 : Int_t nparents = digit->GetNiparent() ;
1089 150 : if ( nparents == 0 ) continue ;
1090 :
1091 : // Loop on all parents in digit
1092 : Int_t jndex ;
1093 178 : for ( jndex = 0 ; jndex < nparents ; jndex++ )
1094 : {
1095 45 : if ( fMulParent > fMaxParent )
1096 : {
1097 0 : fMulTrack = - 1 ;
1098 0 : Error("EvalParents", "increase fMaxParent") ;
1099 0 : break ;
1100 : }
1101 :
1102 45 : Int_t newParent = digit->GetIparent (jndex+1) ;
1103 45 : Float_t newdEParent = digit->GetDEParent(jndex+1) ;
1104 : // Check if parent was not already stored
1105 : // if stored, add its energy deposition.
1106 : Int_t kndex ;
1107 : Bool_t already = kFALSE ;
1108 92 : for ( kndex = 0 ; kndex < fMulParent ; kndex++ )
1109 : {
1110 28 : if ( newParent == parentArray[kndex] )
1111 : {
1112 27 : dEParentArray[kndex] += newdEParent;
1113 : already = kTRUE ;
1114 27 : break ;
1115 : }
1116 : } // end of check
1117 :
1118 : // Store the parent
1119 63 : if ( !already && (fMulParent < fMaxParent) )
1120 : {
1121 18 : parentArray [fMulParent] = newParent ;
1122 18 : dEParentArray[fMulParent] = newdEParent ;
1123 18 : fMulParent++ ;
1124 18 : } // store it
1125 : } // all parents in digit
1126 44 : } // all digits
1127 :
1128 : // Order the parents from highest to lowest energy deposit
1129 33 : if ( fMulParent > 0 )
1130 : {
1131 17 : Int_t *sortIdx = new Int_t[fMulParent];
1132 17 : TMath::Sort(fMulParent,dEParentArray,sortIdx);
1133 :
1134 70 : for(index = 0; index < fMulParent; index++)
1135 : {
1136 18 : fParentsList [index] = parentArray [sortIdx[index]] ;
1137 18 : fDEParentsList[index] = dEParentArray[sortIdx[index]] ;
1138 : }
1139 34 : delete [] sortIdx;
1140 17 : }
1141 :
1142 66 : delete [] parentArray;
1143 66 : delete [] dEParentArray;
1144 33 : }
1145 :
1146 : //____________________________________________________________________________
1147 : void AliEMCALRecPoint::GetLocalPosition(TVector3 & lpos) const
1148 : {
1149 : // returns the position of the cluster in the local reference system
1150 : // of the sub-detector
1151 :
1152 272 : lpos = fLocPos;
1153 136 : }
1154 :
1155 : //____________________________________________________________________________
1156 : void AliEMCALRecPoint::GetGlobalPosition(TVector3 & gpos) const
1157 : {
1158 : // returns the position of the cluster in the global reference system of ALICE
1159 : // These are now the Cartesian X, Y and Z
1160 : // cout<<" geom "<<geom<<endl;
1161 : // fGeomPtr->GetGlobal(fLocPos, gpos, fSuperModuleNumber);
1162 132 : gpos = fGlobPos;
1163 :
1164 66 : }
1165 :
1166 : //____________________________________________________________________________
1167 : //void AliEMCALRecPoint::GetGlobalPosition(TVector3 & gpos, TMatrixF & gmat) const
1168 : //{
1169 : // // returns the position of the cluster in the global reference system of ALICE
1170 : // // These are now the Cartesian X, Y and Z
1171 : // // cout<<" geom "<<geom<<endl;
1172 : //
1173 : // //To be implemented
1174 : // fGeomPtr->GetGlobalEMCAL(this, gpos, gmat);
1175 : //
1176 : //}
1177 :
1178 : //_____________________________________________________________________________
1179 : void AliEMCALRecPoint::EvalLocal2TrackingCSTransform()
1180 : {
1181 : //Evaluates local to "tracking" c.s. transformation (B.P.).
1182 : //All evaluations should be completed before calling for this
1183 : //function.
1184 : //See ALICE PPR Chapter 5 p.18 for "tracking" c.s. definition,
1185 : //or just ask Jouri Belikov. :)
1186 :
1187 66 : SetVolumeId(AliGeomManager::LayerToVolUID(AliGeomManager::kEMCAL,GetSuperModuleNumber()));
1188 :
1189 33 : const TGeoHMatrix* tr2loc = GetTracking2LocalMatrix();
1190 33 : if(!tr2loc) AliFatal(Form("No Tracking2LocalMatrix found."));
1191 :
1192 33 : Double_t lxyz[3] = {fLocPos.X(),fLocPos.Y(),fLocPos.Z()};
1193 33 : Double_t txyz[3] = {0,0,0};
1194 :
1195 33 : tr2loc->MasterToLocal(lxyz,txyz);
1196 33 : SetX(txyz[0]); SetY(txyz[1]); SetZ(txyz[2]);
1197 :
1198 33 : if(AliLog::GetGlobalDebugLevel()>0) {
1199 0 : TVector3 gpos; //TMatrixF gmat;
1200 : //GetGlobalPosition(gpos,gmat); //Not doing anythin special, replace by next line.
1201 0 : fGeomPtr->GetGlobal(fLocPos, gpos, GetSuperModuleNumber());
1202 :
1203 0 : Float_t gxyz[3];
1204 0 : GetGlobalXYZ(gxyz);
1205 0 : AliInfo(Form("lCS-->(%.3f,%.3f,%.3f), tCS-->(%.3f,%.3f,%.3f), gCS-->(%.3f,%.3f,%.3f), gCScalc-\
1206 : ->(%.3f,%.3f,%.3f), supermodule %d",
1207 : fLocPos.X(),fLocPos.Y(),fLocPos.Z(),
1208 : GetX(),GetY(),GetZ(),
1209 : gpos.X(),gpos.Y(),gpos.Z(),
1210 : gxyz[0],gxyz[1],gxyz[2],GetSuperModuleNumber()));
1211 0 : }
1212 :
1213 33 : }
1214 :
1215 : //____________________________________________________________________________
1216 : Float_t AliEMCALRecPoint::GetMaximalEnergy(void) const
1217 : {
1218 : // Finds the maximum energy in the cluster
1219 :
1220 : Float_t menergy = 0. ;
1221 :
1222 : Int_t iDigit;
1223 0 : for(iDigit=0; iDigit<fMulDigit; iDigit++) {
1224 :
1225 0 : if(fEnergyList[iDigit] > menergy)
1226 0 : menergy = fEnergyList[iDigit] ;
1227 : }
1228 0 : return menergy ;
1229 : }
1230 :
1231 : //____________________________________________________________________________
1232 : Int_t AliEMCALRecPoint::GetMaximalEnergyIndex(void) const
1233 : {
1234 : // Finds the maximum energy in the cluster
1235 :
1236 : Float_t menergy = 0. ;
1237 : Int_t mid = 0 ;
1238 : Int_t iDigit;
1239 :
1240 293 : for(iDigit=0; iDigit<fMulDigit; iDigit++) {
1241 :
1242 97 : if(fEnergyList[iDigit] > menergy){
1243 : menergy = fEnergyList[iDigit] ;
1244 : mid = iDigit ;
1245 44 : }
1246 : }//loop on cluster digits
1247 :
1248 33 : return mid ;
1249 : }
1250 :
1251 :
1252 : //____________________________________________________________________________
1253 : Int_t AliEMCALRecPoint::GetMultiplicityAtLevel(Float_t H) const
1254 : {
1255 : // Calculates the multiplicity of digits with energy larger than H*energy
1256 :
1257 : Int_t multipl = 0 ;
1258 : Int_t iDigit ;
1259 0 : for(iDigit=0; iDigit<fMulDigit; iDigit++) {
1260 :
1261 0 : if(fEnergyList[iDigit] > H * fAmp)
1262 0 : multipl++ ;
1263 : }
1264 0 : return multipl ;
1265 : }
1266 :
1267 :
1268 : //____________________________________________________________________________
1269 : ///
1270 : /// Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
1271 : /// energy difference between two local maxima.
1272 : /// Called by AliEMCALClusterizerv1.
1273 : ///
1274 : /// \return Number of local maxima in cluster
1275 : /// \param nMult : Number digits in cluster
1276 : /// \param locMaxCut : Minimum energy difference between local maximum and neighbour towers.
1277 : /// \param digits : List of digits in cluster
1278 : ///
1279 : Int_t AliEMCALRecPoint::GetNumberOfLocalMax(Int_t nMult,
1280 : Float_t locMaxCut,
1281 : TClonesArray * digits) const
1282 : {
1283 66 : AliEMCALDigit ** maxAt = new AliEMCALDigit*[nMult] ;
1284 33 : Float_t * maxAtEnergy = new Float_t [nMult] ;
1285 :
1286 33 : Int_t nMax = GetNumberOfLocalMax(maxAt, maxAtEnergy, locMaxCut, digits);
1287 :
1288 66 : delete[] maxAt ;
1289 66 : delete[] maxAtEnergy ;
1290 :
1291 33 : return nMax;
1292 : }
1293 :
1294 : //____________________________________________________________________________
1295 : ///
1296 : /// Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
1297 : /// energy difference between two local maxima.
1298 : /// Called by AliEMCALUnfolding.
1299 : ///
1300 : /// \return Number of local maxima in cluster
1301 : /// \param maxAt : list of local maxima digits
1302 : /// \param maxAtEnergy : list of local maxima energies
1303 : /// \param locMaxCut : Minimum energy difference between local maximum and neighbour towers.
1304 : /// \param digits : List of digits in cluster
1305 : ///
1306 : Int_t AliEMCALRecPoint::GetNumberOfLocalMax(AliEMCALDigit ** maxAt, Float_t * maxAtEnergy,
1307 : Float_t locMaxCut,TClonesArray * digits) const
1308 : {
1309 : AliEMCALDigit * digit = 0;
1310 : AliEMCALDigit * digitN = 0;
1311 :
1312 : Int_t iDigitN = 0 ;
1313 : Int_t iDigit = 0 ;
1314 :
1315 293 : for(iDigit = 0; iDigit < fMulDigit; iDigit++)
1316 : {
1317 97 : maxAt[iDigit] = (AliEMCALDigit*) digits->At(fDigitsList[iDigit]) ;
1318 : }
1319 :
1320 260 : for(iDigit = 0 ; iDigit < fMulDigit; iDigit++)
1321 : {
1322 97 : if(maxAt[iDigit])
1323 : {
1324 : digit = maxAt[iDigit] ;
1325 :
1326 812 : for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++)
1327 : {
1328 350 : if(iDigitN == iDigit) continue;//the same digit
1329 :
1330 294 : digitN = (AliEMCALDigit *) digits->At(fDigitsList[iDigitN]) ;
1331 :
1332 294 : if ( AreNeighbours(digit, digitN) )
1333 : {
1334 89 : if (fEnergyList[iDigit] > fEnergyList[iDigitN] )
1335 : {
1336 57 : maxAt[iDigitN] = 0 ;
1337 :
1338 : // but may be digit too is not local max ?
1339 57 : if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut)
1340 1 : maxAt[iDigit] = 0 ;
1341 : }
1342 : else
1343 : {
1344 32 : maxAt[iDigit] = 0 ;
1345 :
1346 : // but may be digitN too is not local max ?
1347 32 : if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut)
1348 4 : maxAt[iDigitN] = 0 ;
1349 : }
1350 : } // if Are neighbours
1351 : } // while digitN
1352 : } // slot not empty
1353 : } // while digit
1354 :
1355 : iDigitN = 0 ;
1356 260 : for(iDigit = 0; iDigit < fMulDigit; iDigit++)
1357 : {
1358 97 : if(maxAt[iDigit] )
1359 : {
1360 35 : maxAt[iDigitN] = maxAt[iDigit] ;
1361 35 : maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
1362 35 : iDigitN++ ;
1363 35 : }
1364 : }
1365 :
1366 33 : return iDigitN ;
1367 : }
1368 :
1369 : //____________________________________________________________________________
1370 : Int_t AliEMCALRecPoint::GetPrimaryIndex() const
1371 : {
1372 : // Get the primary track index in TreeK which deposits the most energy
1373 : // in Digits which forms RecPoint.
1374 :
1375 0 : if (fMulTrack)
1376 0 : return fTracksList[0];
1377 0 : return -12345;
1378 0 : }
1379 :
1380 : //____________________________________________________________________________
1381 : void AliEMCALRecPoint::EvalTime(TClonesArray * digits){
1382 : // time is set to the time of the digit with the maximum energy
1383 :
1384 : Float_t maxE = 0;
1385 : Int_t maxAt = 0;
1386 293 : for(Int_t idig=0; idig < fMulDigit; idig++){
1387 97 : if(fEnergyList[idig] > maxE){
1388 : maxE = fEnergyList[idig] ;
1389 : maxAt = idig;
1390 44 : }
1391 : }
1392 33 : fTime = ((AliEMCALDigit*) digits->At(fDigitsList[maxAt]))->GetTime() ;
1393 :
1394 33 : }
1395 :
1396 : //______________________________________________________________________________
1397 : void AliEMCALRecPoint::Paint(Option_t *)
1398 : {
1399 : // Paint this ALiRecPoint as a TMarker with its current attributes
1400 :
1401 0 : TVector3 pos(0.,0.,0.) ;
1402 0 : GetLocalPosition(pos) ;
1403 0 : Coord_t x = pos.X() ;
1404 0 : Coord_t y = pos.Z() ;
1405 : Color_t markercolor = 1 ;
1406 : Size_t markersize = 1.;
1407 : Style_t markerstyle = 5 ;
1408 :
1409 0 : if (!gPad->IsBatch()) {
1410 0 : gVirtualX->SetMarkerColor(markercolor) ;
1411 0 : gVirtualX->SetMarkerSize (markersize) ;
1412 0 : gVirtualX->SetMarkerStyle(markerstyle) ;
1413 : }
1414 0 : gPad->SetAttMarkerPS(markercolor,markerstyle,markersize) ;
1415 0 : gPad->PaintPolyMarker(1,&x,&y,"") ;
1416 0 : }
1417 :
1418 : //_____________________________________________________________________
1419 : Double_t AliEMCALRecPoint::TmaxInCm(const Double_t e , const Int_t key)
1420 : {
1421 : // e energy in GeV)
1422 : // key = 0(gamma, default)
1423 : // != 0(electron)
1424 : const Double_t ca = 4.82; // shower max parameter - first guess; ca=TMath::Log(1000./8.07)
1425 : Double_t tmax = 0.; // position of electromagnetic shower max in cm
1426 :
1427 : Double_t x0 = 1.31; // radiation lenght (cm)
1428 : //If old geometry in use
1429 198 : if(!((fGeomPtr->GetEMCGeometry()->GetGeoName()).Contains("V1"))) x0 = 1.28;
1430 :
1431 66 : if(e>0.1) {
1432 66 : tmax = TMath::Log(e) + ca;
1433 132 : if (key==0) tmax += 0.5;
1434 0 : else tmax -= 0.5;
1435 66 : tmax *= x0; // convert to cm
1436 66 : }
1437 66 : return tmax;
1438 0 : }
1439 :
1440 : //______________________________________________________________________________
1441 : Float_t AliEMCALRecPoint::EtaToTheta(Float_t arg) const
1442 : {
1443 : //Converts Theta (Radians) to Eta(Radians)
1444 0 : return (2.*TMath::ATan(TMath::Exp(-arg)));
1445 : }
1446 :
1447 : //______________________________________________________________________________
1448 : Float_t AliEMCALRecPoint::ThetaToEta(Float_t arg) const
1449 : {
1450 : //Converts Eta (Radians) to Theta(Radians)
1451 0 : return (-1 * TMath::Log(TMath::Tan(0.5 * arg)));
1452 : }
1453 :
1454 : //____________________________________________________________________________
1455 : void AliEMCALRecPoint::Print(Option_t *opt) const
1456 : {
1457 : // Print the list of digits belonging to the cluster
1458 66 : if(strlen(opt)==0) return;
1459 0 : TString message ;
1460 0 : message = "AliEMCALRecPoint:\n" ;
1461 0 : message += " digits # = " ;
1462 0 : AliInfo(message.Data()) ;
1463 :
1464 : Int_t iDigit;
1465 0 : for(iDigit=0; iDigit<fMulDigit; iDigit++)
1466 0 : printf(" %d ", fDigitsList[iDigit] ) ;
1467 0 : printf("\n");
1468 :
1469 0 : AliInfo(" Energies = ") ;
1470 0 : for(iDigit=0; iDigit<fMulDigit; iDigit++)
1471 0 : printf(" %f ", fEnergyList[iDigit] ) ;
1472 0 : printf("\n");
1473 :
1474 0 : AliInfo("\n Abs Ids = ") ;
1475 0 : for(iDigit=0; iDigit<fMulDigit; iDigit++)
1476 0 : printf(" %i ", fAbsIdList[iDigit] ) ;
1477 0 : printf("\n");
1478 :
1479 0 : AliInfo(" Primaries ") ;
1480 0 : for(iDigit = 0;iDigit < fMulTrack; iDigit++)
1481 0 : printf(" %d ", fTracksList[iDigit]) ;
1482 :
1483 0 : printf("\n Local x %6.2f y %7.2f z %7.1f \n", fLocPos[0], fLocPos[1], fLocPos[2]);
1484 :
1485 0 : message = " ClusterType = %d" ;
1486 0 : message += " Multiplicity = %d" ;
1487 0 : message += " Cluster Energy = %f" ;
1488 0 : message += " Core energy = %f" ;
1489 0 : message += " Core radius = %f" ;
1490 0 : message += " Number of primaries %d" ;
1491 0 : message += " Stored at position %d" ;
1492 0 : AliInfo(Form(message.Data(), fClusterType, fMulDigit, fAmp, fCoreEnergy, fCoreRadius, fMulTrack, GetIndexInList()) ) ;
1493 33 : }
1494 :
1495 : //___________________________________________________________
1496 : Double_t AliEMCALRecPoint::GetPointEnergy() const
1497 : {
1498 : //Returns energy ....
1499 : Double_t e=0.0;
1500 0 : for(int ic=0; ic<GetMultiplicity(); ic++) e += double(fEnergyList[ic]);
1501 0 : return e;
1502 : }
|