]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALTowerRecPoint.cxx
switch on special cuts process
[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
305   Float_t d    = 0. ;
306   Float_t wtot = 0. ;
307
308   AliEMCALDigit * digit ;
309  
310   AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ; 
311   AliEMCALGeometry * emcalgeom =  (AliEMCALGeometry*)gime->EMCALGeometry();
312   
313
314  // Calculates the center of gravity in the local EMCAL-module coordinates 
315   
316   Int_t iDigit;
317   Int_t relid[4] ;
318
319   if (!fTheta || !fPhi ) {
320     for(iDigit=0; iDigit<fMulDigit; iDigit++) {
321       digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
322       
323       Float_t thetai ;
324       Float_t phii ;
325       emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
326       emcalgeom->PosInAlice(relid, thetai, phii);
327       Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
328       fTheta = fTheta + thetai * w ;
329       fPhi   += (phii * w );
330       wtot  += w ;
331     }
332     
333     if (wtot > 0 ) {
334       fTheta /= wtot ;
335       fPhi   /= wtot ;
336     } else {
337       fTheta = -1. ;
338       fPhi   = -1. ; 
339     }
340         
341   }
342
343   const Float_t kDeg2Rad = TMath::Pi() / static_cast<Double_t>(180) ; 
344   
345   Float_t cyl_radius = emcalgeom->GetIPDistance()+emcalgeom->GetAirGap() ;
346   Float_t x =  cyl_radius * TMath::Cos(fPhi * kDeg2Rad ) ;
347   Float_t y =  cyl_radius * TMath::Cos(fPhi * kDeg2Rad ) ; 
348   Float_t z =  cyl_radius * TMath::Tan(fTheta * kDeg2Rad ) ; 
349   
350 // Calculates the dispersion in coordinates 
351   wtot = 0.;
352   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
353     digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit])  ;
354     Float_t thetai = 0. ;
355     Float_t phii = 0.;
356     emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
357     emcalgeom->PosInAlice(relid, thetai, phii);
358     Float_t xi =  cyl_radius * TMath::Cos(phii * kDeg2Rad ) ;
359     Float_t yi =  cyl_radius * TMath::Sin(phii * kDeg2Rad ) ; 
360     Float_t zi =  cyl_radius * TMath::Tan(thetai * kDeg2Rad ) ; 
361
362     Float_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
363     d += w*((xi-x)*(xi-x) + (yi-y)*(yi-y)+ (zi-z)*(zi-z) ) ; 
364     wtot+=w ;
365   }
366   
367   if ( wtot > 0 ) 
368     d /= wtot ;
369   else 
370     d = 0. ; 
371
372   fDispersion = TMath::Sqrt(d) ;
373  
374 }
375 //______________________________________________________________________________
376 void AliEMCALTowerRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits)
377 {
378   // This function calculates energy in the core, 
379   // i.e. within a radius rad = 3cm around the center. Beyond this radius
380   // in accordance with shower profile the energy deposition 
381   // should be less than 2%
382
383   Float_t coreRadius = 10. ;
384
385   AliEMCALDigit * digit ;
386   Int_t relid[4] ;
387   Float_t wtot = 0. ;
388
389   AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ; 
390   AliEMCALGeometry * emcalgeom =  (AliEMCALGeometry*)gime->EMCALGeometry();
391     
392   Int_t iDigit;
393
394   if (!fTheta || !fPhi ) {
395     for(iDigit=0; iDigit<fMulDigit; iDigit++) {
396       digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
397       
398       Float_t thetai ;
399       Float_t phii ;
400       emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
401       emcalgeom->PosInAlice(relid, thetai, phii);
402       Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
403       fTheta = fTheta + thetai * w ;
404       fPhi   += (phii * w );
405       wtot  += w ;
406     }
407     
408     if (wtot > 0 ) { 
409       fTheta /= wtot ;
410       fPhi   /= wtot ;
411     } else { 
412       fTheta = -1 ; 
413       fPhi   = -1 ; 
414     }
415   }
416
417   const Float_t kDeg2Rad = TMath::Pi() / static_cast<Double_t>(180) ; 
418
419   Float_t cyl_radius = emcalgeom->GetIPDistance()+emcalgeom->GetAirGap() ;
420   Float_t x =  cyl_radius * TMath::Cos(fPhi * kDeg2Rad ) ;
421   Float_t y =  cyl_radius * TMath::Cos(fPhi * kDeg2Rad ) ; 
422   Float_t z =  cyl_radius * TMath::Tan(fTheta * kDeg2Rad ) ; 
423
424   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
425     digit = (AliEMCALDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
426     Int_t relid[4] ;
427     Float_t thetai = 0. ;
428     Float_t phii = 0. ;
429     emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
430     emcalgeom->PosInAlice(relid, thetai, phii);
431     
432     Float_t xi =  cyl_radius * TMath::Cos(phii * kDeg2Rad ) ;
433     Float_t yi =  cyl_radius * TMath::Sin(phii * kDeg2Rad ) ; 
434     Float_t zi =  cyl_radius * TMath::Tan(thetai * kDeg2Rad ) ; 
435      
436     Float_t distance = TMath::Sqrt((xi-x)*(xi-x)+(yi-y)*(yi-y)+(zi-z)*(zi-z)) ;
437     if(distance < coreRadius)
438       fCoreEnergy += fEnergyList[iDigit] ;
439   }
440
441 }
442
443 //____________________________________________________________________________
444 void  AliEMCALTowerRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
445 {
446   // Calculates the axis of the shower ellipsoid
447
448   Double_t wtot = 0. ;
449   Double_t x    = 0.;
450   Double_t z    = 0.;
451   Double_t dxx  = 0.;
452   Double_t dzz  = 0.;
453   Double_t dxz  = 0.;
454
455   AliEMCALDigit * digit ;
456
457   AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ; 
458   AliEMCALGeometry * emcalgeom =  (AliEMCALGeometry*)gime->EMCALGeometry();
459
460   Int_t iDigit;
461   const Float_t kDeg2Rad = TMath::Pi() / static_cast<Double_t>(180) ; 
462   
463   Float_t cyl_radius = emcalgeom->GetIPDistance()+emcalgeom->GetAirGap() ;
464
465   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
466     digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit])  ;
467     Int_t relid[4] ;
468     Float_t thetai = 0. ;
469     Float_t phii = 0. ; 
470     emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
471     emcalgeom->PosInAlice(relid, thetai, phii);
472     Double_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
473     Float_t xi =  cyl_radius * TMath::Cos(fPhi * kDeg2Rad ) ;
474     Float_t zi =  cyl_radius * TMath::Tan(fTheta * kDeg2Rad ) ; 
475     dxx  += w * xi * xi ;
476     x    += w * xi ;
477     dzz  += w * zi * zi ;
478     z    += w * zi ; 
479     dxz  += w * xi * zi ; 
480     wtot += w ;
481   }
482   if ( wtot > 0 ) { 
483     dxx /= wtot ;
484     x   /= wtot ;
485     dxx -= x * x ;
486     dzz /= wtot ;
487     z   /= wtot ;
488     dzz -= z * z ;
489     dxz /= wtot ;
490     dxz -= x * z ;
491
492     
493     //   //Apply correction due to non-perpendicular incidence
494 //   Double_t CosX ;
495 //   Double_t CosZ ;
496 //   AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ; 
497 //   AliEMCALGeometry * emcalgeom =  (AliEMCALGeometry*)gime->EMCALGeometry();
498   //   Double_t DistanceToIP= (Double_t ) emcalgeom->GetIPDistance() ;
499   
500 //   CosX = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+x*x) ;
501 //   CosZ = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+z*z) ;
502
503 //   dxx = dxx/(CosX*CosX) ;
504 //   dzz = dzz/(CosZ*CosZ) ;
505 //   dxz = dxz/(CosX*CosZ) ;
506
507
508     fLambda[0] =  0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
509     if(fLambda[0] > 0)
510       fLambda[0] = TMath::Sqrt(fLambda[0]) ;
511     
512     fLambda[1] =  0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
513     if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
514       fLambda[1] = TMath::Sqrt(fLambda[1]) ;
515     else
516       fLambda[1]= 0. ;
517   } else { 
518     fLambda[0]= 0. ;
519     fLambda[1]= 0. ;
520   }
521 }
522
523 //____________________________________________________________________________
524 void AliEMCALTowerRecPoint::EvalAll(Float_t logWeight, TClonesArray * digits )
525 {
526   // Evaluates all shower parameters
527
528   AliEMCALRecPoint::EvalAll(logWeight,digits) ;
529   EvalGlobalPosition(logWeight, digits) ;
530   EvalElipsAxis(logWeight, digits) ;
531   EvalDispersion(logWeight, digits) ;
532   EvalCoreEnergy(logWeight, digits);
533   EvalTime(digits) ;
534 }
535
536 //____________________________________________________________________________
537 void AliEMCALTowerRecPoint::EvalGlobalPosition(Float_t logWeight, TClonesArray * digits)
538 {
539   // Calculates the center of gravity in the local EMCAL-module coordinates 
540   Float_t wtot = 0. ;
541  
542   Int_t relid[4] ;
543   
544   AliEMCALDigit * digit ;
545
546   AliEMCALGetter * gime = AliEMCALGetter::GetInstance() ; 
547   AliEMCALGeometry * emcalgeom  =  static_cast<AliEMCALGeometry*>(gime->EMCALGeometry());
548   Int_t iDigit;
549
550   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
551     digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
552
553     Float_t thetai ;
554     Float_t phii ;
555     emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
556     emcalgeom->PosInAlice(relid, thetai, phii);
557     Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
558     fTheta = fTheta + thetai * w ;
559     fPhi   += (phii * w );
560     wtot  += w ;
561   }
562
563   if ( wtot > 0 ) { 
564     fTheta /= wtot ;
565     fPhi   /= wtot ;
566   } else {
567     fTheta = -1 ; 
568     fPhi   = -1.; 
569   }
570   
571   fLocPos.SetX(0.)  ;
572   fLocPos.SetY(0.) ;
573   fLocPos.SetZ(0.)  ;
574
575   fLocPosM = 0 ;
576 }
577
578 //____________________________________________________________________________
579 Float_t AliEMCALTowerRecPoint::GetMaximalEnergy(void) const
580 {
581   // Finds the maximum energy in the cluster
582   
583   Float_t menergy = 0. ;
584
585   Int_t iDigit;
586
587   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
588  
589     if(fEnergyList[iDigit] > menergy) 
590       menergy = fEnergyList[iDigit] ;
591   }
592   return menergy ;
593 }
594
595 //____________________________________________________________________________
596 Int_t AliEMCALTowerRecPoint::GetMultiplicityAtLevel(const Float_t H) const
597 {
598   // Calculates the multiplicity of digits with energy larger than H*energy 
599   
600   Int_t multipl   = 0 ;
601   Int_t iDigit ;
602   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
603
604     if(fEnergyList[iDigit] > H * fAmp) 
605       multipl++ ;
606   }
607   return multipl ;
608 }
609
610 //____________________________________________________________________________
611 Int_t  AliEMCALTowerRecPoint::GetNumberOfLocalMax(AliEMCALDigit **  maxAt, Float_t * maxAtEnergy,
612                                                Float_t locMaxCut,TClonesArray * digits) const
613
614   // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
615   // energy difference between two local maxima
616
617   AliEMCALDigit * digit ;
618   AliEMCALDigit * digitN ;
619   
620
621   Int_t iDigitN ;
622   Int_t iDigit ;
623
624   for(iDigit = 0; iDigit < fMulDigit; iDigit++)
625     maxAt[iDigit] = (AliEMCALDigit*) digits->At(fDigitsList[iDigit])  ;
626
627   
628   for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {   
629     if(maxAt[iDigit]) {
630       digit = maxAt[iDigit] ;
631           
632       for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {        
633         digitN = (AliEMCALDigit *) digits->At(fDigitsList[iDigitN]) ; 
634         
635         if ( AreNeighbours(digit, digitN) ) {
636           if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {    
637             maxAt[iDigitN] = 0 ;
638             // but may be digit too is not local max ?
639             if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut) 
640               maxAt[iDigit] = 0 ;
641           }
642           else {
643             maxAt[iDigit] = 0 ;
644             // but may be digitN too is not local max ?
645             if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut) 
646               maxAt[iDigitN] = 0 ; 
647           } 
648         } // if Areneighbours
649       } // while digitN
650     } // slot not empty
651   } // while digit
652   
653   iDigitN = 0 ;
654   for(iDigit = 0; iDigit < fMulDigit; iDigit++) { 
655     if(maxAt[iDigit] ){
656       maxAt[iDigitN] = maxAt[iDigit] ;
657       maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
658       iDigitN++ ; 
659     }
660   }
661   return iDigitN ;
662 }
663 //____________________________________________________________________________
664 void AliEMCALTowerRecPoint::EvalTime(TClonesArray * digits){
665   
666   Float_t maxE = 0;
667   Int_t maxAt = 0;
668   for(Int_t idig=0; idig < fMulDigit; idig++){
669     if(fEnergyList[idig] > maxE){
670       maxE = fEnergyList[idig] ;
671       maxAt = idig;
672     }
673   }
674   fTime = ((AliEMCALDigit*) digits->At(fDigitsList[maxAt]))->GetTime() ;
675   
676 }
677 //____________________________________________________________________________
678 void AliEMCALTowerRecPoint::Print(Option_t * option) 
679 {
680   // Print the list of digits belonging to the cluster
681   
682   TString message("\n") ; 
683
684   Int_t iDigit;
685   message += "digits # = " ;
686   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
687     message += fDigitsList[iDigit] ; 
688     message += "  " ;
689   } 
690   
691   message += "\nEnergies = " ;
692   for(iDigit=0; iDigit<fMulDigit; iDigit++) { 
693     message += fEnergyList[iDigit] ; 
694     message += "  " ;
695   }
696   
697    message += "\nPrimaries  " ;
698    for(iDigit = 0;iDigit < fMulTrack; iDigit++) {
699      message += fTracksList[iDigit] ;
700      message += " " ;
701    }
702    message += "\n       Multiplicity    = " ; 
703    message += fMulDigit ;
704    message += "\n       Cluster Energy  = " ; 
705    message += fAmp ;
706    message += "\n       Number of primaries " ; 
707    message += fMulTrack ;
708    message += "\n       Stored at position " ;
709    message += GetIndexInList() ; 
710    
711    Info("Print", message.Data() ) ; 
712 }
713