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