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