]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALTowerRecPoint.cxx
default constructor sets datamembers to zero
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALTowerRecPoint.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
16 /* $Id$ */
17
18 //_________________________________________________________________________
19 //  RecPoint implementation for EMCAL-EMC 
20 //  An TowerRecPoint is a cluster of digits   
21 //*--
22 //*-- Author: Dmitri Peressounko (RRC KI & SUBATECH)
23
24
25 // --- ROOT system ---
26 #include "TPad.h"
27 #include "TH2.h"
28 #include "TMath.h" 
29 #include "TCanvas.h" 
30
31 // --- Standard library ---
32
33 #include <iostream.h> 
34
35 // --- AliRoot header files ---
36
37  #include "AliGenerator.h"
38 #include "AliEMCALGeometry.h"
39 #include "AliEMCALTowerRecPoint.h"
40 #include "AliRun.h"
41 #include "AliEMCALGetter.h"
42
43 ClassImp(AliEMCALTowerRecPoint)
44
45 //____________________________________________________________________________
46 AliEMCALTowerRecPoint::AliEMCALTowerRecPoint() : AliEMCALRecPoint()
47 {
48   // ctor
49
50   fMulDigit   = 0 ;  
51   fAmp   = 0. ;   
52   fCoreEnergy = 0 ; 
53   fEnergyList = 0 ;
54   fTime = 0. ;
55   fLocPos.SetX(0.)  ;      //Local position should be evaluated
56    
57 }
58
59 //____________________________________________________________________________
60 AliEMCALTowerRecPoint::AliEMCALTowerRecPoint(const char * opt) : AliEMCALRecPoint(opt)
61 {
62   // ctor
63   
64   fMulDigit   = 0 ;  
65   fAmp   = 0. ;   
66   fCoreEnergy = 0 ; 
67   fEnergyList = 0 ;
68   fTime = -1. ;
69   fLocPos.SetX(1000000.)  ;      //Local position should be evaluated
70   
71 }
72
73 //____________________________________________________________________________
74 AliEMCALTowerRecPoint::~AliEMCALTowerRecPoint()
75 {
76   // dtor
77
78   if ( fEnergyList )
79     delete[] fEnergyList ; 
80 }
81
82 //____________________________________________________________________________
83 void AliEMCALTowerRecPoint::AddDigit(AliEMCALDigit & digit, Float_t Energy)
84 {
85   // Adds a digit to the RecPoint
86   // and accumulates the total amplitude and the multiplicity 
87   
88   if(fEnergyList == 0)
89     fEnergyList =  new Float_t[fMaxDigit]; 
90
91   if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists 
92     fMaxDigit*=2 ; 
93     Int_t * tempo = new ( Int_t[fMaxDigit] ) ; 
94     Float_t * tempoE =  new ( Float_t[fMaxDigit] ) ;
95
96     Int_t index ;     
97     for ( index = 0 ; index < fMulDigit ; index++ ){
98       tempo[index]  = fDigitsList[index] ;
99       tempoE[index] = fEnergyList[index] ; 
100     }
101     
102     delete [] fDigitsList ; 
103     fDigitsList =  new ( Int_t[fMaxDigit] ) ;
104  
105     delete [] fEnergyList ;
106     fEnergyList =  new ( Float_t[fMaxDigit] ) ;
107
108     for ( index = 0 ; index < fMulDigit ; index++ ){
109       fDigitsList[index] = tempo[index] ;
110       fEnergyList[index] = tempoE[index] ; 
111     }
112  
113     delete [] tempo ;
114     delete [] tempoE ; 
115   } // if
116   
117   fDigitsList[fMulDigit]   = digit.GetIndexInList()  ; 
118   fEnergyList[fMulDigit]   = Energy ;
119   fMulDigit++ ; 
120   fAmp += Energy ; 
121
122   // EvalEMCALMod(&digit) ;
123 }
124
125 //____________________________________________________________________________
126 Bool_t AliEMCALTowerRecPoint::AreNeighbours(AliEMCALDigit * digit1, AliEMCALDigit * digit2 ) const
127 {
128   // Tells if (true) or not (false) two digits are neighbors
129   
130   Bool_t aren = kFALSE ;
131   
132   AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ; 
133   AliEMCALGeometry * phosgeom =  (AliEMCALGeometry*)gime->EMCALGeometry();
134
135   Int_t relid1[4] ; 
136   phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ; 
137
138   Int_t relid2[4] ; 
139   phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ; 
140   
141   Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;  
142   Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;  
143
144   if (( coldiff <= 1 )  && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0)) 
145     aren = kTRUE ;
146   
147   return aren ;
148 }
149
150 //____________________________________________________________________________
151 Int_t AliEMCALTowerRecPoint::Compare(const TObject * obj) const
152 {
153   // Compares two RecPoints according to their position in the EMCAL modules
154
155   Float_t delta = 1 ; //Width of "Sorting row". If you changibg this 
156                       //value (what is senseless) change as vell delta in
157                       //AliEMCALTrackSegmentMakerv* and other RecPoints...
158   Int_t rv ; 
159
160   AliEMCALTowerRecPoint * clu = (AliEMCALTowerRecPoint *)obj ; 
161
162  
163   Int_t phosmod1 = GetEMCALArm() ;
164   Int_t phosmod2 = clu->GetEMCALArm() ;
165
166   TVector3 locpos1; 
167   GetLocalPosition(locpos1) ;
168   TVector3 locpos2;  
169   clu->GetLocalPosition(locpos2) ;  
170
171   if(phosmod1 == phosmod2 ) {
172     Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
173     if (rowdif> 0) 
174       rv = 1 ;
175      else if(rowdif < 0) 
176        rv = -1 ;
177     else if(locpos1.Z()>locpos2.Z()) 
178       rv = -1 ;
179     else 
180       rv = 1 ; 
181      }
182
183   else {
184     if(phosmod1 < phosmod2 ) 
185       rv = -1 ;
186     else 
187       rv = 1 ;
188   }
189
190   return rv ; 
191 }
192 //______________________________________________________________________________
193 void AliEMCALTowerRecPoint::ExecuteEvent(Int_t event, Int_t px, Int_t py) const
194 {
195   
196   // Execute action corresponding to one event
197   //  This member function is called when a AliEMCALRecPoint is clicked with the locator
198   //
199   //  If Left button is clicked on AliEMCALRecPoint, the digits are switched on    
200   //  and switched off when the mouse button is released.
201   
202     
203   // AliEMCALGetter * gime =  AliEMCALGetter::GetInstance() ; 
204 //   if(!gime) return ;
205 //   AliEMCALGeometry * emcalgeom =  (AliEMCALGeometry*)gime->EMCALGeometry();
206   
207 //   static TGraph *  digitgraph = 0 ;
208   
209 //   if (!gPad->IsEditable()) return;
210   
211 //   TH2F * histo = 0 ;
212 //   TCanvas * histocanvas ; 
213
214 //   const TClonesArray * digits = gime->Digits() ;
215   
216 //   switch (event) {
217     
218 //   case kButton1Down: {
219 //     AliEMCALDigit * digit ;
220 //     Int_t iDigit;
221 //     Int_t relid[4] ;
222     
223 //     const Int_t kMulDigit = AliEMCALTowerRecPoint::GetDigitsMultiplicity() ; 
224 //     Float_t * xi = new Float_t[kMulDigit] ; 
225 //     Float_t * zi = new Float_t[kMulDigit] ; 
226     
227 //     // create the histogram for the single cluster 
228 //     // 1. gets histogram boundaries
229 //     Float_t ximax = -999. ; 
230 //     Float_t zimax = -999. ; 
231 //     Float_t ximin = 999. ; 
232 //     Float_t zimin = 999. ;
233     
234 //     for(iDigit=0; iDigit<kMulDigit; iDigit++) {
235 //       digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit])  ;
236 //       emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
237 //       emcalgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
238 //       if ( xi[iDigit] > ximax )
239 //      ximax = xi[iDigit] ; 
240 //       if ( xi[iDigit] < ximin )
241 //      ximin = xi[iDigit] ; 
242 //       if ( zi[iDigit] > zimax )
243 //      zimax = zi[iDigit] ; 
244 //       if ( zi[iDigit] < zimin )
245 //      zimin = zi[iDigit] ;     
246 //     }
247 //     ximax += emcalgeom->GetCrystalSize(0) / 2. ;
248 //     zimax += emcalgeom->GetCrystalSize(2) / 2. ;
249 //     ximin -= emcalgeom->GetCrystalSize(0) / 2. ;
250 //     zimin -= emcalgeom->GetCrystalSize(2) / 2. ;
251 //     Int_t xdim = (int)( (ximax - ximin ) / emcalgeom->GetCrystalSize(0) + 0.5  ) ; 
252 //     Int_t zdim = (int)( (zimax - zimin ) / emcalgeom->GetCrystalSize(2) + 0.5 ) ;
253     
254 //     // 2. gets the histogram title
255     
256 //     Text_t title[100] ; 
257 //     sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
258     
259 //     if (!histo) {
260 //       delete histo ; 
261 //       histo = 0 ; 
262 //     }
263 //     histo = new TH2F("cluster3D", title,  xdim, ximin, ximax, zdim, zimin, zimax)  ;
264     
265 //     Float_t x, z ; 
266 //     for(iDigit=0; iDigit<kMulDigit; iDigit++) {
267 //       digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit])  ;
268 //       emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
269 //       emcalgeom->RelPosInModule(relid, x, z);
270 //       histo->Fill(x, z, fEnergyList[iDigit] ) ;
271 //     }
272     
273 //     if (!digitgraph) {
274 //       digitgraph = new TGraph(kMulDigit,xi,zi);
275 //       digitgraph-> SetMarkerStyle(5) ; 
276 //       digitgraph-> SetMarkerSize(1.) ;
277 //       digitgraph-> SetMarkerColor(1) ;
278 //       digitgraph-> Paint("P") ;
279 //     }
280     
281 //     //    Print() ;
282 //     histocanvas = new TCanvas("cluster", "a single cluster", 600, 500) ; 
283 //     histocanvas->Draw() ; 
284 //     histo->Draw("lego1") ; 
285     
286 //     delete[] xi ; 
287 //     delete[] zi ; 
288     
289 //     break;
290 //   }
291   
292 //   case kButton1Up: 
293 //     if (digitgraph) {
294 //       delete digitgraph  ;
295 //       digitgraph = 0 ;
296 //     }
297 //     break;
298   
299 //    }
300 }
301
302 //____________________________________________________________________________
303 void  AliEMCALTowerRecPoint::EvalDispersion(Float_t logWeight,TClonesArray * digits)
304 {
305   // Calculates the dispersion of the shower at the origine of the RecPoint
306
307   Float_t d    = 0. ;
308   Float_t wtot = 0. ;
309
310   AliEMCALDigit * digit ;
311  
312   AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ; 
313   AliEMCALGeometry * emcalgeom =  (AliEMCALGeometry*)gime->EMCALGeometry();
314   
315
316  // Calculates the center of gravity in the local EMCAL-module coordinates 
317   
318   Int_t iDigit;
319   Int_t relid[4] ;
320
321   if (!fTheta || !fPhi ) {
322     for(iDigit=0; iDigit<fMulDigit; iDigit++) {
323       digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
324       
325       Float_t thetai ;
326       Float_t phii ;
327       emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
328       emcalgeom->PosInAlice(relid, thetai, phii);
329       Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
330       fTheta = fTheta + thetai * w ;
331       fPhi   += (phii * w );
332       wtot  += w ;
333     }
334     
335     if (wtot > 0 ) {
336       fTheta /= wtot ;
337       fPhi   /= wtot ;
338     } else {
339       fTheta = -1. ;
340       fPhi   = -1. ; 
341     }
342         
343   }
344
345   const Float_t kDeg2Rad = TMath::Pi() / static_cast<Double_t>(180) ; 
346   
347   Float_t cyl_radius = emcalgeom->GetIPDistance()+emcalgeom->GetAirGap() ;
348   Float_t x =  cyl_radius * TMath::Cos(fPhi * kDeg2Rad ) ;
349   Float_t y =  cyl_radius * TMath::Cos(fPhi * kDeg2Rad ) ; 
350   Float_t z =  cyl_radius * TMath::Tan(fTheta * kDeg2Rad ) ; 
351   
352 // Calculates the dispersion in coordinates 
353   wtot = 0.;
354   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
355     digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit])  ;
356     Float_t thetai = 0. ;
357     Float_t phii = 0.;
358     emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
359     emcalgeom->PosInAlice(relid, thetai, phii);
360     Float_t xi =  cyl_radius * TMath::Cos(phii * kDeg2Rad ) ;
361     Float_t yi =  cyl_radius * TMath::Sin(phii * kDeg2Rad ) ; 
362     Float_t zi =  cyl_radius * TMath::Tan(thetai * kDeg2Rad ) ; 
363
364     Float_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
365     d += w*((xi-x)*(xi-x) + (yi-y)*(yi-y)+ (zi-z)*(zi-z) ) ; 
366     wtot+=w ;
367   }
368   
369   if ( wtot > 0 ) 
370     d /= wtot ;
371   else 
372     d = 0. ; 
373
374   fDispersion = TMath::Sqrt(d) ;
375  
376 }
377 //______________________________________________________________________________
378 void AliEMCALTowerRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits)
379 {
380   // This function calculates energy in the core, 
381   // i.e. within a radius rad = 3cm around the center. Beyond this radius
382   // in accordance with shower profile the energy deposition 
383   // should be less than 2%
384
385   Float_t coreRadius = 10. ;
386
387   AliEMCALDigit * digit ;
388   Int_t relid[4] ;
389   Float_t wtot = 0. ;
390
391   AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ; 
392   AliEMCALGeometry * emcalgeom =  (AliEMCALGeometry*)gime->EMCALGeometry();
393     
394   Int_t iDigit;
395
396   if (!fTheta || !fPhi ) {
397     for(iDigit=0; iDigit<fMulDigit; iDigit++) {
398       digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
399       
400       Float_t thetai ;
401       Float_t phii ;
402       emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
403       emcalgeom->PosInAlice(relid, thetai, phii);
404       Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
405       fTheta = fTheta + thetai * w ;
406       fPhi   += (phii * w );
407       wtot  += w ;
408     }
409     
410     if (wtot > 0 ) { 
411       fTheta /= wtot ;
412       fPhi   /= wtot ;
413     } else { 
414       fTheta = -1 ; 
415       fPhi   = -1 ; 
416     }
417   }
418
419   const Float_t kDeg2Rad = TMath::Pi() / static_cast<Double_t>(180) ; 
420
421   Float_t cyl_radius = emcalgeom->GetIPDistance()+emcalgeom->GetAirGap() ;
422   Float_t x =  cyl_radius * TMath::Cos(fPhi * kDeg2Rad ) ;
423   Float_t y =  cyl_radius * TMath::Cos(fPhi * kDeg2Rad ) ; 
424   Float_t z =  cyl_radius * TMath::Tan(fTheta * kDeg2Rad ) ; 
425
426   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
427     digit = (AliEMCALDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
428     Int_t relid[4] ;
429     Float_t thetai = 0. ;
430     Float_t phii = 0. ;
431     emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
432     emcalgeom->PosInAlice(relid, thetai, phii);
433     
434     Float_t xi =  cyl_radius * TMath::Cos(phii * kDeg2Rad ) ;
435     Float_t yi =  cyl_radius * TMath::Sin(phii * kDeg2Rad ) ; 
436     Float_t zi =  cyl_radius * TMath::Tan(thetai * kDeg2Rad ) ; 
437      
438     Float_t distance = TMath::Sqrt((xi-x)*(xi-x)+(yi-y)*(yi-y)+(zi-z)*(zi-z)) ;
439     if(distance < coreRadius)
440       fCoreEnergy += fEnergyList[iDigit] ;
441   }
442
443 }
444
445 //____________________________________________________________________________
446 void  AliEMCALTowerRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
447 {
448   // Calculates the axis of the shower ellipsoid
449
450   Double_t wtot = 0. ;
451   Double_t x    = 0.;
452   Double_t z    = 0.;
453   Double_t dxx  = 0.;
454   Double_t dzz  = 0.;
455   Double_t dxz  = 0.;
456
457   AliEMCALDigit * digit ;
458
459   AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ; 
460   AliEMCALGeometry * emcalgeom =  (AliEMCALGeometry*)gime->EMCALGeometry();
461
462   Int_t iDigit;
463   const Float_t kDeg2Rad = TMath::Pi() / static_cast<Double_t>(180) ; 
464   
465   Float_t cyl_radius = emcalgeom->GetIPDistance()+emcalgeom->GetAirGap() ;
466
467   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
468     digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit])  ;
469     Int_t relid[4] ;
470     Float_t thetai = 0. ;
471     Float_t phii = 0. ; 
472     emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
473     emcalgeom->PosInAlice(relid, thetai, phii);
474     Double_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
475     Float_t xi =  cyl_radius * TMath::Cos(fPhi * kDeg2Rad ) ;
476     Float_t zi =  cyl_radius * TMath::Tan(fTheta * kDeg2Rad ) ; 
477     dxx  += w * xi * xi ;
478     x    += w * xi ;
479     dzz  += w * zi * zi ;
480     z    += w * zi ; 
481     dxz  += w * xi * zi ; 
482     wtot += w ;
483   }
484   if ( wtot > 0 ) { 
485     dxx /= wtot ;
486     x   /= wtot ;
487     dxx -= x * x ;
488     dzz /= wtot ;
489     z   /= wtot ;
490     dzz -= z * z ;
491     dxz /= wtot ;
492     dxz -= x * z ;
493
494     
495     //   //Apply correction due to non-perpendicular incidence
496 //   Double_t CosX ;
497 //   Double_t CosZ ;
498 //   AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ; 
499 //   AliEMCALGeometry * emcalgeom =  (AliEMCALGeometry*)gime->EMCALGeometry();
500   //   Double_t DistanceToIP= (Double_t ) emcalgeom->GetIPDistance() ;
501   
502 //   CosX = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+x*x) ;
503 //   CosZ = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+z*z) ;
504
505 //   dxx = dxx/(CosX*CosX) ;
506 //   dzz = dzz/(CosZ*CosZ) ;
507 //   dxz = dxz/(CosX*CosZ) ;
508
509
510     fLambda[0] =  0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
511     if(fLambda[0] > 0)
512       fLambda[0] = TMath::Sqrt(fLambda[0]) ;
513     
514     fLambda[1] =  0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
515     if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
516       fLambda[1] = TMath::Sqrt(fLambda[1]) ;
517     else
518       fLambda[1]= 0. ;
519   } else { 
520     fLambda[0]= 0. ;
521     fLambda[1]= 0. ;
522   }
523 }
524
525 //____________________________________________________________________________
526 void AliEMCALTowerRecPoint::EvalAll(Float_t logWeight, TClonesArray * digits )
527 {
528   // Evaluates all shower parameters
529
530   AliEMCALRecPoint::EvalAll(logWeight,digits) ;
531   EvalGlobalPosition(logWeight, digits) ;
532   EvalElipsAxis(logWeight, digits) ;
533   EvalDispersion(logWeight, digits) ;
534   EvalCoreEnergy(logWeight, digits);
535   EvalTime(digits) ;
536 }
537
538 //____________________________________________________________________________
539 void AliEMCALTowerRecPoint::EvalGlobalPosition(Float_t logWeight, TClonesArray * digits)
540 {
541   // Calculates the center of gravity in the local EMCAL-module coordinates 
542   Float_t wtot = 0. ;
543  
544   Int_t relid[4] ;
545   
546   AliEMCALDigit * digit ;
547
548   AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ; 
549   AliEMCALGeometry * emcalgeom  =  static_cast<AliEMCALGeometry*>(gime->EMCALGeometry());
550   Int_t iDigit;
551
552   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
553     digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
554
555     Float_t thetai ;
556     Float_t phii ;
557     emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
558     emcalgeom->PosInAlice(relid, thetai, phii);
559     Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
560     fTheta = fTheta + thetai * w ;
561     fPhi   += (phii * w );
562     wtot  += w ;
563   }
564
565   if ( wtot > 0 ) { 
566     fTheta /= wtot ;
567     fPhi   /= wtot ;
568   } else {
569     fTheta = -1 ; 
570     fPhi   = -1.; 
571   }
572   
573   fLocPos.SetX(0.)  ;
574   fLocPos.SetY(0.) ;
575   fLocPos.SetZ(0.)  ;
576
577   fLocPosM = 0 ;
578 }
579
580 //____________________________________________________________________________
581 Float_t AliEMCALTowerRecPoint::GetMaximalEnergy(void) const
582 {
583   // Finds the maximum energy in the cluster
584   
585   Float_t menergy = 0. ;
586
587   Int_t iDigit;
588
589   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
590  
591     if(fEnergyList[iDigit] > menergy) 
592       menergy = fEnergyList[iDigit] ;
593   }
594   return menergy ;
595 }
596
597 //____________________________________________________________________________
598 Int_t AliEMCALTowerRecPoint::GetMultiplicityAtLevel(const Float_t H) const
599 {
600   // Calculates the multiplicity of digits with energy larger than H*energy 
601   
602   Int_t multipl   = 0 ;
603   Int_t iDigit ;
604   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
605
606     if(fEnergyList[iDigit] > H * fAmp) 
607       multipl++ ;
608   }
609   return multipl ;
610 }
611
612 //____________________________________________________________________________
613 Int_t  AliEMCALTowerRecPoint::GetNumberOfLocalMax(Int_t *  maxAt, Float_t * maxAtEnergy,
614                                                Float_t locMaxCut,TClonesArray * digits) const
615
616   // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
617   // energy difference between two local maxima
618
619   AliEMCALDigit * digit ;
620   AliEMCALDigit * digitN ;
621   
622
623   Int_t iDigitN ;
624   Int_t iDigit ;
625
626   for(iDigit = 0; iDigit < fMulDigit; iDigit++)
627     maxAt[iDigit] = (Int_t) digits->At(fDigitsList[iDigit])  ;
628
629   
630   for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {   
631     if(maxAt[iDigit] != -1) {
632       digit = (AliEMCALDigit *) maxAt[iDigit] ;
633           
634       for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {        
635         digitN = (AliEMCALDigit *) digits->At(fDigitsList[iDigitN]) ; 
636         
637         if ( AreNeighbours(digit, digitN) ) {
638           if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {    
639             maxAt[iDigitN] = -1 ;
640             // but may be digit too is not local max ?
641             if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut) 
642               maxAt[iDigit] = -1 ;
643           }
644           else {
645             maxAt[iDigit] = -1 ;
646             // but may be digitN too is not local max ?
647             if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut) 
648               maxAt[iDigitN] = -1 ; 
649           } 
650         } // if Areneighbours
651       } // while digitN
652     } // slot not empty
653   } // while digit
654   
655   iDigitN = 0 ;
656   for(iDigit = 0; iDigit < fMulDigit; iDigit++) { 
657     if(maxAt[iDigit] != -1){
658       maxAt[iDigitN] = maxAt[iDigit] ;
659       maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
660       iDigitN++ ; 
661     }
662   }
663   return iDigitN ;
664 }
665 //____________________________________________________________________________
666 void AliEMCALTowerRecPoint::EvalTime(TClonesArray * digits){
667   
668   Float_t maxE = 0;
669   Int_t maxAt = 0;
670   for(Int_t idig=0; idig < fMulDigit; idig++){
671     if(fEnergyList[idig] > maxE){
672       maxE = fEnergyList[idig] ;
673       maxAt = idig;
674     }
675   }
676   fTime = ((AliEMCALDigit*) digits->At(fDigitsList[maxAt]))->GetTime() ;
677   
678 }
679 //____________________________________________________________________________
680 void AliEMCALTowerRecPoint::Print(Option_t * option) 
681 {
682   // Print the list of digits belonging to the cluster
683   
684   cout << "AliEMCALTowerRecPoint: " << endl ;
685
686   Int_t iDigit;
687   cout << " digits # = " ;
688   for(iDigit=0; iDigit<fMulDigit; iDigit++)
689     cout << fDigitsList[iDigit] << "  " ;  
690   cout << endl ;
691   
692   cout << " Energies = " ;
693   for(iDigit=0; iDigit<fMulDigit; iDigit++) 
694     cout  << fEnergyList[iDigit] << "  ";
695   cout << endl ;
696   
697   cout << " Primaries  " ;
698   for(iDigit = 0;iDigit < fMulTrack; iDigit++)
699     cout << fTracksList[iDigit] << " " << endl ;
700         
701   cout << "       Multiplicity    = " << fMulDigit  << endl ;
702   cout << "       Cluster Energy  = " << fAmp << endl ;
703   cout << "       Number of primaries " << fMulTrack << endl ;
704   cout << "       Stored at position " << GetIndexInList() << endl ; 
705  
706 }
707