]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSEmcRecPoint.cxx
Threshold for digits in RecPoint introduced
[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 // --- AliRoot header files ---
34
35  #include "AliGenerator.h"
36 #include "AliPHOSGeometry.h"
37 #include "AliPHOSEmcRecPoint.h"
38 #include "AliRun.h"
39 #include "AliPHOSGetter.h"
40
41 ClassImp(AliPHOSEmcRecPoint)
42
43 //____________________________________________________________________________
44 AliPHOSEmcRecPoint::AliPHOSEmcRecPoint() : AliPHOSRecPoint()
45 {
46   // ctor
47
48   fMulDigit   = 0 ;  
49   fAmp   = 0. ;   
50   fCoreEnergy = 0 ; 
51   fEnergyList = 0 ;
52   fTime = -1. ;
53   fLocPos.SetX(1000000.)  ;      //Local position should be evaluated
54    
55 }
56
57 //____________________________________________________________________________
58 AliPHOSEmcRecPoint::AliPHOSEmcRecPoint(const char * opt) : AliPHOSRecPoint(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 AliPHOSEmcRecPoint::~AliPHOSEmcRecPoint()
73 {
74   // dtor
75
76   if ( fEnergyList )
77     delete[] fEnergyList ; 
78 }
79
80 //____________________________________________________________________________
81 void AliPHOSEmcRecPoint::AddDigit(AliPHOSDigit & 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   EvalPHOSMod(&digit) ;
121 }
122
123 //____________________________________________________________________________
124 Bool_t AliPHOSEmcRecPoint::AreNeighbours(AliPHOSDigit * digit1, AliPHOSDigit * digit2 ) const
125 {
126   // Tells if (true) or not (false) two digits are neighbors
127   
128   Bool_t aren = kFALSE ;
129   
130   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
131   AliPHOSGeometry * phosgeom =  (AliPHOSGeometry*)gime->PHOSGeometry();
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 AliPHOSEmcRecPoint::Compare(const TObject * obj) const
150 {
151   // Compares two RecPoints according to their position in the PHOS 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                       //AliPHOSTrackSegmentMakerv* and other RecPoints...
156   Int_t rv ; 
157
158   AliPHOSEmcRecPoint * clu = (AliPHOSEmcRecPoint *)obj ; 
159
160  
161   Int_t phosmod1 = GetPHOSMod() ;
162   Int_t phosmod2 = clu->GetPHOSMod() ;
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 AliPHOSEmcRecPoint::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 AliPHOSRecPoint is clicked with the locator
196   //
197   //  If Left button is clicked on AliPHOSRecPoint, the digits are switched on    
198   //  and switched off when the mouse button is released.
199   
200     
201   AliPHOSGetter * gime =  AliPHOSGetter::GetInstance() ; 
202   if(!gime) return ;
203   AliPHOSGeometry * phosgeom =  (AliPHOSGeometry*)gime->PHOSGeometry();
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     AliPHOSDigit * digit ;
218     Int_t iDigit;
219     Int_t relid[4] ;
220     
221     const Int_t kMulDigit = AliPHOSEmcRecPoint::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 = (AliPHOSDigit *) digits->At(fDigitsList[iDigit])  ;
234       phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
235       phosgeom->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 += phosgeom->GetCrystalSize(0) / 2. ;
246     zimax += phosgeom->GetCrystalSize(2) / 2. ;
247     ximin -= phosgeom->GetCrystalSize(0) / 2. ;
248     zimin -= phosgeom->GetCrystalSize(2) / 2. ;
249     Int_t xdim = (int)( (ximax - ximin ) / phosgeom->GetCrystalSize(0) + 0.5  ) ; 
250     Int_t zdim = (int)( (zimax - zimin ) / phosgeom->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 = (AliPHOSDigit *) digits->At(fDigitsList[iDigit])  ;
266       phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
267       phosgeom->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  AliPHOSEmcRecPoint::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   Float_t x = 0.;
309   Float_t z = 0.;
310
311   AliPHOSDigit * digit ;
312  
313   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
314   AliPHOSGeometry * phosgeom =  (AliPHOSGeometry*)gime->PHOSGeometry();
315   
316
317  // Calculates the center of gravity in the local PHOS-module coordinates 
318   
319   Int_t iDigit;
320
321   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
322     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
323     Int_t relid[4] ;
324     Float_t xi ;
325     Float_t zi ;
326     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
327     phosgeom->RelPosInModule(relid, xi, zi);
328     Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
329     x    += xi * w ;
330     z    += zi * w ;
331     wtot += w ;
332   }
333   x /= wtot ;
334   z /= wtot ;
335
336
337 // Calculates the dispersion in coordinates 
338   wtot = 0.;
339   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
340     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit])  ;
341     Int_t relid[4] ;
342     Float_t xi ;
343     Float_t zi ;
344     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
345     phosgeom->RelPosInModule(relid, xi, zi);
346     Float_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
347     d += w*((xi-x)*(xi-x) + (zi-z)*(zi-z) ) ; 
348     wtot+=w ;
349
350   
351   }
352   
353
354   d /= wtot ;
355
356   fDispersion = TMath::Sqrt(d) ;
357  
358 }
359 //______________________________________________________________________________
360 void AliPHOSEmcRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits)
361 {
362   // This function calculates energy in the core, 
363   // i.e. within a radius rad = 3cm around the center. Beyond this radius
364   // in accordance with shower profile the energy deposition 
365   // should be less than 2%
366
367   Float_t coreRadius = 3 ;
368
369   Float_t x = 0 ;
370   Float_t z = 0 ;
371
372   AliPHOSDigit * digit ;
373
374   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
375   AliPHOSGeometry * phosgeom =  (AliPHOSGeometry*)gime->PHOSGeometry();
376     
377   Int_t iDigit;
378
379 // Calculates the center of gravity in the local PHOS-module coordinates 
380   Float_t wtot = 0;
381   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
382     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
383     Int_t relid[4] ;
384     Float_t xi ;
385     Float_t zi ;
386     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
387     phosgeom->RelPosInModule(relid, xi, zi);
388     Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
389     x    += xi * w ;
390     z    += zi * w ;
391     wtot += w ;
392   }
393   x /= wtot ;
394   z /= wtot ;
395
396
397   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
398     digit = (AliPHOSDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
399     Int_t relid[4] ;
400     Float_t xi ;
401     Float_t zi ;
402     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
403     phosgeom->RelPosInModule(relid, xi, zi);    
404     Float_t distance = TMath::Sqrt((xi-x)*(xi-x)+(zi-z)*(zi-z)) ;
405     if(distance < coreRadius)
406       fCoreEnergy += fEnergyList[iDigit] ;
407   }
408
409 }
410
411 //____________________________________________________________________________
412 void  AliPHOSEmcRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
413 {
414   // Calculates the axis of the shower ellipsoid
415
416   Double_t wtot = 0. ;
417   Double_t x    = 0.;
418   Double_t z    = 0.;
419   Double_t dxx  = 0.;
420   Double_t dzz  = 0.;
421   Double_t dxz  = 0.;
422
423   AliPHOSDigit * digit ;
424
425   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
426   AliPHOSGeometry * phosgeom =  (AliPHOSGeometry*)gime->PHOSGeometry();
427
428   Int_t iDigit;
429
430
431   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
432     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit])  ;
433     Int_t relid[4] ;
434     Float_t xi ;
435     Float_t zi ;
436     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
437     phosgeom->RelPosInModule(relid, xi, zi);
438     Double_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
439     dxx  += w * xi * xi ;
440     x    += w * xi ;
441     dzz  += w * zi * zi ;
442     z    += w * zi ; 
443     dxz  += w * xi * zi ; 
444     wtot += w ;
445   }
446   dxx /= wtot ;
447   x   /= wtot ;
448   dxx -= x * x ;
449   dzz /= wtot ;
450   z   /= wtot ;
451   dzz -= z * z ;
452   dxz /= wtot ;
453   dxz -= x * z ;
454
455 //   //Apply correction due to non-perpendicular incidence
456 //   Double_t CosX ;
457 //   Double_t CosZ ;
458 //   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
459 //   AliPHOSGeometry * phosgeom =  (AliPHOSGeometry*)gime->PHOSGeometry();
460   //   Double_t DistanceToIP= (Double_t ) phosgeom->GetIPtoCrystalSurface() ;
461   
462 //   CosX = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+x*x) ;
463 //   CosZ = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+z*z) ;
464
465 //   dxx = dxx/(CosX*CosX) ;
466 //   dzz = dzz/(CosZ*CosZ) ;
467 //   dxz = dxz/(CosX*CosZ) ;
468
469
470   fLambda[0] =  0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
471   if(fLambda[0] > 0)
472     fLambda[0] = TMath::Sqrt(fLambda[0]) ;
473
474   fLambda[1] =  0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
475   if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
476     fLambda[1] = TMath::Sqrt(fLambda[1]) ;
477   else
478     fLambda[1]= 0. ;
479 }
480
481 //____________________________________________________________________________
482 void AliPHOSEmcRecPoint::EvalAll(Float_t logWeight, TClonesArray * digits )
483 {
484   // Evaluates all shower parameters
485
486   AliPHOSRecPoint::EvalAll(logWeight,digits) ;
487   EvalLocalPosition(logWeight, digits) ;
488   EvalElipsAxis(logWeight, digits) ;
489   EvalDispersion(logWeight, digits) ;
490   EvalCoreEnergy(logWeight, digits);
491   EvalTime(digits) ;
492 }
493 //____________________________________________________________________________
494 void AliPHOSEmcRecPoint::EvalLocalPosition(Float_t logWeight, TClonesArray * digits)
495 {
496   // Calculates the center of gravity in the local PHOS-module coordinates 
497   Float_t wtot = 0. ;
498  
499   Int_t relid[4] ;
500
501   Float_t x = 0. ;
502   Float_t z = 0. ;
503   
504   AliPHOSDigit * digit ;
505
506   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
507   AliPHOSGeometry * phosgeom =  (AliPHOSGeometry*)gime->PHOSGeometry();
508
509   Int_t iDigit;
510
511   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
512     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
513
514     Float_t xi ;
515     Float_t zi ;
516     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
517     phosgeom->RelPosInModule(relid, xi, zi);
518     Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
519     x    += xi * w ;
520     z    += zi * w ;
521     wtot += w ;
522
523   }
524
525   x /= wtot ;
526   z /= wtot ;
527
528   // Correction for the depth of the shower starting point (TDR p 127)  
529   Float_t para = 0.925 ; 
530   Float_t parb = 6.52 ; 
531
532   Float_t xo,yo,zo ; //Coordinates of the origin
533   gAlice->Generator()->GetOrigin(xo,yo,zo) ;
534
535   Float_t phi = phosgeom->GetPHOSAngle(relid[0]) ;
536
537   //Transform to the local ref.frame
538   Float_t xoL,yoL ;
539   xoL = xo*TMath::Cos(phi)-yo*TMath::Sin(phi) ;
540   yoL = xo*TMath::Sin(phi)+yo*TMath::Cos(phi) ;
541   
542   Float_t radius = phosgeom->GetIPtoCrystalSurface()-yoL;
543  
544   Float_t incidencephi = TMath::ATan((x-xoL ) / radius) ; 
545   Float_t incidencetheta = TMath::ATan((z-zo) / radius) ;
546  
547   Float_t depthx =  ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencephi) ; 
548   Float_t depthz =  ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencetheta) ; 
549
550   fLocPos.SetX(x - depthx)  ;
551   fLocPos.SetY(0.) ;
552   fLocPos.SetZ(z - depthz)  ;
553
554   fLocPosM = 0 ;
555 }
556
557 //____________________________________________________________________________
558 Float_t AliPHOSEmcRecPoint::GetMaximalEnergy(void) const
559 {
560   // Finds the maximum energy in the cluster
561   
562   Float_t menergy = 0. ;
563
564   Int_t iDigit;
565
566   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
567  
568     if(fEnergyList[iDigit] > menergy) 
569       menergy = fEnergyList[iDigit] ;
570   }
571   return menergy ;
572 }
573
574 //____________________________________________________________________________
575 Int_t AliPHOSEmcRecPoint::GetMultiplicityAtLevel(const Float_t H) const
576 {
577   // Calculates the multiplicity of digits with energy larger than H*energy 
578   
579   Int_t multipl   = 0 ;
580   Int_t iDigit ;
581   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
582
583     if(fEnergyList[iDigit] > H * fAmp) 
584       multipl++ ;
585   }
586   return multipl ;
587 }
588
589 //____________________________________________________________________________
590 Int_t  AliPHOSEmcRecPoint::GetNumberOfLocalMax( AliPHOSDigit **  maxAt, Float_t * maxAtEnergy,
591                                                Float_t locMaxCut,TClonesArray * digits) const
592
593   // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
594   // energy difference between two local maxima
595
596   AliPHOSDigit * digit ;
597   AliPHOSDigit * digitN ;
598   
599
600   Int_t iDigitN ;
601   Int_t iDigit ;
602
603   for(iDigit = 0; iDigit < fMulDigit; iDigit++)
604     maxAt[iDigit] = (AliPHOSDigit*) digits->At(fDigitsList[iDigit])  ;
605
606   
607   for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {   
608     if(maxAt[iDigit]) {
609       digit = maxAt[iDigit] ;
610           
611       for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {        
612         digitN = (AliPHOSDigit *) digits->At(fDigitsList[iDigitN]) ; 
613         
614         if ( AreNeighbours(digit, digitN) ) {
615           if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {    
616             maxAt[iDigitN] = 0 ;
617             // but may be digit too is not local max ?
618             if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut) 
619               maxAt[iDigit] = 0 ;
620           }
621           else {
622             maxAt[iDigit] = 0 ;
623             // but may be digitN too is not local max ?
624             if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut) 
625               maxAt[iDigitN] = 0 ; 
626           } 
627         } // if Areneighbours
628       } // while digitN
629     } // slot not empty
630   } // while digit
631   
632   iDigitN = 0 ;
633   for(iDigit = 0; iDigit < fMulDigit; iDigit++) { 
634     if(maxAt[iDigit]){
635       maxAt[iDigitN] = maxAt[iDigit] ;
636       maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
637       iDigitN++ ; 
638     }
639   }
640   return iDigitN ;
641 }
642 //____________________________________________________________________________
643 void AliPHOSEmcRecPoint::EvalTime(TClonesArray * digits){
644   
645   Float_t maxE = 0;
646   Int_t maxAt = 0;
647   for(Int_t idig=0; idig < fMulDigit; idig++){
648     if(fEnergyList[idig] > maxE){
649       maxE = fEnergyList[idig] ;
650       maxAt = idig;
651     }
652   }
653   fTime = ((AliPHOSDigit*) digits->At(fDigitsList[maxAt]))->GetTime() ;
654   
655 }
656 //____________________________________________________________________________
657 void AliPHOSEmcRecPoint::Purify(Float_t threshold){
658   //Removes digits below threshold
659
660   Int_t * tempo = new ( Int_t[fMaxDigit] ) ; 
661   Float_t * tempoE =  new ( Float_t[fMaxDigit] ) ;
662
663   Int_t mult = 0 ;
664   for(Int_t iDigit=0;iDigit< fMulDigit ;iDigit++){
665     if(fEnergyList[iDigit] > threshold){
666       tempo[mult] = fDigitsList[iDigit] ;
667       tempoE[mult] = fEnergyList[iDigit] ;
668       mult++ ;
669     }
670   }
671   
672   fMulDigit = mult ;
673   delete [] fDigitsList ;
674   delete [] fEnergyList ;
675   fDigitsList = new (Int_t[fMulDigit]) ;
676   fEnergyList = new ( Float_t[fMulDigit]) ;
677
678   for(Int_t iDigit=0;iDigit< fMulDigit ;iDigit++){
679       fDigitsList[iDigit] = tempo[iDigit];
680       fEnergyList[iDigit] = tempoE[iDigit] ;
681   }
682       
683   delete [] tempo ;
684   delete [] tempoE ;
685
686 }
687 //____________________________________________________________________________
688 void AliPHOSEmcRecPoint::Print(Option_t * option) const
689 {
690   // Print the list of digits belonging to the cluster
691   
692   TString message ; 
693   message  = "AliPHOSEmcRecPoint:\n" ;
694   message +=  " digits # = " ; 
695   Info("Print", message.Data()) ; 
696
697   Int_t iDigit;
698   for(iDigit=0; iDigit<fMulDigit; iDigit++)
699     Info("Print", " %d ", fDigitsList[iDigit] ) ;  
700   
701   Info("Print", " Energies = ") ;
702   for(iDigit=0; iDigit<fMulDigit; iDigit++) 
703     Info("Print", " %f ", fEnergyList[iDigit] ) ;
704   
705    Info("Print", " Primaries  ") ;
706   for(iDigit = 0;iDigit < fMulTrack; iDigit++)
707     Info("Print", " %d ", fTracksList[iDigit]) ;
708         
709   message  = "       Multiplicity    = %d" ;
710   message += "       Cluster Energy  = %f" ; 
711   message += "       Number of primaries %d" ; 
712   message += "       Stored at position %d" ; 
713  
714   Info("Print", message.Data(), fMulDigit, fAmp, fMulTrack,GetIndexInList() ) ;  
715 }
716  
717