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