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