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