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