corrected a bug in dispersion calculations (replaced Sin by Cos)
[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 // --- AliRoot header files ---
34
35  #include "AliGenerator.h"
36 #include "AliEMCALGeometry.h"
37 #include "AliEMCALTowerRecPoint.h"
38 #include "AliRun.h"
39 #include "AliEMCALGetter.h"
40
41 ClassImp(AliEMCALTowerRecPoint)
42
43 //____________________________________________________________________________
44 AliEMCALTowerRecPoint::AliEMCALTowerRecPoint() : AliEMCALRecPoint()
45 {
46   // ctor
47
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    
55 }
56
57 //____________________________________________________________________________
58 AliEMCALTowerRecPoint::AliEMCALTowerRecPoint(const char * opt) : AliEMCALRecPoint(opt)
59 {
60   // ctor
61   
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   
69 }
70
71 //____________________________________________________________________________
72 AliEMCALTowerRecPoint::~AliEMCALTowerRecPoint()
73 {
74   // dtor
75
76   if ( fEnergyList )
77     delete[] fEnergyList ; 
78 }
79
80 //____________________________________________________________________________
81 void AliEMCALTowerRecPoint::AddDigit(AliEMCALDigit & digit, Float_t Energy)
82 {
83   // Adds a digit to the RecPoint
84   // and accumulates the total amplitude and the multiplicity 
85   
86   if(fEnergyList == 0)
87     fEnergyList =  new Float_t[fMaxDigit]; 
88
89   if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists 
90     fMaxDigit*=2 ; 
91     Int_t * tempo = new ( Int_t[fMaxDigit] ) ; 
92     Float_t * tempoE =  new ( Float_t[fMaxDigit] ) ;
93
94     Int_t index ;     
95     for ( index = 0 ; index < fMulDigit ; index++ ){
96       tempo[index]  = fDigitsList[index] ;
97       tempoE[index] = fEnergyList[index] ; 
98     }
99     
100     delete [] fDigitsList ; 
101     fDigitsList =  new ( Int_t[fMaxDigit] ) ;
102  
103     delete [] fEnergyList ;
104     fEnergyList =  new ( Float_t[fMaxDigit] ) ;
105
106     for ( index = 0 ; index < fMulDigit ; index++ ){
107       fDigitsList[index] = tempo[index] ;
108       fEnergyList[index] = tempoE[index] ; 
109     }
110  
111     delete [] tempo ;
112     delete [] tempoE ; 
113   } // if
114   
115   fDigitsList[fMulDigit]   = digit.GetIndexInList()  ; 
116   fEnergyList[fMulDigit]   = Energy ;
117   fMulDigit++ ; 
118   fAmp += Energy ; 
119
120   // EvalEMCALMod(&digit) ;
121 }
122
123 //____________________________________________________________________________
124 Bool_t AliEMCALTowerRecPoint::AreNeighbours(AliEMCALDigit * digit1, AliEMCALDigit * digit2 ) const
125 {
126   // Tells if (true) or not (false) two digits are neighbors
127   
128   Bool_t aren = kFALSE ;
129   
130   AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ; 
131   AliEMCALGeometry * phosgeom =  (AliEMCALGeometry*)gime->EMCALGeometry();
132
133   Int_t relid1[4] ; 
134   phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ; 
135
136   Int_t relid2[4] ; 
137   phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ; 
138   
139   Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;  
140   Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;  
141
142   if (( coldiff <= 1 )  && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0)) 
143     aren = kTRUE ;
144   
145   return aren ;
146 }
147
148 //____________________________________________________________________________
149 Int_t AliEMCALTowerRecPoint::Compare(const TObject * obj) const
150 {
151   // Compares two RecPoints according to their position in the EMCAL modules
152
153   Float_t delta = 1 ; //Width of "Sorting row". If you changibg this 
154                       //value (what is senseless) change as vell delta in
155                       //AliEMCALTrackSegmentMakerv* and other RecPoints...
156   Int_t rv ; 
157
158   AliEMCALTowerRecPoint * clu = (AliEMCALTowerRecPoint *)obj ; 
159
160  
161   Int_t phosmod1 = GetEMCALArm() ;
162   Int_t phosmod2 = clu->GetEMCALArm() ;
163
164   TVector3 locpos1; 
165   GetLocalPosition(locpos1) ;
166   TVector3 locpos2;  
167   clu->GetLocalPosition(locpos2) ;  
168
169   if(phosmod1 == phosmod2 ) {
170     Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
171     if (rowdif> 0) 
172       rv = 1 ;
173      else if(rowdif < 0) 
174        rv = -1 ;
175     else if(locpos1.Z()>locpos2.Z()) 
176       rv = -1 ;
177     else 
178       rv = 1 ; 
179      }
180
181   else {
182     if(phosmod1 < phosmod2 ) 
183       rv = -1 ;
184     else 
185       rv = 1 ;
186   }
187
188   return rv ; 
189 }
190 //______________________________________________________________________________
191 void AliEMCALTowerRecPoint::ExecuteEvent(Int_t event, Int_t px, Int_t py) const
192 {
193   
194   // Execute action corresponding to one event
195   //  This member function is called when a AliEMCALRecPoint is clicked with the locator
196   //
197   //  If Left button is clicked on AliEMCALRecPoint, the digits are switched on    
198   //  and switched off when the mouse button is released.
199   
200     
201   // AliEMCALGetter * gime =  AliEMCALGetter::GetInstance() ; 
202 //   if(!gime) return ;
203 //   AliEMCALGeometry * emcalgeom =  (AliEMCALGeometry*)gime->EMCALGeometry();
204   
205 //   static TGraph *  digitgraph = 0 ;
206   
207 //   if (!gPad->IsEditable()) return;
208   
209 //   TH2F * histo = 0 ;
210 //   TCanvas * histocanvas ; 
211
212 //   const TClonesArray * digits = gime->Digits() ;
213   
214 //   switch (event) {
215     
216 //   case kButton1Down: {
217 //     AliEMCALDigit * digit ;
218 //     Int_t iDigit;
219 //     Int_t relid[4] ;
220     
221 //     const Int_t kMulDigit = AliEMCALTowerRecPoint::GetDigitsMultiplicity() ; 
222 //     Float_t * xi = new Float_t[kMulDigit] ; 
223 //     Float_t * zi = new Float_t[kMulDigit] ; 
224     
225 //     // create the histogram for the single cluster 
226 //     // 1. gets histogram boundaries
227 //     Float_t ximax = -999. ; 
228 //     Float_t zimax = -999. ; 
229 //     Float_t ximin = 999. ; 
230 //     Float_t zimin = 999. ;
231     
232 //     for(iDigit=0; iDigit<kMulDigit; iDigit++) {
233 //       digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit])  ;
234 //       emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
235 //       emcalgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
236 //       if ( xi[iDigit] > ximax )
237 //      ximax = xi[iDigit] ; 
238 //       if ( xi[iDigit] < ximin )
239 //      ximin = xi[iDigit] ; 
240 //       if ( zi[iDigit] > zimax )
241 //      zimax = zi[iDigit] ; 
242 //       if ( zi[iDigit] < zimin )
243 //      zimin = zi[iDigit] ;     
244 //     }
245 //     ximax += emcalgeom->GetCrystalSize(0) / 2. ;
246 //     zimax += emcalgeom->GetCrystalSize(2) / 2. ;
247 //     ximin -= emcalgeom->GetCrystalSize(0) / 2. ;
248 //     zimin -= emcalgeom->GetCrystalSize(2) / 2. ;
249 //     Int_t xdim = (int)( (ximax - ximin ) / emcalgeom->GetCrystalSize(0) + 0.5  ) ; 
250 //     Int_t zdim = (int)( (zimax - zimin ) / emcalgeom->GetCrystalSize(2) + 0.5 ) ;
251     
252 //     // 2. gets the histogram title
253     
254 //     Text_t title[100] ; 
255 //     sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
256     
257 //     if (!histo) {
258 //       delete histo ; 
259 //       histo = 0 ; 
260 //     }
261 //     histo = new TH2F("cluster3D", title,  xdim, ximin, ximax, zdim, zimin, zimax)  ;
262     
263 //     Float_t x, z ; 
264 //     for(iDigit=0; iDigit<kMulDigit; iDigit++) {
265 //       digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit])  ;
266 //       emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
267 //       emcalgeom->RelPosInModule(relid, x, z);
268 //       histo->Fill(x, z, fEnergyList[iDigit] ) ;
269 //     }
270     
271 //     if (!digitgraph) {
272 //       digitgraph = new TGraph(kMulDigit,xi,zi);
273 //       digitgraph-> SetMarkerStyle(5) ; 
274 //       digitgraph-> SetMarkerSize(1.) ;
275 //       digitgraph-> SetMarkerColor(1) ;
276 //       digitgraph-> Paint("P") ;
277 //     }
278     
279 //     //    Print() ;
280 //     histocanvas = new TCanvas("cluster", "a single cluster", 600, 500) ; 
281 //     histocanvas->Draw() ; 
282 //     histo->Draw("lego1") ; 
283     
284 //     delete[] xi ; 
285 //     delete[] zi ; 
286     
287 //     break;
288 //   }
289   
290 //   case kButton1Up: 
291 //     if (digitgraph) {
292 //       delete digitgraph  ;
293 //       digitgraph = 0 ;
294 //     }
295 //     break;
296   
297 //    }
298 }
299
300 //____________________________________________________________________________
301 void  AliEMCALTowerRecPoint::EvalDispersion(Float_t logWeight,TClonesArray * digits)
302 {
303   // Calculates the dispersion of the shower at the origine of the RecPoint
304   printf("**************** EVAL Dispersion *****************") ; 
305
306   Float_t d    = 0. ;
307   Float_t wtot = 0. ;
308
309   AliEMCALDigit * digit ;
310  
311   AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ; 
312   AliEMCALGeometry * emcalgeom =  (AliEMCALGeometry*)gime->EMCALGeometry();
313   
314
315  // Calculates the center of gravity in the local EMCAL-module coordinates 
316   
317   Int_t iDigit;
318   Int_t relid[4] ;
319
320   if (!fTheta || !fPhi ) {
321     for(iDigit=0; iDigit<fMulDigit; iDigit++) {
322       digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
323       
324       Float_t thetai ;
325       Float_t phii ;
326       emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
327       emcalgeom->PosInAlice(relid, thetai, phii);
328       Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
329       fTheta = fTheta + thetai * w ;
330       fPhi   += (phii * w );
331       wtot  += w ;
332     }
333     
334     if (wtot > 0 ) {
335       fTheta /= wtot ;
336       fPhi   /= wtot ;
337     } else {
338       fTheta = -1. ;
339       fPhi   = -1. ; 
340     }
341         
342   }
343
344   const Float_t kDeg2Rad = TMath::Pi() / static_cast<Double_t>(180) ; 
345   
346   Float_t cyl_radius = emcalgeom->GetIPDistance()+emcalgeom->GetAirGap() ;
347   Float_t x =  cyl_radius * TMath::Cos(fPhi * kDeg2Rad ) ;
348   Float_t y =  cyl_radius * TMath::Sin(fPhi * kDeg2Rad ) ; 
349   Float_t z =  cyl_radius * TMath::Tan(fTheta * kDeg2Rad ) ; 
350   
351 // Calculates the dispersion in coordinates 
352   wtot = 0.;
353   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
354     digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit])  ;
355     Float_t thetai = 0. ;
356     Float_t phii = 0.;
357     emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
358     emcalgeom->PosInAlice(relid, thetai, phii);
359     Float_t xi =  cyl_radius * TMath::Cos(phii * kDeg2Rad ) ;
360     Float_t yi =  cyl_radius * TMath::Sin(phii * kDeg2Rad ) ; 
361     Float_t zi =  cyl_radius * TMath::Tan(thetai * kDeg2Rad ) ; 
362
363     Float_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
364     d += w*((xi-x)*(xi-x) + (yi-y)*(yi-y)+ (zi-z)*(zi-z) ) ; 
365     wtot+=w ;
366     printf("xi=%f, x=%f\n yi=%f, y=%f\n zi=%f, z=%f\n phi=%f phii=%f theta=%f thetaii=%f\n\n",xi,x,yi,y,zi,z, fPhi, phii, fTheta, thetai) ;  
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(AliEMCALDigit **  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] = (AliEMCALDigit*) digits->At(fDigitsList[iDigit])  ;
628
629   
630   for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {   
631     if(maxAt[iDigit]) {
632       digit = 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] = 0 ;
640             // but may be digit too is not local max ?
641             if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut) 
642               maxAt[iDigit] = 0 ;
643           }
644           else {
645             maxAt[iDigit] = 0 ;
646             // but may be digitN too is not local max ?
647             if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut) 
648               maxAt[iDigitN] = 0 ; 
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] ){
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   TString message("\n") ; 
685
686   Int_t iDigit;
687   message += "digits # = " ;
688   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
689     message += fDigitsList[iDigit] ; 
690     message += "  " ;
691   } 
692   
693   message += "\nEnergies = " ;
694   for(iDigit=0; iDigit<fMulDigit; iDigit++) { 
695     message += fEnergyList[iDigit] ; 
696     message += "  " ;
697   }
698   
699    message += "\nPrimaries  " ;
700    for(iDigit = 0;iDigit < fMulTrack; iDigit++) {
701      message += fTracksList[iDigit] ;
702      message += " " ;
703    }
704    message += "\n       Multiplicity    = " ; 
705    message += fMulDigit ;
706    message += "\n       Cluster Energy  = " ; 
707    message += fAmp ;
708    message += "\n       Number of primaries " ; 
709    message += fMulTrack ;
710    message += "\n       Stored at position " ;
711    message += GetIndexInList() ; 
712    
713    Info("Print", message.Data() ) ; 
714 }
715