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