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