]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSEmcRecPoint.cxx
Reading Trees branches directly to the TFolders added
[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   fLocPos.SetX(1000000.)  ;      //Local position should be evaluated
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 //______________________________________________________________________________
179 void AliPHOSEmcRecPoint::ExecuteEvent(Int_t event, Int_t px, Int_t py) const
180 {
181 //   Commented by Dmitri Peressounko: there is no possibility to ensure, 
182 //   that AliPHOSGetter keeps the correct information.
183
184 //   // Execute action corresponding to one event
185 //   //  This member function is called when a AliPHOSRecPoint is clicked with the locator
186 //   //
187 //   //  If Left button is clicked on AliPHOSRecPoint, the digits are switched on    
188 //   //  and switched off when the mouse button is released.
189 //   //
190
191 //   //   static Int_t pxold, pyold;
192
193 //   AliPHOSGetter * gime =  AliPHOSGetter::GetInstance() ; 
194   
195 //   static TGraph *  digitgraph = 0 ;
196   
197 //   if (!gPad->IsEditable()) return;
198   
199 //   TH2F * histo = 0 ;
200 //   TCanvas * histocanvas ; 
201   
202 //   switch (event) {
203     
204 //   case kButton1Down: {
205 //     AliPHOSDigit * digit ;
206 //     AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
207 //     AliPHOSGeometry * phosgeom =  (AliPHOSGeometry*)gime->PHOSGeometry();
208 //     Int_t iDigit;
209 //     Int_t relid[4] ;
210     
211 //     const Int_t kMulDigit = AliPHOSEmcRecPoint::GetDigitsMultiplicity() ; 
212 //     Float_t * xi = new Float_t[kMulDigit] ; 
213 //     Float_t * zi = new Float_t[kMulDigit] ; 
214     
215 //     // create the histogram for the single cluster 
216 //     // 1. gets histogram boundaries
217 //     Float_t ximax = -999. ; 
218 //     Float_t zimax = -999. ; 
219 //     Float_t ximin = 999. ; 
220 //     Float_t zimin = 999. ;
221     
222 //     for(iDigit=0; iDigit<kMulDigit; iDigit++) {
223 //       digit = (AliPHOSDigit *) ( gime->Digit(fDigitsList[iDigit]) ) ;
224 //       phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
225 //       phosgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
226 //       if ( xi[iDigit] > ximax )
227 //      ximax = xi[iDigit] ; 
228 //       if ( xi[iDigit] < ximin )
229 //      ximin = xi[iDigit] ; 
230 //       if ( zi[iDigit] > zimax )
231 //      zimax = zi[iDigit] ; 
232 //       if ( zi[iDigit] < zimin )
233 //      zimin = zi[iDigit] ;     
234 //     }
235 //     ximax += phosgeom->GetCrystalSize(0) / 2. ;
236 //     zimax += phosgeom->GetCrystalSize(2) / 2. ;
237 //     ximin -= phosgeom->GetCrystalSize(0) / 2. ;
238 //     zimin -= phosgeom->GetCrystalSize(2) / 2. ;
239 //     Int_t xdim = (int)( (ximax - ximin ) / phosgeom->GetCrystalSize(0) + 0.5  ) ; 
240 //     Int_t zdim = (int)( (zimax - zimin ) / phosgeom->GetCrystalSize(2) + 0.5 ) ;
241     
242 //     // 2. gets the histogram title
243     
244 //     Text_t title[100] ; 
245 //     sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
246     
247 //     if (!histo) {
248 //       delete histo ; 
249 //       histo = 0 ; 
250 //     }
251 //     histo = new TH2F("cluster3D", title,  xdim, ximin, ximax, zdim, zimin, zimax)  ;
252     
253 //     Float_t x, z ; 
254 //     for(iDigit=0; iDigit<kMulDigit; iDigit++) {
255 //       digit = (AliPHOSDigit *) ( gime->Digit(fDigitsList[iDigit]) ) ;
256 //       phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
257 //       phosgeom->RelPosInModule(relid, x, z);
258 //       histo->Fill(x, z, fEnergyList[iDigit] ) ;
259 //     }
260     
261 //     if (!digitgraph) {
262 //       digitgraph = new TGraph(kMulDigit,xi,zi);
263 //       digitgraph-> SetMarkerStyle(5) ; 
264 //       digitgraph-> SetMarkerSize(1.) ;
265 //       digitgraph-> SetMarkerColor(1) ;
266 //       digitgraph-> Paint("P") ;
267 //     }
268     
269 //     Print() ;
270 //     histocanvas = new TCanvas("cluster", "a single cluster", 600, 500) ; 
271 //     histocanvas->Draw() ; 
272 //     histo->Draw("lego1") ; 
273     
274 //     delete[] xi ; 
275 //     delete[] zi ; 
276     
277 //     break;
278 //   }
279   
280 //   case kButton1Up: 
281 //     if (digitgraph) {
282 //       delete digitgraph  ;
283 //       digitgraph = 0 ;
284 //     }
285 //     break;
286   
287 //    }
288 }
289
290 //____________________________________________________________________________
291 void  AliPHOSEmcRecPoint::EvalDispersion(Float_t logWeight,TClonesArray * digits)
292 {
293   // Calculates the dispersion of the shower at the origine of the RecPoint
294
295   Float_t d    = 0. ;
296   Float_t wtot = 0. ;
297
298   Float_t x = 0.;
299   Float_t z = 0.;
300
301   AliPHOSDigit * digit ;
302  
303   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
304   AliPHOSGeometry * phosgeom =  (AliPHOSGeometry*)gime->PHOSGeometry();
305   
306
307  // Calculates the center of gravity in the local PHOS-module coordinates 
308   
309   Int_t iDigit;
310
311   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
312     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
313     Int_t relid[4] ;
314     Float_t xi ;
315     Float_t zi ;
316     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
317     phosgeom->RelPosInModule(relid, xi, zi);
318     Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
319     x    += xi * w ;
320     z    += zi * w ;
321     wtot += w ;
322   }
323   x /= wtot ;
324   z /= wtot ;
325
326
327 // Calculates the dispersion in coordinates 
328   wtot = 0.;
329   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
330     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit])  ;
331     Int_t relid[4] ;
332     Float_t xi ;
333     Float_t zi ;
334     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
335     phosgeom->RelPosInModule(relid, xi, zi);
336     Float_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
337     d += w*((xi-x)*(xi-x) + (zi-z)*(zi-z) ) ; 
338     wtot+=w ;
339
340   
341   }
342   
343
344   d /= wtot ;
345
346   fDispersion = TMath::Sqrt(d) ;
347  
348 }
349 //______________________________________________________________________________
350 void AliPHOSEmcRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits)
351 {
352   // This function calculates energy in the core, 
353   // i.e. within a radius rad = 3cm around the center. Beyond this radius
354   // in accordance with shower profile the energy deposition 
355   // should be less than 2%
356
357   Float_t coreRadius = 3 ;
358
359   Float_t x = 0 ;
360   Float_t z = 0 ;
361
362   AliPHOSDigit * digit ;
363
364   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
365   AliPHOSGeometry * phosgeom =  (AliPHOSGeometry*)gime->PHOSGeometry();
366     
367   Int_t iDigit;
368
369 // Calculates the center of gravity in the local PHOS-module coordinates 
370   Float_t wtot = 0;
371   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
372     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
373     Int_t relid[4] ;
374     Float_t xi ;
375     Float_t zi ;
376     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
377     phosgeom->RelPosInModule(relid, xi, zi);
378     Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
379     x    += xi * w ;
380     z    += zi * w ;
381     wtot += w ;
382   }
383   x /= wtot ;
384   z /= wtot ;
385
386
387   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
388     digit = (AliPHOSDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
389     Int_t relid[4] ;
390     Float_t xi ;
391     Float_t zi ;
392     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
393     phosgeom->RelPosInModule(relid, xi, zi);    
394     Float_t distance = TMath::Sqrt((xi-x)*(xi-x)+(zi-z)*(zi-z)) ;
395     if(distance < coreRadius)
396       fCoreEnergy += fEnergyList[iDigit] ;
397   }
398
399 }
400
401 //____________________________________________________________________________
402 void  AliPHOSEmcRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
403 {
404   // Calculates the axis of the shower ellipsoid
405
406   Double_t wtot = 0. ;
407   Double_t x    = 0.;
408   Double_t z    = 0.;
409   Double_t dxx  = 0.;
410   Double_t dzz  = 0.;
411   Double_t dxz  = 0.;
412
413   AliPHOSDigit * digit ;
414
415   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
416   AliPHOSGeometry * phosgeom =  (AliPHOSGeometry*)gime->PHOSGeometry();
417
418   Int_t iDigit;
419
420
421   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
422     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit])  ;
423     Int_t relid[4] ;
424     Float_t xi ;
425     Float_t zi ;
426     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
427     phosgeom->RelPosInModule(relid, xi, zi);
428     Double_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
429     dxx  += w * xi * xi ;
430     x    += w * xi ;
431     dzz  += w * zi * zi ;
432     z    += w * zi ; 
433     dxz  += w * xi * zi ; 
434     wtot += w ;
435   }
436   dxx /= wtot ;
437   x   /= wtot ;
438   dxx -= x * x ;
439   dzz /= wtot ;
440   z   /= wtot ;
441   dzz -= z * z ;
442   dxz /= wtot ;
443   dxz -= x * z ;
444
445 //   //Apply correction due to non-perpendicular incidence
446 //   Double_t CosX ;
447 //   Double_t CosZ ;
448 //   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
449 //   AliPHOSGeometry * phosgeom =  (AliPHOSGeometry*)gime->PHOSGeometry();
450   //   Double_t DistanceToIP= (Double_t ) phosgeom->GetIPtoCrystalSurface() ;
451   
452 //   CosX = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+x*x) ;
453 //   CosZ = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+z*z) ;
454
455 //   dxx = dxx/(CosX*CosX) ;
456 //   dzz = dzz/(CosZ*CosZ) ;
457 //   dxz = dxz/(CosX*CosZ) ;
458
459
460   fLambda[0] =  0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
461   if(fLambda[0] > 0)
462     fLambda[0] = TMath::Sqrt(fLambda[0]) ;
463
464   fLambda[1] =  0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
465   if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
466     fLambda[1] = TMath::Sqrt(fLambda[1]) ;
467   else
468     fLambda[1]= 0. ;
469 }
470
471 //____________________________________________________________________________
472 void AliPHOSEmcRecPoint::EvalAll(Float_t logWeight, TClonesArray * digits )
473 {
474   // Evaluates all shower parameters
475
476   AliPHOSRecPoint::EvalAll(logWeight,digits) ;
477   EvalLocalPosition(logWeight, digits) ;
478   EvalElipsAxis(logWeight, digits) ;
479   EvalDispersion(logWeight, digits) ;
480   EvalCoreEnergy(logWeight, digits);
481 }
482 //____________________________________________________________________________
483 void AliPHOSEmcRecPoint::EvalLocalPosition(Float_t logWeight, TClonesArray * digits)
484 {
485   // Calculates the center of gravity in the local PHOS-module coordinates 
486   Float_t wtot = 0. ;
487  
488   Int_t relid[4] ;
489
490   Float_t x = 0. ;
491   Float_t z = 0. ;
492   
493   AliPHOSDigit * digit ;
494
495   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
496   AliPHOSGeometry * phosgeom =  (AliPHOSGeometry*)gime->PHOSGeometry();
497
498   Int_t iDigit;
499
500   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
501     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
502
503     Float_t xi ;
504     Float_t zi ;
505     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
506     phosgeom->RelPosInModule(relid, xi, zi);
507     Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
508     x    += xi * w ;
509     z    += zi * w ;
510     wtot += w ;
511
512   }
513
514   x /= wtot ;
515   z /= wtot ;
516
517   // Correction for the depth of the shower starting point (TDR p 127)  
518   Float_t para = 0.925 ; 
519   Float_t parb = 6.52 ; 
520
521   Float_t xo,yo,zo ; //Coordinates of the origin
522   gAlice->Generator()->GetOrigin(xo,yo,zo) ;
523
524   Float_t phi = phosgeom->GetPHOSAngle(relid[0]) ;
525
526   //Transform to the local ref.frame
527   Float_t xoL,yoL ;
528   xoL = xo*TMath::Cos(phi)-yo*TMath::Sin(phi) ;
529   yoL = xo*TMath::Sin(phi)+yo*TMath::Cos(phi) ;
530   
531   Float_t radius = TMath::Sqrt((xoL-x)*(xoL-x)+
532                                (phosgeom->GetIPtoCrystalSurface()-yoL)*(phosgeom->GetIPtoCrystalSurface()-yoL)+
533                                (zo-z)*(zo-z));
534  
535   Float_t incidencephi = TMath::ATan((x-xoL ) / radius) ; 
536   Float_t incidencetheta = TMath::ATan((z-zo) / radius) ;
537  
538   Float_t depthx =  ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencephi) ; 
539   Float_t depthz =  ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencetheta) ; 
540   
541
542   fLocPos.SetX(x - depthx)  ;
543   fLocPos.SetY(0.) ;
544   fLocPos.SetZ(z - depthz)  ;
545
546   fLocPosM = 0 ;
547 }
548
549 //____________________________________________________________________________
550 Float_t AliPHOSEmcRecPoint::GetMaximalEnergy(void) const
551 {
552   // Finds the maximum energy in the cluster
553   
554   Float_t menergy = 0. ;
555
556   Int_t iDigit;
557
558   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
559  
560     if(fEnergyList[iDigit] > menergy) 
561       menergy = fEnergyList[iDigit] ;
562   }
563   return menergy ;
564 }
565
566 //____________________________________________________________________________
567 Int_t AliPHOSEmcRecPoint::GetMultiplicityAtLevel(const Float_t H) const
568 {
569   // Calculates the multiplicity of digits with energy larger than H*energy 
570   
571   Int_t multipl   = 0 ;
572   Int_t iDigit ;
573   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
574
575     if(fEnergyList[iDigit] > H * fAmp) 
576       multipl++ ;
577   }
578   return multipl ;
579 }
580
581 //____________________________________________________________________________
582 Int_t  AliPHOSEmcRecPoint::GetNumberOfLocalMax(Int_t *  maxAt, Float_t * maxAtEnergy,
583                                                Float_t locMaxCut,TClonesArray * digits) const
584
585   // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
586   // energy difference between two local maxima
587
588   AliPHOSDigit * digit ;
589   AliPHOSDigit * digitN ;
590   
591
592   Int_t iDigitN ;
593   Int_t iDigit ;
594
595   for(iDigit = 0; iDigit < fMulDigit; iDigit++)
596     maxAt[iDigit] = (Int_t) digits->At(fDigitsList[iDigit])  ;
597
598   
599   for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {   
600     if(maxAt[iDigit] != -1) {
601       digit = (AliPHOSDigit *) maxAt[iDigit] ;
602           
603       for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {        
604         digitN = (AliPHOSDigit *) digits->At(fDigitsList[iDigitN]) ; 
605         
606         if ( AreNeighbours(digit, digitN) ) {
607           if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {    
608             maxAt[iDigitN] = -1 ;
609             // but may be digit too is not local max ?
610             if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut) 
611               maxAt[iDigit] = -1 ;
612           }
613           else {
614             maxAt[iDigit] = -1 ;
615             // but may be digitN too is not local max ?
616             if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut) 
617               maxAt[iDigitN] = -1 ; 
618           } 
619         } // if Areneighbours
620       } // while digitN
621     } // slot not empty
622   } // while digit
623   
624   iDigitN = 0 ;
625   for(iDigit = 0; iDigit < fMulDigit; iDigit++) { 
626     if(maxAt[iDigit] != -1){
627       maxAt[iDigitN] = maxAt[iDigit] ;
628       maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
629       iDigitN++ ; 
630     }
631   }
632   return iDigitN ;
633 }
634
635
636
637 //____________________________________________________________________________
638 void AliPHOSEmcRecPoint::Print(Option_t * option) 
639 {
640   // Print the list of digits belonging to the cluster
641   
642   cout << "AliPHOSEmcRecPoint: " << endl ;
643
644   Int_t iDigit;
645   cout << " digits # = " ;
646   for(iDigit=0; iDigit<fMulDigit; iDigit++)
647     cout << fDigitsList[iDigit] << "  " ;  
648   cout << endl ;
649   
650   cout << " Energies = " ;
651   for(iDigit=0; iDigit<fMulDigit; iDigit++) 
652     cout  << fEnergyList[iDigit] << "  ";
653   cout << endl ;
654   
655   cout << " Primaries  " ;
656   for(iDigit = 0;iDigit < fMulTrack; iDigit++)
657     cout << fTracksList[iDigit] << " " << endl ;
658         
659   cout << "       Multiplicity    = " << fMulDigit  << endl ;
660   cout << "       Cluster Energy  = " << fAmp << endl ;
661   cout << "       Number of primaries " << fMulTrack << endl ;
662   cout << "       Stored at position " << GetIndexInList() << endl ; 
663  
664 }
665  
666