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