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