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