]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSEmcRecPoint.cxx
Coding convention
[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 "TH2.h"
27 #include "TMath.h" 
28 #include "TCanvas.h" 
29 #include "TGraph.h"
30
31 // --- Standard library ---
32
33 // --- AliRoot header files ---
34 #include "AliPHOSLoader.h"
35 #include "AliGenerator.h"
36 #include "AliPHOSGeometry.h"
37 #include "AliPHOSDigit.h"
38 #include "AliPHOSEmcRecPoint.h"
39  
40 ClassImp(AliPHOSEmcRecPoint)
41
42 //____________________________________________________________________________
43 AliPHOSEmcRecPoint::AliPHOSEmcRecPoint() : AliPHOSRecPoint()
44 {
45   // ctor
46
47   fMulDigit   = 0 ;  
48   fAmp   = 0. ;   
49   fCoreEnergy = 0 ; 
50   fEnergyList = 0 ;
51   fNExMax     = 0 ;   //Not unfolded yet
52   fTime = -1. ;
53   fLocPos.SetX(1000000.)  ;      //Local position should be evaluated
54   fDebug=0;
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   fDebug=0;
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   AliPHOSGeometry * phosgeom =  AliPHOSLoader::GetPHOSGeometry();
150
151   Int_t relid1[4] ; 
152   phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ; 
153
154   Int_t relid2[4] ; 
155   phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ; 
156   
157   Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;  
158   Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;  
159
160   if (( coldiff <= 1 )  && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0)) 
161     aren = kTRUE ;
162   
163   return aren ;
164 }
165
166 //____________________________________________________________________________
167 Int_t AliPHOSEmcRecPoint::Compare(const TObject * obj) const
168 {
169   // Compares two RecPoints according to their position in the PHOS modules
170
171   Float_t delta = 1 ; //Width of "Sorting row". If you changibg this 
172                       //value (what is senseless) change as vell delta in
173                       //AliPHOSTrackSegmentMakerv* and other RecPoints...
174   Int_t rv ; 
175
176   AliPHOSEmcRecPoint * clu = (AliPHOSEmcRecPoint *)obj ; 
177
178  
179   Int_t phosmod1 = GetPHOSMod() ;
180   Int_t phosmod2 = clu->GetPHOSMod() ;
181
182   TVector3 locpos1; 
183   GetLocalPosition(locpos1) ;
184   TVector3 locpos2;  
185   clu->GetLocalPosition(locpos2) ;  
186
187   if(phosmod1 == phosmod2 ) {
188     Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
189     if (rowdif> 0) 
190       rv = 1 ;
191      else if(rowdif < 0) 
192        rv = -1 ;
193     else if(locpos1.Z()>locpos2.Z()) 
194       rv = -1 ;
195     else 
196       rv = 1 ; 
197      }
198
199   else {
200     if(phosmod1 < phosmod2 ) 
201       rv = -1 ;
202     else 
203       rv = 1 ;
204   }
205
206   return rv ; 
207 }
208 //______________________________________________________________________________
209 void AliPHOSEmcRecPoint::ExecuteEvent(Int_t event, Int_t, Int_t) const
210 {
211   
212   // Execute action corresponding to one event
213   //  This member function is called when a AliPHOSRecPoint is clicked with the locator
214   //
215   //  If Left button is clicked on AliPHOSRecPoint, the digits are switched on    
216   //  and switched off when the mouse button is released.
217   
218   
219   AliPHOSGeometry * phosgeom =  AliPHOSLoader::GetPHOSGeometry();
220   
221   static TGraph *  digitgraph = 0 ;
222   
223   if (!gPad->IsEditable()) return;
224   
225   TH2F * histo = 0 ;
226   TCanvas * histocanvas ; 
227
228   
229   //try to get run loader from default event folder
230   AliRunLoader* rn = AliRunLoader::GetRunLoader(AliConfig::GetDefaultEventFolderName());
231   if (rn == 0x0) 
232     {
233       Error("ExecuteEvent","Can not find Run Loader in Default Event Folder");
234       return;
235     }
236   AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(rn->GetLoader("PHOSLoader"));
237   if (gime == 0x0) 
238     {
239       Error("ExecuteEvent","Can not find PHOS Loader from Run Loader");
240       return;
241     }
242   
243   
244   const TClonesArray * digits = gime->Digits() ;
245   
246   switch (event) {
247     
248   case kButton1Down: {
249     AliPHOSDigit * digit ;
250     Int_t iDigit;
251     Int_t relid[4] ;
252     
253     const Int_t kMulDigit = AliPHOSEmcRecPoint::GetDigitsMultiplicity() ; 
254     Float_t * xi = new Float_t[kMulDigit] ; 
255     Float_t * zi = new Float_t[kMulDigit] ; 
256     
257     // create the histogram for the single cluster 
258     // 1. gets histogram boundaries
259     Float_t ximax = -999. ; 
260     Float_t zimax = -999. ; 
261     Float_t ximin = 999. ; 
262     Float_t zimin = 999. ;
263     
264     for(iDigit=0; iDigit<kMulDigit; iDigit++) {
265       digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit])  ;
266       phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
267       phosgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
268       if ( xi[iDigit] > ximax )
269         ximax = xi[iDigit] ; 
270       if ( xi[iDigit] < ximin )
271         ximin = xi[iDigit] ; 
272       if ( zi[iDigit] > zimax )
273         zimax = zi[iDigit] ; 
274       if ( zi[iDigit] < zimin )
275         zimin = zi[iDigit] ;     
276     }
277     ximax += phosgeom->GetCrystalSize(0) / 2. ;
278     zimax += phosgeom->GetCrystalSize(2) / 2. ;
279     ximin -= phosgeom->GetCrystalSize(0) / 2. ;
280     zimin -= phosgeom->GetCrystalSize(2) / 2. ;
281     Int_t xdim = (int)( (ximax - ximin ) / phosgeom->GetCrystalSize(0) + 0.5  ) ; 
282     Int_t zdim = (int)( (zimax - zimin ) / phosgeom->GetCrystalSize(2) + 0.5 ) ;
283     
284     // 2. gets the histogram title
285     
286     Text_t title[100] ; 
287     sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
288     
289     if (!histo) {
290       delete histo ; 
291       histo = 0 ; 
292     }
293     histo = new TH2F("cluster3D", title,  xdim, ximin, ximax, zdim, zimin, zimax)  ;
294     
295     Float_t x, z ; 
296     for(iDigit=0; iDigit<kMulDigit; iDigit++) {
297       digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit])  ;
298       phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
299       phosgeom->RelPosInModule(relid, x, z);
300       histo->Fill(x, z, fEnergyList[iDigit] ) ;
301     }
302     
303     if (!digitgraph) {
304       digitgraph = new TGraph(kMulDigit,xi,zi);
305       digitgraph-> SetMarkerStyle(5) ; 
306       digitgraph-> SetMarkerSize(1.) ;
307       digitgraph-> SetMarkerColor(1) ;
308       digitgraph-> Paint("P") ;
309     }
310     
311     //    Print() ;
312     histocanvas = new TCanvas("cluster", "a single cluster", 600, 500) ; 
313     histocanvas->Draw() ; 
314     histo->Draw("lego1") ; 
315     
316     delete[] xi ; 
317     delete[] zi ; 
318     
319     break;
320   }
321   
322   case kButton1Up: 
323     if (digitgraph) {
324       delete digitgraph  ;
325       digitgraph = 0 ;
326     }
327     break;
328   
329    }
330 }
331
332 //____________________________________________________________________________
333 void  AliPHOSEmcRecPoint::EvalDispersion(Float_t logWeight,TClonesArray * digits)
334 {
335   // Calculates the dispersion of the shower at the origine of the RecPoint
336
337   Float_t d    = 0. ;
338   Float_t wtot = 0. ;
339
340   Float_t x = 0.;
341   Float_t z = 0.;
342
343   AliPHOSDigit * digit ;
344  
345   AliPHOSGeometry * phosgeom =  AliPHOSLoader::GetPHOSGeometry();
346
347  // Calculates the center of gravity in the local PHOS-module coordinates 
348   
349   Int_t iDigit;
350
351   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
352     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
353     Int_t relid[4] ;
354     Float_t xi ;
355     Float_t zi ;
356     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
357     phosgeom->RelPosInModule(relid, xi, zi);
358     Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
359     x    += xi * w ;
360     z    += zi * w ;
361     wtot += w ;
362   }
363   x /= wtot ;
364   z /= wtot ;
365
366
367 // Calculates the dispersion in coordinates 
368   wtot = 0.;
369   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
370     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit])  ;
371     Int_t relid[4] ;
372     Float_t xi ;
373     Float_t zi ;
374     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
375     phosgeom->RelPosInModule(relid, xi, zi);
376     Float_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
377     d += w*((xi-x)*(xi-x) + (zi-z)*(zi-z) ) ; 
378     wtot+=w ;
379
380   
381   }
382   
383
384   d /= wtot ;
385
386   fDispersion = TMath::Sqrt(d) ;
387  
388 }
389 //______________________________________________________________________________
390 void AliPHOSEmcRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits)
391 {
392   // This function calculates energy in the core, 
393   // i.e. within a radius rad = 3cm around the center. Beyond this radius
394   // in accordance with shower profile the energy deposition 
395   // should be less than 2%
396
397   Float_t coreRadius = 3 ;
398
399   Float_t x = 0 ;
400   Float_t z = 0 ;
401
402   AliPHOSDigit * digit ;
403
404   AliPHOSGeometry * phosgeom =  AliPHOSLoader::GetPHOSGeometry();
405
406   Int_t iDigit;
407
408 // Calculates the center of gravity in the local PHOS-module coordinates 
409   Float_t wtot = 0;
410   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
411     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
412     Int_t relid[4] ;
413     Float_t xi ;
414     Float_t zi ;
415     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
416     phosgeom->RelPosInModule(relid, xi, zi);
417     Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
418     x    += xi * w ;
419     z    += zi * w ;
420     wtot += w ;
421   }
422   x /= wtot ;
423   z /= wtot ;
424
425
426   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
427     digit = (AliPHOSDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
428     Int_t relid[4] ;
429     Float_t xi ;
430     Float_t zi ;
431     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
432     phosgeom->RelPosInModule(relid, xi, zi);    
433     Float_t distance = TMath::Sqrt((xi-x)*(xi-x)+(zi-z)*(zi-z)) ;
434     if(distance < coreRadius)
435       fCoreEnergy += fEnergyList[iDigit] ;
436   }
437
438 }
439
440 //____________________________________________________________________________
441 void  AliPHOSEmcRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
442 {
443   // Calculates the axis of the shower ellipsoid
444
445   Double_t wtot = 0. ;
446   Double_t x    = 0.;
447   Double_t z    = 0.;
448   Double_t dxx  = 0.;
449   Double_t dzz  = 0.;
450   Double_t dxz  = 0.;
451
452   AliPHOSDigit * digit ;
453
454   AliPHOSGeometry * phosgeom =  AliPHOSLoader::GetPHOSGeometry();
455
456   Int_t iDigit;
457
458
459   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
460     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit])  ;
461     Int_t relid[4] ;
462     Float_t xi ;
463     Float_t zi ;
464     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
465     phosgeom->RelPosInModule(relid, xi, zi);
466     Double_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
467     dxx  += w * xi * xi ;
468     x    += w * xi ;
469     dzz  += w * zi * zi ;
470     z    += w * zi ; 
471     dxz  += w * xi * zi ; 
472     wtot += w ;
473   }
474   dxx /= wtot ;
475   x   /= wtot ;
476   dxx -= x * x ;
477   dzz /= wtot ;
478   z   /= wtot ;
479   dzz -= z * z ;
480   dxz /= wtot ;
481   dxz -= x * z ;
482
483 //   //Apply correction due to non-perpendicular incidence
484 //   Double_t CosX ;
485 //   Double_t CosZ ;
486 //   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
487 //   AliPHOSGeometry * phosgeom =  (AliPHOSGeometry*)gime->PHOSGeometry();
488   //   Double_t DistanceToIP= (Double_t ) phosgeom->GetIPtoCrystalSurface() ;
489   
490 //   CosX = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+x*x) ;
491 //   CosZ = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+z*z) ;
492
493 //   dxx = dxx/(CosX*CosX) ;
494 //   dzz = dzz/(CosZ*CosZ) ;
495 //   dxz = dxz/(CosX*CosZ) ;
496
497
498   fLambda[0] =  0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
499   if(fLambda[0] > 0)
500     fLambda[0] = TMath::Sqrt(fLambda[0]) ;
501
502   fLambda[1] =  0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
503   if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
504     fLambda[1] = TMath::Sqrt(fLambda[1]) ;
505   else
506     fLambda[1]= 0. ;
507 }
508
509 //____________________________________________________________________________
510 void  AliPHOSEmcRecPoint::EvalMoments(Float_t logWeight,TClonesArray * digits)
511 {
512   // Calculate the shower moments in the eigen reference system
513   // M2x, M2z, M3x, M4z
514   // Calculate the angle between the shower position vector and the eigen vector
515
516   Double_t wtot = 0. ;
517   Double_t x    = 0.;
518   Double_t z    = 0.;
519   Double_t dxx  = 0.;
520   Double_t dzz  = 0.;
521   Double_t dxz  = 0.;
522   Double_t lambda0=0, lambda1=0;
523
524   AliPHOSDigit * digit ;
525
526   AliPHOSGeometry * phosgeom =  AliPHOSLoader::GetPHOSGeometry();
527
528   Int_t iDigit;
529
530   // 1) Find covariance matrix elements:
531   //    || dxx dxz ||
532   //    || dxz dzz ||
533
534   Float_t xi ;
535   Float_t zi ;
536   Int_t relid[4] ;
537   Double_t w;
538   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
539     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit])  ;
540     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
541     phosgeom->RelPosInModule(relid, xi, zi);
542     w     = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
543     x    += w * xi ;
544     z    += w * zi ; 
545     dxx  += w * xi * xi ;
546     dzz  += w * zi * zi ;
547     dxz  += w * xi * zi ; 
548     wtot += w ;
549   }
550   x   /= wtot ;
551   z   /= wtot ;
552   dxx /= wtot ;
553   dzz /= wtot ;
554   dxz /= wtot ;
555   dxx -= x * x ;
556   dzz -= z * z ;
557   dxz -= x * z ;
558
559   // 2) Find covariance matrix eigen values lambda0 and lambda1
560
561   lambda0 =  0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
562   lambda1 =  0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
563
564   // 3) Find covariance matrix eigen vectors e0 and e1
565
566   TVector2 e0,e1;
567   if (dxz != 0)
568     e0.Set(1.,(lambda0-dxx)/dxz);
569   else
570     e0.Set(0.,1.);
571
572   e0 = e0.Unit();
573   e1.Set(-e0.Y(),e0.X());
574
575   // 4) Rotate cluster tensor from (x,z) to (e0,e1) system
576   //    and calculate moments M3x and M4z
577
578   Float_t cosPhi = e0.X();
579   Float_t sinPhi = e0.Y();
580
581   Float_t xiPHOS ;
582   Float_t ziPHOS ;
583   Double_t dx3, dz3, dz4;
584   wtot = 0.;
585   x    = 0.;
586   z    = 0.;
587   dxx  = 0.;
588   dzz  = 0.;
589   dxz  = 0.;
590   dx3  = 0.;
591   dz3  = 0.;
592   dz4  = 0.;
593   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
594     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit])  ;
595     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
596     phosgeom->RelPosInModule(relid, xiPHOS, ziPHOS);
597     xi    = xiPHOS*cosPhi + ziPHOS*sinPhi;
598     zi    = ziPHOS*cosPhi - xiPHOS*sinPhi;
599     w     = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
600     x    += w * xi ;
601     z    += w * zi ; 
602     dxx  += w * xi * xi ;
603     dzz  += w * zi * zi ;
604     dxz  += w * xi * zi ; 
605     dx3  += w * xi * xi * xi;
606     dz3  += w * zi * zi * zi ;
607     dz4  += w * zi * zi * zi * zi ;
608     wtot += w ;
609   }
610   x   /= wtot ;
611   z   /= wtot ;
612   dxx /= wtot ;
613   dzz /= wtot ;
614   dxz /= wtot ;
615   dx3 /= wtot ;
616   dz3 /= wtot ;
617   dz4 /= wtot ;
618   dx3 += -3*dxx*x + 2*x*x*x;
619   dz4 += -4*dz3*z + 6*dzz*z*z -3*z*z*z*z;
620   dxx -= x * x ;
621   dzz -= z * z ;
622   dxz -= x * z ;
623
624   // 5) Find an angle between cluster center vector and eigen vector e0
625
626   Float_t phi = TMath::ACos ((x*e0.X() + z*e0.Y()) / sqrt(x*x + z*z));
627
628   fM2x   = lambda0;
629   fM2z   = lambda1;
630   fM3x   = dx3;
631   fM4z   = dz4;
632   fPhixe = phi;
633
634 }
635 //____________________________________________________________________________
636 void AliPHOSEmcRecPoint::EvalAll(Float_t logWeight, TClonesArray * digits )
637 {
638   // Evaluates all shower parameters
639
640   EvalLocalPosition(logWeight, digits) ;
641   EvalElipsAxis(logWeight, digits) ;
642   EvalMoments(logWeight, digits) ;
643   EvalDispersion(logWeight, digits) ;
644   EvalCoreEnergy(logWeight, digits);
645   EvalTime(digits) ;
646   AliPHOSRecPoint::EvalAll(digits) ;
647 }
648 //____________________________________________________________________________
649 void AliPHOSEmcRecPoint::EvalLocalPosition(Float_t logWeight, TClonesArray * digits)
650 {
651   // Calculates the center of gravity in the local PHOS-module coordinates 
652   Float_t wtot = 0. ;
653  
654   Int_t relid[4] ;
655
656   Float_t x = 0. ;
657   Float_t z = 0. ;
658   
659   AliPHOSDigit * digit ;
660
661   AliPHOSGeometry * phosgeom =  AliPHOSLoader::GetPHOSGeometry();
662
663   Int_t iDigit;
664
665   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
666     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
667
668     Float_t xi ;
669     Float_t zi ;
670     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
671     phosgeom->RelPosInModule(relid, xi, zi);
672     Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
673     x    += xi * w ;
674     z    += zi * w ;
675     wtot += w ;
676
677   }
678
679   x /= wtot ;
680   z /= wtot ;
681
682   // Correction for the depth of the shower starting point (TDR p 127)  
683   Float_t para = 0.925 ; 
684   Float_t parb = 6.52 ; 
685
686   Float_t xo,yo,zo ; //Coordinates of the origin
687   gAlice->Generator()->GetOrigin(xo,yo,zo) ;
688
689   Float_t phi = phosgeom->GetPHOSAngle(relid[0]) ;
690
691   //Transform to the local ref.frame
692   Float_t xoL,yoL ;
693   xoL = xo*TMath::Cos(phi)-yo*TMath::Sin(phi) ;
694   yoL = xo*TMath::Sin(phi)+yo*TMath::Cos(phi) ;
695   
696   Float_t radius = phosgeom->GetIPtoCrystalSurface()-yoL;
697  
698   Float_t incidencephi = TMath::ATan((x-xoL ) / radius) ; 
699   Float_t incidencetheta = TMath::ATan((z-zo) / radius) ;
700  
701   Float_t depthx =  ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencephi) ; 
702   Float_t depthz =  ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencetheta) ; 
703
704   fLocPos.SetX(x - depthx)  ;
705   fLocPos.SetY(0.) ;
706   fLocPos.SetZ(z - depthz)  ;
707
708   fLocPosM = 0 ;
709 }
710
711 //____________________________________________________________________________
712 Float_t AliPHOSEmcRecPoint::GetMaximalEnergy(void) const
713 {
714   // Finds the maximum energy in the cluster
715   
716   Float_t menergy = 0. ;
717
718   Int_t iDigit;
719
720   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
721  
722     if(fEnergyList[iDigit] > menergy) 
723       menergy = fEnergyList[iDigit] ;
724   }
725   return menergy ;
726 }
727
728 //____________________________________________________________________________
729 Int_t AliPHOSEmcRecPoint::GetMultiplicityAtLevel(Float_t H) const
730 {
731   // Calculates the multiplicity of digits with energy larger than H*energy 
732   
733   Int_t multipl   = 0 ;
734   Int_t iDigit ;
735   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
736
737     if(fEnergyList[iDigit] > H * fAmp) 
738       multipl++ ;
739   }
740   return multipl ;
741 }
742
743 //____________________________________________________________________________
744 Int_t  AliPHOSEmcRecPoint::GetNumberOfLocalMax( AliPHOSDigit **  maxAt, Float_t * maxAtEnergy,
745                                                Float_t locMaxCut,TClonesArray * digits) const
746
747   // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
748   // energy difference between two local maxima
749
750   AliPHOSDigit * digit ;
751   AliPHOSDigit * digitN ;
752   
753
754   Int_t iDigitN ;
755   Int_t iDigit ;
756
757   for(iDigit = 0; iDigit < fMulDigit; iDigit++)
758     maxAt[iDigit] = (AliPHOSDigit*) digits->At(fDigitsList[iDigit])  ;
759
760   
761   for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {   
762     if(maxAt[iDigit]) {
763       digit = maxAt[iDigit] ;
764           
765       for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {        
766         digitN = (AliPHOSDigit *) digits->At(fDigitsList[iDigitN]) ; 
767         
768         if ( AreNeighbours(digit, digitN) ) {
769           if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {    
770             maxAt[iDigitN] = 0 ;
771             // but may be digit too is not local max ?
772             if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut) 
773               maxAt[iDigit] = 0 ;
774           }
775           else {
776             maxAt[iDigit] = 0 ;
777             // but may be digitN too is not local max ?
778             if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut) 
779               maxAt[iDigitN] = 0 ; 
780           } 
781         } // if Areneighbours
782       } // while digitN
783     } // slot not empty
784   } // while digit
785   
786   iDigitN = 0 ;
787   for(iDigit = 0; iDigit < fMulDigit; iDigit++) { 
788     if(maxAt[iDigit]){
789       maxAt[iDigitN] = maxAt[iDigit] ;
790       maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
791       iDigitN++ ; 
792     }
793   }
794   return iDigitN ;
795 }
796 //____________________________________________________________________________
797 void AliPHOSEmcRecPoint::EvalTime(TClonesArray * digits)
798 {
799   // Define a rec.point time as a time in the cell with the maximum energy
800
801   Float_t maxE = 0;
802   Int_t maxAt = 0;
803   for(Int_t idig=0; idig < fMulDigit; idig++){
804     if(fEnergyList[idig] > maxE){
805       maxE = fEnergyList[idig] ;
806       maxAt = idig;
807     }
808   }
809   fTime = ((AliPHOSDigit*) digits->At(fDigitsList[maxAt]))->GetTime() ;
810   
811 }
812 //____________________________________________________________________________
813 void AliPHOSEmcRecPoint::Purify(Float_t threshold){
814   //Removes digits below threshold
815
816   Int_t * tempo = new Int_t[fMaxDigit]; 
817   Float_t * tempoE =  new Float_t[fMaxDigit];
818
819   Int_t mult = 0 ;
820   for(Int_t iDigit=0;iDigit< fMulDigit ;iDigit++){
821     if(fEnergyList[iDigit] > threshold){
822       tempo[mult] = fDigitsList[iDigit] ;
823       tempoE[mult] = fEnergyList[iDigit] ;
824       mult++ ;
825     }
826   }
827   
828   fMulDigit = mult ;
829   delete [] fDigitsList ;
830   delete [] fEnergyList ;
831   fDigitsList = new Int_t[fMulDigit];
832   fEnergyList = new Float_t[fMulDigit];
833
834   fAmp = 0.; //Recalculate total energy
835   for(Int_t iDigit=0;iDigit< fMulDigit ;iDigit++){
836      fDigitsList[iDigit] = tempo[iDigit];
837      fEnergyList[iDigit] = tempoE[iDigit] ;
838      fAmp+=tempoE[iDigit];
839   }
840       
841   delete [] tempo ;
842   delete [] tempoE ;
843
844 }
845 //____________________________________________________________________________
846 void AliPHOSEmcRecPoint::Print(Option_t *) const
847 {
848   // Print the list of digits belonging to the cluster
849   
850   TString message ; 
851   message  = "AliPHOSEmcRecPoint:\n" ;
852   message +=  " digits # = " ; 
853   Info("Print", message.Data()) ; 
854
855   Int_t iDigit;
856   for(iDigit=0; iDigit<fMulDigit; iDigit++)
857     printf(" %d ", fDigitsList[iDigit] ) ;  
858   
859   Info("Print", " Energies = ") ;
860   for(iDigit=0; iDigit<fMulDigit; iDigit++) 
861     printf(" %f ", fEnergyList[iDigit] ) ;
862   printf("\n") ; 
863    Info("Print", " Primaries  ") ;
864   for(iDigit = 0;iDigit < fMulTrack; iDigit++)
865     printf(" %d ", fTracksList[iDigit]) ;
866   printf("\n") ;        
867   message  = "       Multiplicity    = %d" ;
868   message += "       Cluster Energy  = %f" ; 
869   message += "       Number of primaries %d" ; 
870   message += "       Stored at position %d" ; 
871  
872   Info("Print", message.Data(), fMulDigit, fAmp, fMulTrack,GetIndexInList() ) ;  
873 }
874  
875