]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALRecPoint.cxx
Small syntax corrections for IBM xlc on the Mac
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALRecPoint.cxx
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 //*-- Author: Yves Schutz (SUBATECH)
20 //*-- Author: Dmitri Peressounko (RRC KI & SUBATECH)
21 //*-- Author: Heather Gray (LBL) merged AliEMCALRecPoint and AliEMCALTowerRecPoint 02/04
22
23 // --- ROOT system ---
24 #include "TPad.h"
25 #include "TGraph.h"
26 #include "TPaveText.h"
27 #include "TClonesArray.h"
28 #include "TMath.h" 
29
30 // --- Standard library ---
31
32 // --- AliRoot header files ---
33 #include "AliGenerator.h"
34 #include "AliEMCALGeometry.h"
35 #include "AliEMCALDigit.h"
36 #include "AliEMCALRecPoint.h"
37 #include "AliEMCALGetter.h"
38
39 ClassImp(AliEMCALRecPoint)
40
41
42 //____________________________________________________________________________
43 AliEMCALRecPoint::AliEMCALRecPoint()
44   : AliRecPoint()
45 {
46   // ctor
47   fMaxTrack = 0 ;
48   fMulDigit   = 0 ;  
49   fAmp   = 0. ;   
50   fCoreEnergy = 0 ; 
51   fEnergyList = 0 ;
52   fTime = 0. ;
53   fLocPos.SetX(0.)  ;      //Local position should be evaluated
54   fCoreRadius = 10;        //HG Check this
55 }
56
57 //____________________________________________________________________________
58 AliEMCALRecPoint::AliEMCALRecPoint(const char * opt) : AliRecPoint(opt)
59 {
60   // ctor
61   fMaxTrack = 200 ;
62   fMulDigit   = 0 ;  
63   fAmp   = 0. ;   
64   fCoreEnergy = 0 ; 
65   fEnergyList = 0 ;
66   fTime = -1. ;
67   fLocPos.SetX(1000000.)  ;      //Local position should be evaluated
68   fCoreRadius = 10;        //HG Check this
69 }
70 //____________________________________________________________________________
71 AliEMCALRecPoint::~AliEMCALRecPoint()
72 {
73   // dtor
74   if ( fEnergyList )
75     delete[] fEnergyList ; 
76 }
77
78 //____________________________________________________________________________
79 void AliEMCALRecPoint::AddDigit(AliEMCALDigit & digit, Float_t Energy)
80 {
81   // Adds a digit to the RecPoint
82   // and accumulates the total amplitude and the multiplicity 
83   
84   if(fEnergyList == 0)
85     fEnergyList =  new Float_t[fMaxDigit]; 
86
87   if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists 
88     fMaxDigit*=2 ; 
89     Int_t * tempo = new Int_t[fMaxDigit]; 
90     Float_t * tempoE =  new Float_t[fMaxDigit];
91
92     Int_t index ;     
93     for ( index = 0 ; index < fMulDigit ; index++ ){
94       tempo[index]  = fDigitsList[index] ;
95       tempoE[index] = fEnergyList[index] ; 
96     }
97     
98     delete [] fDigitsList ; 
99     fDigitsList =  new Int_t[fMaxDigit];
100  
101     delete [] fEnergyList ;
102     fEnergyList =  new Float_t[fMaxDigit];
103
104     for ( index = 0 ; index < fMulDigit ; index++ ){
105       fDigitsList[index] = tempo[index] ;
106       fEnergyList[index] = tempoE[index] ; 
107     }
108  
109     delete [] tempo ;
110     delete [] tempoE ; 
111   } // if
112   
113   fDigitsList[fMulDigit]   = digit.GetIndexInList()  ; 
114   fEnergyList[fMulDigit]   = Energy ;
115   fMulDigit++ ; 
116   fAmp += Energy ; 
117
118 }
119 //____________________________________________________________________________
120 Bool_t AliEMCALRecPoint::AreNeighbours(AliEMCALDigit * digit1, AliEMCALDigit * digit2 ) const
121 {
122   // Tells if (true) or not (false) two digits are neighbours
123   // A neighbour is defined as being two digits which share a corner
124   
125   Bool_t areNeighbours = kFALSE ;
126   
127   AliEMCALGeometry * geom =  (AliEMCALGetter::Instance())->EMCALGeometry();
128
129   Int_t relid1[2] ; 
130   geom->AbsToRelNumbering(digit1->GetId(), relid1) ; 
131
132   Int_t relid2[2] ; 
133   geom->AbsToRelNumbering(digit2->GetId(), relid2) ; 
134   
135   Int_t rowdiff = TMath::Abs( relid1[0] - relid2[0] ) ;  
136   Int_t coldiff = TMath::Abs( relid1[1] - relid2[1] ) ;  
137
138   if (( coldiff <= 1 )  && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0)) 
139     areNeighbours = kTRUE ;
140   
141   return areNeighbours;
142 }
143
144 //____________________________________________________________________________
145 Int_t AliEMCALRecPoint::Compare(const TObject * obj) const
146 {
147   // Compares two RecPoints according to their position in the EMCAL modules
148
149   Float_t delta = 1 ; //Width of "Sorting row". If you change this 
150                       //value (what is senseless) change as well delta in
151                       //AliEMCALTrackSegmentMakerv* and other RecPoints...
152   Int_t rv ; 
153
154   AliEMCALRecPoint * clu = (AliEMCALRecPoint *)obj ; 
155
156   TVector3 locpos1; 
157   GetLocalPosition(locpos1);
158   TVector3 locpos2;  
159   clu->GetLocalPosition(locpos2);  
160
161   Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
162   if (rowdif> 0) 
163     rv = 1 ;
164   else if(rowdif < 0) 
165     rv = -1 ;
166   else if(locpos1.Y()>locpos2.Y()) 
167     rv = -1 ;
168   else 
169     rv = 1 ; 
170
171   return rv ; 
172 }
173
174 //____________________________________________________________________________
175 Int_t AliEMCALRecPoint::DistancetoPrimitive(Int_t px, Int_t py)
176 {
177   // Compute distance from point px,py to  a AliEMCALRecPoint considered as a Tmarker
178   // Compute the closest distance of approach from point px,py to this marker.
179   // The distance is computed in pixels units.
180   // HG Still need to update -> Not sure what this should achieve
181
182   TVector3 pos(0.,0.,0.) ;
183   GetLocalPosition(pos) ;
184   Float_t x =  pos.X() ;
185   Float_t y =  pos.Y() ;
186   const Int_t kMaxDiff = 10;
187   Int_t pxm  = gPad->XtoAbsPixel(x);
188   Int_t pym  = gPad->YtoAbsPixel(y);
189   Int_t dist = (px-pxm)*(px-pxm) + (py-pym)*(py-pym);
190   
191   if (dist > kMaxDiff) return 9999;
192   return dist;
193 }
194
195 //___________________________________________________________________________
196  void AliEMCALRecPoint::Draw(Option_t *option)
197  {
198    // Draw this AliEMCALRecPoint with its current attributes
199    
200    AppendPad(option);
201  }
202
203 //______________________________________________________________________________
204 void AliEMCALRecPoint::ExecuteEvent(Int_t /*event*/, Int_t, Int_t)
205 {
206   // Execute action corresponding to one event
207   // This member function is called when a AliEMCALRecPoint is clicked with the locator
208   //
209   // If Left button is clicked on AliEMCALRecPoint, the digits are switched on    
210   // and switched off when the mouse button is released.
211
212   //  static Int_t pxold, pyold;
213
214   /* static TGraph *  digitgraph = 0 ;
215   static TPaveText* clustertext = 0 ;
216   
217   if (!gPad->IsEditable()) return;
218   
219   switch (event) {
220     
221     
222   case kButton1Down:{
223     AliEMCALDigit * digit ;
224     AliEMCALGeometry * emcalgeom =  (AliEMCALGetter::Instance())->EMCALGeometry() ;
225
226     Int_t iDigit;
227     Int_t relid[2] ;
228   
229     const Int_t kMulDigit=AliEMCALRecPoint::GetDigitsMultiplicity() ;
230     Float_t * xi = new Float_t [kMulDigit] ; 
231     Float_t * zi = new Float_t [kMulDigit] ;
232     
233     for(iDigit = 0; iDigit < kMulDigit; iDigit++) {
234       Fatal("AliEMCALRecPoint::ExecuteEvent", " -> Something wrong with the code"); 
235       digit = 0 ; //dynamic_cast<AliEMCALDigit *>((fDigitsList)[iDigit]);
236       emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
237       emcalgeom->PosInAlice(relid, xi[iDigit], zi[iDigit]) ;
238     }
239     
240     if (!digitgraph) {
241       digitgraph = new TGraph(fMulDigit,xi,zi);
242       digitgraph-> SetMarkerStyle(5) ; 
243       digitgraph-> SetMarkerSize(1.) ;
244       digitgraph-> SetMarkerColor(1) ;
245       digitgraph-> Draw("P") ;
246     }
247     if (!clustertext) {
248       
249       TVector3 pos(0.,0.,0.) ;
250       GetLocalPosition(pos) ;
251       clustertext = new TPaveText(pos.X()-10,pos.Z()+10,pos.X()+50,pos.Z()+35,"") ;
252       Text_t  line1[40] ;
253       Text_t  line2[40] ;
254       sprintf(line1,"Energy=%1.2f GeV",GetEnergy()) ;
255       sprintf(line2,"%d Digits",GetDigitsMultiplicity()) ;
256       clustertext ->AddText(line1) ;
257       clustertext ->AddText(line2) ;
258       clustertext ->Draw("");
259     }
260     gPad->Update() ; 
261     Print("") ;
262     delete[] xi ; 
263     delete[] zi ; 
264    }
265   
266 break;
267   
268   case kButton1Up:
269     if (digitgraph) {
270       delete digitgraph  ;
271       digitgraph = 0 ;
272     }
273     if (clustertext) {
274       delete clustertext ;
275       clustertext = 0 ;
276     }
277     
278     break;
279     
280     }*/
281 }
282 //____________________________________________________________________________
283 void AliEMCALRecPoint::EvalAll(Float_t logWeight,TClonesArray * digits) 
284 {
285   // Evaluates all shower parameters
286
287   EvalLocalPosition(logWeight, digits) ;
288   EvalElipsAxis(logWeight, digits) ;
289   EvalDispersion(logWeight, digits) ;
290   EvalCoreEnergy(logWeight, digits);
291   EvalTime(digits) ;
292
293   //EvalPrimaries(digits) ;
294 }
295
296 //____________________________________________________________________________
297 void  AliEMCALRecPoint::EvalDispersion(Float_t logWeight, TClonesArray * digits)
298 {
299   // Calculates the dispersion of the shower at the origin of the RecPoint
300
301   Float_t d    = 0. ;
302   Float_t wtot = 0. ;
303
304   AliEMCALDigit * digit ;
305  
306   AliEMCALGeometry * geom = (AliEMCALGetter::Instance())->EMCALGeometry();
307   
308   // Calculates the centre of gravity in the local EMCAL-module coordinates   
309   Int_t iDigit;
310
311   if (!fLocPos.X() || !fLocPos.Y()) 
312     EvalLocalPosition(logWeight, digits) ;
313   
314   const Float_t kDeg2Rad = TMath::DegToRad() ; 
315     
316   Float_t cluEta =  fLocPos.X() ; 
317   Float_t cluPhi =  fLocPos.Y() ; 
318   Float_t cluR =  fLocPos.Z() ; 
319   
320   if (gDebug == 2) 
321     printf("EvalDispersion: eta,phi,r = %f,%f,%f", cluEta, cluPhi, cluR) ;
322   
323   // Calculates the dispersion in coordinates 
324   wtot = 0.;
325   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
326     digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit])  ;
327     Float_t etai = 0.;
328     Float_t phii = 0.;
329     geom->EtaPhiFromIndex(digit->GetId(), etai, phii);
330     phii = phii * kDeg2Rad;
331     if (gDebug == 2) 
332       printf("EvalDispersion: id = %d, etai,phii = %f,%f", digit->GetId(), etai, phii) ;
333
334     Float_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
335     d += w * ( (etai-cluEta)*(etai-cluEta) + (phii-cluPhi)*(phii-cluPhi) ) ; 
336     wtot+=w ;
337   }
338   
339   if ( wtot > 0 ) 
340     d /= wtot ;
341   else 
342     d = 0. ; 
343
344   fDispersion = TMath::Sqrt(d) ;
345  
346 }
347
348 //____________________________________________________________________________
349 void AliEMCALRecPoint::EvalLocalPosition(Float_t logWeight, TClonesArray * digits)
350 {
351   // Calculates the center of gravity in the local EMCAL-module coordinates 
352   Float_t wtot = 0. ;
353  
354   //  Int_t relid[3] ;
355   
356   AliEMCALDigit * digit ;
357   AliEMCALGeometry * geom  =  (AliEMCALGetter::Instance())->EMCALGeometry();
358   Int_t iDigit;
359   Float_t cluEta = 0;
360   Float_t cluPhi = 0;
361   const Float_t kDeg2Rad = TMath::DegToRad();
362
363   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
364     digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
365
366     Float_t etai ;
367     Float_t phii ;
368     geom->EtaPhiFromIndex(digit->GetId(), etai, phii);
369     phii = phii * kDeg2Rad;
370     Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
371     cluEta += (etai * w) ;
372     cluPhi += (phii * w );
373     wtot += w ;
374   }
375
376   if ( wtot > 0 ) { 
377     cluEta /= wtot ;
378     cluPhi /= wtot ;
379   } else {
380     cluEta = -1 ; 
381     cluPhi = -1.; 
382   }
383   
384   fLocPos.SetX(cluEta);
385   fLocPos.SetY(cluPhi);
386   fLocPos.SetZ(geom->GetIP2ECASection());
387     
388   if (gDebug==2)
389     printf("EvalLocalPosition: eta,phi,r = %f,%f,%f", fLocPos.X(), fLocPos.Y(), fLocPos.Z()) ; 
390   fLocPosM = 0 ;
391 }
392
393 //______________________________________________________________________________
394 void AliEMCALRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits)
395 {
396   // This function calculates energy in the core, 
397   // i.e. within a radius rad = 3cm around the center. Beyond this radius
398   // in accordance with shower profile the energy deposition 
399   // should be less than 2%
400
401   AliEMCALDigit * digit ;
402   const Float_t kDeg2Rad = TMath::DegToRad() ;
403   AliEMCALGeometry * geom = (AliEMCALGetter::Instance())->EMCALGeometry();    
404   Int_t iDigit;
405
406   if (!fLocPos.X() || !fLocPos.Y() ) {
407     EvalLocalPosition(logWeight, digits);
408   }
409   
410   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
411     digit = (AliEMCALDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
412     Float_t etai = 0. ;
413     Float_t phii = 0. ;
414     geom->PosInAlice(digit->GetId(), etai, phii);
415     phii = phii * kDeg2Rad;
416   
417     Float_t distance = TMath::Sqrt((etai-fLocPos.X())*(etai-fLocPos.X())+(phii-fLocPos.Y())*(phii-fLocPos.Y())) ;
418     if(distance < fCoreRadius)
419       fCoreEnergy += fEnergyList[iDigit] ;
420   }
421   
422 }
423 //____________________________________________________________________________
424 void  AliEMCALRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
425 {
426   // Calculates the axis of the shower ellipsoid in eta and phi
427
428   Double_t wtot = 0. ;
429   Double_t x    = 0.;
430   Double_t z    = 0.;
431   Double_t dxx  = 0.;
432   Double_t dzz  = 0.;
433   Double_t dxz  = 0.;
434
435   AliEMCALDigit * digit ;
436
437   AliEMCALGeometry * geom = (AliEMCALGetter::Instance())->EMCALGeometry();
438
439   Int_t iDigit;
440   
441   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
442     digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit])  ;
443     Float_t etai = 0. ;
444     Float_t phii = 0. ; 
445     geom->EtaPhiFromIndex(digit->GetId(), etai, phii);
446     Double_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
447     dxx  += w * etai * etai ;
448     x    += w * etai ;
449     dzz  += w * phii * phii ;
450     z    += w * phii ; 
451     dxz  += w * etai * etai ; 
452     wtot += w ;
453   }
454   if ( wtot > 0 ) { 
455     dxx /= wtot ;
456     x   /= wtot ;
457     dxx -= x * x ;
458     dzz /= wtot ;
459     z   /= wtot ;
460     dzz -= z * z ;
461     dxz /= wtot ;
462     dxz -= x * z ;
463
464     fLambda[0] =  0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
465     if(fLambda[0] > 0)
466       fLambda[0] = TMath::Sqrt(fLambda[0]) ;
467     else
468       fLambda[0] = 0;
469     
470     fLambda[1] =  0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
471     if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
472       fLambda[1] = TMath::Sqrt(fLambda[1]) ;
473     else
474       fLambda[1]= 0. ;
475   } else { 
476     fLambda[0]= 0. ;
477     fLambda[1]= 0. ;
478   }
479 }
480
481 //______________________________________________________________________________
482 void  AliEMCALRecPoint::EvalPrimaries(TClonesArray * digits)
483 {
484   // Constructs the list of primary particles (tracks) which have contributed to this RecPoint
485   
486   AliEMCALDigit * digit ;
487   Int_t * tempo    = new Int_t[fMaxTrack] ;
488
489   Int_t index ;  
490   for ( index = 0 ; index < GetDigitsMultiplicity() ; index++ ) { // all digits
491     digit = dynamic_cast<AliEMCALDigit *>(digits->At( fDigitsList[index] )) ; 
492     Int_t nprimaries = digit->GetNprimary() ;
493     Int_t * newprimaryarray = new Int_t[nprimaries] ;
494     Int_t ii ; 
495     for ( ii = 0 ; ii < nprimaries ; ii++)
496       newprimaryarray[ii] = digit->GetPrimary(ii+1) ; 
497
498     Int_t jndex ;
499     for ( jndex = 0 ; jndex < nprimaries ; jndex++ ) { // all primaries in digit
500       if ( fMulTrack > fMaxTrack ) {
501         fMulTrack = - 1 ;
502         Error("GetNprimaries", "increase fMaxTrack ")  ;
503         break ;
504       }
505       Int_t newprimary = newprimaryarray[jndex] ;
506       Int_t kndex ;
507       Bool_t already = kFALSE ;
508       for ( kndex = 0 ; kndex < fMulTrack ; kndex++ ) { //check if not already stored
509         if ( newprimary == tempo[kndex] ){
510           already = kTRUE ;
511           break ;
512         }
513       } // end of check
514       if ( !already) { // store it
515         tempo[fMulTrack] = newprimary ; 
516         fMulTrack++ ;
517       } // store it
518     } // all primaries in digit
519     delete newprimaryarray ; 
520   } // all digits
521
522   
523   fTracksList = new Int_t[fMulTrack] ;
524   for(index = 0; index < fMulTrack; index++)
525    fTracksList[index] = tempo[index] ;
526  
527   delete tempo ;
528
529 }
530
531 //____________________________________________________________________________
532 void AliEMCALRecPoint::GetLocalPosition(TVector3 & lpos) const
533 {
534   // returns the position of the cluster in the local reference system of ALICE
535   // X = eta, Y = phi, Z = r (a constant for the EMCAL)
536   
537   lpos.SetX(fLocPos.X()) ;
538   lpos.SetY(fLocPos.Y()) ;
539   lpos.SetZ(fLocPos.Z()) ;
540 }
541
542 //____________________________________________________________________________
543 void AliEMCALRecPoint::GetGlobalPosition(TVector3 & gpos) const
544 {
545   // returns the position of the cluster in the global reference system of ALICE
546   // These are now the Cartesian X, Y and Z
547
548   AliEMCALGeometry * geom = (AliEMCALGetter::Instance())->EMCALGeometry();
549   Int_t absid = geom->TowerIndexFromEtaPhi(fLocPos.X(), TMath::RadToDeg()*fLocPos.Y());
550   geom->XYZFromIndex(absid, gpos);
551 }
552
553 //____________________________________________________________________________
554 Float_t AliEMCALRecPoint::GetMaximalEnergy(void) const
555 {
556   // Finds the maximum energy in the cluster
557   
558   Float_t menergy = 0. ;
559
560   Int_t iDigit;
561
562   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
563  
564     if(fEnergyList[iDigit] > menergy) 
565       menergy = fEnergyList[iDigit] ;
566   }
567   return menergy ;
568 }
569
570 //____________________________________________________________________________
571 Int_t AliEMCALRecPoint::GetMultiplicityAtLevel(Float_t H) const
572 {
573   // Calculates the multiplicity of digits with energy larger than H*energy 
574   
575   Int_t multipl   = 0 ;
576   Int_t iDigit ;
577   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
578
579     if(fEnergyList[iDigit] > H * fAmp) 
580       multipl++ ;
581   }
582   return multipl ;
583 }
584
585 //____________________________________________________________________________
586 Int_t  AliEMCALRecPoint::GetNumberOfLocalMax(AliEMCALDigit **  maxAt, Float_t * maxAtEnergy,
587                                                Float_t locMaxCut,TClonesArray * digits) const
588
589   // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
590   // energy difference between two local maxima
591
592   AliEMCALDigit * digit ;
593   AliEMCALDigit * digitN ;
594   
595   Int_t iDigitN ;
596   Int_t iDigit ;
597
598   for(iDigit = 0; iDigit < fMulDigit; iDigit++)
599     maxAt[iDigit] = (AliEMCALDigit*) digits->At(fDigitsList[iDigit])  ;
600   
601   for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {   
602     if(maxAt[iDigit]) {
603       digit = maxAt[iDigit] ;
604           
605       for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {        
606         digitN = (AliEMCALDigit *) digits->At(fDigitsList[iDigitN]) ; 
607         
608         if ( AreNeighbours(digit, digitN) ) {
609           if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {    
610             maxAt[iDigitN] = 0 ;
611             // but may be digit too is not local max ?
612             if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut) 
613               maxAt[iDigit] = 0 ;
614           }
615           else {
616             maxAt[iDigit] = 0 ;
617             // but may be digitN too is not local max ?
618             if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut) 
619               maxAt[iDigitN] = 0 ; 
620           } 
621         } // if Areneighbours
622       } // while digitN
623     } // slot not empty
624   } // while digit
625   
626   iDigitN = 0 ;
627   for(iDigit = 0; iDigit < fMulDigit; iDigit++) { 
628     if(maxAt[iDigit] ){
629       maxAt[iDigitN] = maxAt[iDigit] ;
630       maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
631       iDigitN++ ; 
632     }
633   }
634   return iDigitN ;
635 }
636 //____________________________________________________________________________
637 void AliEMCALRecPoint::EvalTime(TClonesArray * digits){
638   // time is set to the time of the digit with the maximum energy
639
640   Float_t maxE = 0;
641   Int_t maxAt = 0;
642   for(Int_t idig=0; idig < fMulDigit; idig++){
643     if(fEnergyList[idig] > maxE){
644       maxE = fEnergyList[idig] ;
645       maxAt = idig;
646     }
647   }
648   fTime = ((AliEMCALDigit*) digits->At(fDigitsList[maxAt]))->GetTime() ;
649   
650 }
651
652 //______________________________________________________________________________
653 void AliEMCALRecPoint::Paint(Option_t *)
654 {
655   // Paint this ALiRecPoint as a TMarker  with its current attributes
656   
657   TVector3 pos(0.,0.,0.)  ;
658   GetLocalPosition(pos)   ;
659   Coord_t x = pos.X()     ;
660   Coord_t y = pos.Z()     ;
661   Color_t markercolor = 1 ;
662   Size_t  markersize = 1. ;
663   Style_t markerstyle = 5 ;
664   
665   if (!gPad->IsBatch()) {
666     gVirtualX->SetMarkerColor(markercolor) ;
667     gVirtualX->SetMarkerSize (markersize)  ;
668     gVirtualX->SetMarkerStyle(markerstyle) ;
669   }
670   gPad->SetAttMarkerPS(markercolor,markerstyle,markersize) ;
671   gPad->PaintPolyMarker(1,&x,&y,"") ;
672 }
673
674 //______________________________________________________________________________
675 Float_t AliEMCALRecPoint::EtaToTheta(Float_t arg) const
676 {
677   //Converts Theta (Radians) to Eta(Radians)
678   return (2.*TMath::ATan(TMath::Exp(-arg)));
679 }
680
681 //______________________________________________________________________________
682 Float_t AliEMCALRecPoint::ThetaToEta(Float_t arg) const
683 {
684   //Converts Eta (Radians) to Theta(Radians)
685   return (-1 * TMath::Log(TMath::Tan(0.5 * arg)));
686 }