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