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