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