]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSEmcRecPoint.cxx
Radius for core energy is set via AliPHOSRecoParam.
[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, 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 coreRadius = AliPHOSReconstructor::GetRecoParamEmc()->GetEcoreRadius() ;
433
434   Float_t x = 0 ;
435   Float_t z = 0 ;
436
437   AliPHOSDigit * digit ;
438
439   AliPHOSGeometry * phosgeom =  AliPHOSGeometry::GetInstance();
440
441   Int_t iDigit;
442
443 // Calculates the center of gravity in the local PHOS-module coordinates 
444   Float_t wtot = 0;
445   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
446     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
447     Int_t relid[4] ;
448     Float_t xi ;
449     Float_t zi ;
450     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
451     phosgeom->RelPosInModule(relid, xi, zi);
452     if (fAmp>0 && fEnergyList[iDigit]>0) {
453       Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
454       x    += xi * w ;
455       z    += zi * w ;
456       wtot += w ;
457     }
458     else
459       AliError(Form("Wrong energy %f and/or amplitude %f\n", fEnergyList[iDigit], fAmp));
460   }
461   if (wtot>0) {
462     x /= wtot ;
463     z /= wtot ;
464   }
465   else
466     AliError(Form("Wrong weight %f\n", wtot));
467
468
469   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
470     digit = (AliPHOSDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
471     Int_t relid[4] ;
472     Float_t xi ;
473     Float_t zi ;
474     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
475     phosgeom->RelPosInModule(relid, xi, zi);    
476     Float_t distance = TMath::Sqrt((xi-x)*(xi-x)+(zi-z)*(zi-z)) ;
477     if(distance < coreRadius)
478       fCoreEnergy += fEnergyList[iDigit] ;
479   }
480
481
482 }
483
484 //____________________________________________________________________________
485 void  AliPHOSEmcRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits, TVector3 & /* vInc */)
486 {
487   // Calculates the axis of the shower ellipsoid
488
489   Double_t wtot = 0. ;
490   Double_t x    = 0.;
491   Double_t z    = 0.;
492   Double_t dxx  = 0.;
493   Double_t dzz  = 0.;
494   Double_t dxz  = 0.;
495
496   AliPHOSDigit * digit ;
497
498   AliPHOSGeometry * phosgeom =  AliPHOSGeometry::GetInstance();
499
500   Int_t iDigit;
501
502
503   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
504     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit])  ;
505     Int_t relid[4] ;
506     Float_t xi ;
507     Float_t zi ;
508     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
509     phosgeom->RelPosInModule(relid, xi, zi);
510     if (fAmp>0 && fEnergyList[iDigit]>0) {
511       Double_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
512       dxx  += w * xi * xi ;
513       x    += w * xi ;
514       dzz  += w * zi * zi ;
515       z    += w * zi ; 
516       dxz  += w * xi * zi ; 
517       wtot += w ;
518     }
519     else
520       AliError(Form("Wrong energy %f and/or amplitude %f\n", fEnergyList[iDigit], fAmp));
521   }
522   if (wtot>0) {
523     dxx /= wtot ;
524     x   /= wtot ;
525     dxx -= x * x ;
526     dzz /= wtot ;
527     z   /= wtot ;
528     dzz -= z * z ;
529     dxz /= wtot ;
530     dxz -= x * z ;
531
532 //   //Apply correction due to non-perpendicular incidence
533 //   Double_t CosX ;
534 //   Double_t CosZ ;
535 //   AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance() ;
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 }
850
851 //____________________________________________________________________________
852 Float_t AliPHOSEmcRecPoint::GetMaximalEnergy(void) const
853 {
854   // Finds the maximum energy in the cluster
855   
856   Float_t menergy = 0. ;
857
858   Int_t iDigit;
859
860   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
861  
862     if(fEnergyList[iDigit] > menergy) 
863       menergy = fEnergyList[iDigit] ;
864   }
865   return menergy ;
866 }
867
868 //____________________________________________________________________________
869 Int_t AliPHOSEmcRecPoint::GetMultiplicityAtLevel(Float_t H) const
870 {
871   // Calculates the multiplicity of digits with energy larger than H*energy 
872   
873   Int_t multipl   = 0 ;
874   Int_t iDigit ;
875   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
876
877     if(fEnergyList[iDigit] > H * fAmp) 
878       multipl++ ;
879   }
880   return multipl ;
881 }
882
883 //____________________________________________________________________________
884 Int_t  AliPHOSEmcRecPoint::GetNumberOfLocalMax( AliPHOSDigit **  maxAt, Float_t * maxAtEnergy,
885                                                Float_t locMaxCut,TClonesArray * digits) const
886
887   // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
888   // energy difference between two local maxima
889
890   AliPHOSDigit * digit ;
891   AliPHOSDigit * digitN ;
892   
893
894   Int_t iDigitN ;
895   Int_t iDigit ;
896
897   for(iDigit = 0; iDigit < fMulDigit; iDigit++)
898     maxAt[iDigit] = (AliPHOSDigit*) digits->At(fDigitsList[iDigit])  ;
899
900   
901   for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {   
902     if(maxAt[iDigit]) {
903       digit = maxAt[iDigit] ;
904           
905       for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {        
906         if(iDigit == iDigitN)
907           continue ;
908         
909         digitN = (AliPHOSDigit *) digits->At(fDigitsList[iDigitN]) ; 
910         
911         if ( AreNeighbours(digit, digitN) ) {
912           if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {    
913             maxAt[iDigitN] = 0 ;
914             // but may be digit too is not local max ?
915             if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut) 
916               maxAt[iDigit] = 0 ;
917           }
918           else {
919             maxAt[iDigit] = 0 ;
920             // but may be digitN too is not local max ?
921             if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut) 
922               maxAt[iDigitN] = 0 ; 
923           } 
924         } // if Areneighbours
925       } // while digitN
926     } // slot not empty
927   } // while digit
928   
929   iDigitN = 0 ;
930   for(iDigit = 0; iDigit < fMulDigit; iDigit++) { 
931     if(maxAt[iDigit]){
932       maxAt[iDigitN] = maxAt[iDigit] ;
933       maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
934       iDigitN++ ; 
935     }
936   }
937   return iDigitN ;
938 }
939 //____________________________________________________________________________
940 void AliPHOSEmcRecPoint::EvalTime(TClonesArray * digits)
941 {
942   // Define a rec.point time as a time in the cell with the maximum energy
943
944   Float_t maxE = 0;
945   Int_t maxAt = 0;
946   for(Int_t idig=0; idig < fMulDigit; idig++){
947     if(fEnergyList[idig] > maxE){
948       maxE = fEnergyList[idig] ;
949       maxAt = idig;
950     }
951   }
952   fTime = ((AliPHOSDigit*) digits->At(fDigitsList[maxAt]))->GetTime() ;
953   
954 }
955 //____________________________________________________________________________
956 void AliPHOSEmcRecPoint::Purify(Float_t threshold){
957   //Removes digits below threshold
958
959   Int_t * tempo = new Int_t[fMaxDigit]; 
960   Float_t * tempoE =  new Float_t[fMaxDigit];
961
962   Int_t mult = 0 ;
963   for(Int_t iDigit=0;iDigit< fMulDigit ;iDigit++){
964     if(fEnergyList[iDigit] > threshold){
965       tempo[mult] = fDigitsList[iDigit] ;
966       tempoE[mult] = fEnergyList[iDigit] ;
967       mult++ ;
968     }
969   }
970   
971   fMulDigit = mult ;
972   delete [] fDigitsList ;
973   delete [] fEnergyList ;
974   fDigitsList = new Int_t[fMulDigit];
975   fEnergyList = new Float_t[fMulDigit];
976
977   fAmp = 0.; //Recalculate total energy
978   for(Int_t iDigit=0;iDigit< fMulDigit ;iDigit++){
979      fDigitsList[iDigit] = tempo[iDigit];
980      fEnergyList[iDigit] = tempoE[iDigit] ;
981      fAmp+=tempoE[iDigit];
982   }
983       
984   delete [] tempo ;
985   delete [] tempoE ;
986
987 }
988 //____________________________________________________________________________
989 void AliPHOSEmcRecPoint::Print(Option_t *) const
990 {
991   // Print the list of digits belonging to the cluster
992   
993   TString message ; 
994   message  = "AliPHOSEmcRecPoint:\n" ;
995   message += "Digit multiplicity = %d" ;
996   message += ", cluster Energy = %f" ; 
997   message += ", number of primaries = %d" ; 
998   message += ", rec.point index = %d \n" ; 
999   printf(message.Data(), fMulDigit, fAmp, fMulTrack,GetIndexInList() ) ;  
1000
1001   Int_t iDigit;
1002   printf(" digits ids = ") ; 
1003   for(iDigit=0; iDigit<fMulDigit; iDigit++)
1004     printf(" %d ", fDigitsList[iDigit] ) ;  
1005   
1006   printf("\n digit energies = ") ;
1007   for(iDigit=0; iDigit<fMulDigit; iDigit++) 
1008     printf(" %f ", fEnergyList[iDigit] ) ;
1009
1010   printf("\n digit primaries  ") ;
1011   if (fMulTrack<1) printf("... no primaries");
1012   for(iDigit = 0;iDigit < fMulTrack; iDigit++)
1013     printf(" %d ", fTracksList[iDigit]) ;
1014   printf("\n") ;        
1015
1016 }
1017  
1018