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