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