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