]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALRecPoint.cxx
TRD1, 12 SM, 110 degree as default
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALRecPoint.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 /* $Id$ */
16 //_________________________________________________________________________
17 //  Reconstructed Points for the EMCAL
18 //  A RecPoint is a cluster of digits
19 //*-- Author: Yves Schutz (SUBATECH)
20 //*-- Author: Dmitri Peressounko (RRC KI & SUBATECH)
21 //*-- Author: Heather Gray (LBL) merged AliEMCALRecPoint and AliEMCALTowerRecPoint 02/04
22
23 // --- ROOT system ---
24 #include "TPad.h"
25 #include "TGraph.h"
26 #include "TPaveText.h"
27 #include "TClonesArray.h"
28 #include "TMath.h"
29
30 // --- Standard library ---
31
32 // --- AliRoot header files ---
33 #include "AliGenerator.h"
34 #include "AliEMCALGeometry.h"
35 #include "AliEMCALDigit.h"
36 #include "AliEMCALRecPoint.h"
37 #include "AliEMCALGetter.h"
38
39 ClassImp(AliEMCALRecPoint)
40
41 //____________________________________________________________________________
42 AliEMCALRecPoint::AliEMCALRecPoint()
43   : AliRecPoint()
44 {
45   // ctor
46   fMaxTrack = 0 ;
47   fMulDigit   = 0 ;  
48   fMaxParent = 0;
49   fMulParent = 0;
50   fAmp   = 0. ;   
51   fCoreEnergy = 0 ; 
52   fEnergyList = 0 ;
53   fParentsList = 0;
54   fTime = 0. ;
55   fLocPos.SetX(0.)  ;      //Local position should be evaluated
56   fCoreRadius = 10;        //HG Check this
57 }
58
59 //____________________________________________________________________________
60 AliEMCALRecPoint::AliEMCALRecPoint(const char * opt) : AliRecPoint(opt)
61 {
62   // ctor
63   fMaxTrack = 1000 ;
64   fMaxParent = 1000;
65   fMulDigit   = 0 ; 
66   fMulParent = 0; 
67   fAmp   = 0. ;   
68   fCoreEnergy = 0 ; 
69   fEnergyList = 0 ;
70   fParentsList = new Int_t[fMaxParent];
71   fTime = -1. ;
72   fLocPos.SetX(1000000.)  ;      //Local position should be evaluated
73   fCoreRadius = 10;        //HG Check this
74 }
75 //____________________________________________________________________________
76 AliEMCALRecPoint::~AliEMCALRecPoint()
77 {
78   // dtor
79   if ( fEnergyList )
80     delete[] fEnergyList ; 
81    if ( fParentsList)
82     delete[] fParentsList;
83 }
84
85 //____________________________________________________________________________
86 void AliEMCALRecPoint::AddDigit(AliEMCALDigit & digit, Float_t Energy)
87 {
88   // Adds a digit to the RecPoint
89   // and accumulates the total amplitude and the multiplicity 
90   
91   if(fEnergyList == 0)
92     fEnergyList =  new Float_t[fMaxDigit]; 
93
94   if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists 
95     fMaxDigit*=2 ; 
96     Int_t * tempo = new Int_t[fMaxDigit]; 
97     Float_t * tempoE =  new Float_t[fMaxDigit];
98
99     Int_t index ;     
100     for ( index = 0 ; index < fMulDigit ; index++ ){
101       tempo[index]  = fDigitsList[index] ;
102       tempoE[index] = fEnergyList[index] ; 
103     }
104     
105     delete [] fDigitsList ; 
106     fDigitsList =  new Int_t[fMaxDigit];
107  
108     delete [] fEnergyList ;
109     fEnergyList =  new Float_t[fMaxDigit];
110
111     for ( index = 0 ; index < fMulDigit ; index++ ){
112       fDigitsList[index] = tempo[index] ;
113       fEnergyList[index] = tempoE[index] ; 
114     }
115  
116     delete [] tempo ;
117     delete [] tempoE ; 
118   } // if
119   
120   fDigitsList[fMulDigit]   = digit.GetIndexInList()  ; 
121   fEnergyList[fMulDigit]   = Energy ;
122   fMulDigit++ ; 
123   fAmp += Energy ; 
124
125 }
126 //____________________________________________________________________________
127 Bool_t AliEMCALRecPoint::AreNeighbours(AliEMCALDigit * digit1, AliEMCALDigit * digit2 ) const
128 {
129   // Tells if (true) or not (false) two digits are neighbours
130   // A neighbour is defined as being two digits which share a corner
131   
132   Bool_t areNeighbours = kFALSE ;
133   
134   AliEMCALGeometry * geom =  (AliEMCALGetter::Instance())->EMCALGeometry();
135
136   Int_t relid1[2] ; 
137     //copied for shish-kebab geometry, ieta,iphi is cast as float as eta,phi conversion
138     // for this geometry does not exist
139     int nSupMod=0, nTower=0, nIphi=0, nIeta=0;
140     int iphi=0, ieta=0;
141        geom->GetCellIndex(digit1->GetId(), nSupMod,nTower,nIphi,nIeta);
142        geom->GetCellPhiEtaIndexInSModule(nSupMod,nTower,nIphi,nIeta, iphi,ieta);
143        relid1[0]=ieta;
144        relid1[1]=iphi;
145 //  geom->AbsToRelNumbering(digit1->GetId(), relid1) ; 
146
147   Int_t relid2[2] ; 
148     //copied for shish-kebab geometry, ieta,iphi is cast as float as eta,phi conversion
149     // for this geometry does not exist
150     int nSupMod1=0, nTower1=0, nIphi1=0, nIeta1=0;
151     int iphi1=0, ieta1=0;
152        geom->GetCellIndex(digit2->GetId(), nSupMod1,nTower1,nIphi1,nIeta1);
153        geom->GetCellPhiEtaIndexInSModule(nSupMod1,nTower1,nIphi1,nIeta1, iphi1,ieta1);
154        relid2[0]=ieta1;
155        relid2[1]=iphi1;
156 //  geom->AbsToRelNumbering(digit2->GetId(), relid2) ; 
157   
158   Int_t rowdiff = TMath::Abs( relid1[0] - relid2[0] ) ;  
159   Int_t coldiff = TMath::Abs( relid1[1] - relid2[1] ) ;  
160
161   if (( coldiff <= 1 )  && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0)) 
162     areNeighbours = kTRUE ;
163   
164   return areNeighbours;
165 }
166
167 //____________________________________________________________________________
168 Int_t AliEMCALRecPoint::Compare(const TObject * obj) const
169 {
170   // Compares two RecPoints according to their position in the EMCAL modules
171
172   Float_t delta = 1 ; //Width of "Sorting row". If you change this 
173                       //value (what is senseless) change as well delta in
174                       //AliEMCALTrackSegmentMakerv* and other RecPoints...
175   Int_t rv ; 
176
177   AliEMCALRecPoint * clu = (AliEMCALRecPoint *)obj ; 
178
179   TVector3 locpos1; 
180   GetLocalPosition(locpos1);
181   TVector3 locpos2;  
182   clu->GetLocalPosition(locpos2);  
183
184   Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
185   if (rowdif> 0) 
186     rv = 1 ;
187   else if(rowdif < 0) 
188     rv = -1 ;
189   else if(locpos1.Y()>locpos2.Y()) 
190     rv = -1 ;
191   else 
192     rv = 1 ; 
193
194   return rv ; 
195 }
196
197 //____________________________________________________________________________
198 Int_t AliEMCALRecPoint::DistancetoPrimitive(Int_t px, Int_t py)
199 {
200   // Compute distance from point px,py to  a AliEMCALRecPoint considered as a Tmarker
201   // Compute the closest distance of approach from point px,py to this marker.
202   // The distance is computed in pixels units.
203   // HG Still need to update -> Not sure what this should achieve
204
205   TVector3 pos(0.,0.,0.) ;
206   GetLocalPosition(pos) ;
207   Float_t x =  pos.X() ;
208   Float_t y =  pos.Y() ;
209   const Int_t kMaxDiff = 10;
210   Int_t pxm  = gPad->XtoAbsPixel(x);
211   Int_t pym  = gPad->YtoAbsPixel(y);
212   Int_t dist = (px-pxm)*(px-pxm) + (py-pym)*(py-pym);
213   
214   if (dist > kMaxDiff) return 9999;
215   return dist;
216 }
217
218 //___________________________________________________________________________
219  void AliEMCALRecPoint::Draw(Option_t *option)
220  {
221    // Draw this AliEMCALRecPoint with its current attributes
222    
223    AppendPad(option);
224  }
225
226 //______________________________________________________________________________
227 void AliEMCALRecPoint::ExecuteEvent(Int_t /*event*/, Int_t, Int_t)
228 {
229   // Execute action corresponding to one event
230   // This member function is called when a AliEMCALRecPoint is clicked with the locator
231   //
232   // If Left button is clicked on AliEMCALRecPoint, the digits are switched on    
233   // and switched off when the mouse button is released.
234
235   //  static Int_t pxold, pyold;
236
237   /* static TGraph *  digitgraph = 0 ;
238   static TPaveText* clustertext = 0 ;
239   
240   if (!gPad->IsEditable()) return;
241   
242   switch (event) {
243     
244     
245   case kButton1Down:{
246     AliEMCALDigit * digit ;
247     AliEMCALGeometry * emcalgeom =  (AliEMCALGetter::Instance())->EMCALGeometry() ;
248
249     Int_t iDigit;
250     Int_t relid[2] ;
251   
252     const Int_t kMulDigit=AliEMCALRecPoint::GetDigitsMultiplicity() ;
253     Float_t * xi = new Float_t [kMulDigit] ; 
254     Float_t * zi = new Float_t [kMulDigit] ;
255     
256     for(iDigit = 0; iDigit < kMulDigit; iDigit++) {
257       Fatal("AliEMCALRecPoint::ExecuteEvent", " -> Something wrong with the code"); 
258       digit = 0 ; //dynamic_cast<AliEMCALDigit *>((fDigitsList)[iDigit]);
259       emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
260       emcalgeom->PosInAlice(relid, xi[iDigit], zi[iDigit]) ;
261     }
262     
263     if (!digitgraph) {
264       digitgraph = new TGraph(fMulDigit,xi,zi);
265       digitgraph-> SetMarkerStyle(5) ; 
266       digitgraph-> SetMarkerSize(1.) ;
267       digitgraph-> SetMarkerColor(1) ;
268       digitgraph-> Draw("P") ;
269     }
270     if (!clustertext) {
271       
272       TVector3 pos(0.,0.,0.) ;
273       GetLocalPosition(pos) ;
274       clustertext = new TPaveText(pos.X()-10,pos.Z()+10,pos.X()+50,pos.Z()+35,"") ;
275       Text_t  line1[40] ;
276       Text_t  line2[40] ;
277       sprintf(line1,"Energy=%1.2f GeV",GetEnergy()) ;
278       sprintf(line2,"%d Digits",GetDigitsMultiplicity()) ;
279       clustertext ->AddText(line1) ;
280       clustertext ->AddText(line2) ;
281       clustertext ->Draw("");
282     }
283     gPad->Update() ; 
284     Print("") ;
285     delete[] xi ; 
286     delete[] zi ; 
287    }
288   
289 break;
290   
291   case kButton1Up:
292     if (digitgraph) {
293       delete digitgraph  ;
294       digitgraph = 0 ;
295     }
296     if (clustertext) {
297       delete clustertext ;
298       clustertext = 0 ;
299     }
300     
301     break;
302     
303     }*/
304 }
305 //____________________________________________________________________________
306 void AliEMCALRecPoint::EvalAll(Float_t logWeight,TClonesArray * digits) 
307 {
308   // Evaluates all shower parameters
309
310   EvalLocalPosition(logWeight, digits) ;
311   //  printf("eval position done\n");
312   EvalElipsAxis(logWeight, digits) ;
313   //  printf("eval axis done\n");
314   EvalDispersion(logWeight, digits) ;
315   //  printf("eval dispersion done\n");
316  // EvalCoreEnergy(logWeight, digits);
317  // printf("eval energy done\n");
318   EvalTime(digits) ;
319   //  printf("eval time done\n");
320
321   EvalPrimaries(digits) ;
322   //  printf("eval pri done\n");
323   EvalParents(digits);
324   //  printf("eval parent done\n");
325 }
326
327 //____________________________________________________________________________
328 void  AliEMCALRecPoint::EvalDispersion(Float_t logWeight, TClonesArray * digits)
329 {
330   // Calculates the dispersion of the shower at the origin of the RecPoint
331
332   Float_t d    = 0. ;
333   Float_t wtot = 0. ;
334
335   AliEMCALDigit * digit ;
336  
337   AliEMCALGeometry * geom = (AliEMCALGetter::Instance())->EMCALGeometry();
338   
339   // Calculates the centre of gravity in the local EMCAL-module coordinates   
340   Int_t iDigit;
341
342   if (!fLocPos.X() || !fLocPos.Y()) 
343     EvalLocalPosition(logWeight, digits) ;
344   
345   //Sub const Float_t kDeg2Rad = TMath::DegToRad() ; 
346     
347   Float_t cluEta =  fLocPos.X() ; 
348   Float_t cluPhi =  fLocPos.Y() ; 
349   Float_t cluR =  fLocPos.Z() ; 
350   
351   if (gDebug == 2) 
352     printf("EvalDispersion: eta,phi,r = %f,%f,%f", cluEta, cluPhi, cluR) ;
353   
354   // Calculates the dispersion in coordinates 
355   wtot = 0.;
356   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
357     digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit])  ;
358     Float_t etai = 0.;
359     Float_t phii = 0.;
360     //copied for shish-kebab geometry, ieta,iphi is cast as float as eta,phi conversion
361     // for this geometry does not exist
362     int nSupMod=0, nTower=0, nIphi=0, nIeta=0;
363     int iphi=0, ieta=0;
364        geom->GetCellIndex(digit->GetId(), nSupMod,nTower,nIphi,nIeta);
365        geom->GetCellPhiEtaIndexInSModule(nSupMod,nTower,nIphi,nIeta, iphi,ieta);
366         etai=(Float_t)ieta;
367         phii=(Float_t)iphi;
368         //        printf("%f,%d,%d \n", fAmp, ieta, iphi) ;
369
370 /* Sub
371   geom->EtaPhiFromIndex(digit->GetId(), etai, phii);
372     phii = phii * kDeg2Rad;
373 */
374 ///////////////////////////
375 //  if(fAmp>0)printf("%f %d %f", fAmp,iDigit,fEnergyList[iDigit]) ;
376 /////////////////////////
377
378     if (gDebug == 2) 
379       printf("EvalDispersion: id = %d, etai,phii = %f,%f", digit->GetId(), etai, phii) ;
380
381     Float_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
382     d += w * ( (etai-cluEta)*(etai-cluEta) + (phii-cluPhi)*(phii-cluPhi) ) ; 
383     wtot+=w ;
384   }
385   
386   if ( wtot > 0 ) 
387     d /= wtot ;
388   else 
389     d = 0. ; 
390
391   fDispersion = TMath::Sqrt(d) ;
392   //      printf("Dispersion: = %f", fDispersion) ;
393 }
394
395 //____________________________________________________________________________
396 void AliEMCALRecPoint::EvalLocalPosition(Float_t logWeight, TClonesArray * digits)
397 {
398   // Calculates the center of gravity in the local EMCAL-module coordinates 
399   Float_t wtot = 0. ;
400  
401   //  Int_t relid[3] ;
402   
403   AliEMCALDigit * digit ;
404   AliEMCALGeometry * geom  =  (AliEMCALGetter::Instance())->EMCALGeometry();
405   Int_t iDigit;
406   Float_t cluEta = 0;
407   Float_t cluPhi = 0;
408  //Sub const Float_t kDeg2Rad = TMath::DegToRad();
409
410   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
411     digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
412
413     Float_t etai ;
414     Float_t phii ;
415     //copied for shish-kebab geometry, ieta,iphi is cast as float as eta,phi conversion
416     // for this geometry does not exist
417     int nSupMod=0, nTower=0, nIphi=0, nIeta=0;
418     int iphi=0, ieta=0;
419        geom->GetCellIndex(digit->GetId(), nSupMod,nTower,nIphi,nIeta);
420        geom->GetCellPhiEtaIndexInSModule(nSupMod, nTower,nIphi,nIeta, iphi,ieta); //19-oct-05
421         etai=(Float_t)ieta;
422         phii=(Float_t)iphi;
423 //Sub    geom->EtaPhiFromIndex(digit->GetId(), etai, phii);
424 //Sub    phii = phii * kDeg2Rad;
425     Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
426     cluEta += (etai * w) ;
427     cluPhi += (phii * w );
428     wtot += w ;
429   }
430
431   if ( wtot > 0 ) { 
432     cluEta /= wtot ;
433     cluPhi /= wtot ;
434   } else {
435     cluEta = -1 ; 
436     cluPhi = -1.; 
437   }
438   
439   fLocPos.SetX(cluEta);
440   fLocPos.SetY(cluPhi);
441   fLocPos.SetZ(geom->GetIP2ECASection());
442     
443 //  if (gDebug==2)
444 //    printf("EvalLocalPosition: eta,phi,r = %f,%f,%f", fLocPos.X(), fLocPos.Y(), fLocPos.Z()) ; 
445   fLocPosM = 0 ;
446 }
447
448 //______________________________________________________________________________
449 void AliEMCALRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits)
450 {
451   // This function calculates energy in the core, 
452   // i.e. within a radius rad = 3cm around the center. Beyond this radius
453   // in accordance with shower profile the energy deposition 
454   // should be less than 2%
455
456   AliEMCALDigit * digit ;
457   const Float_t kDeg2Rad = TMath::DegToRad() ;
458   AliEMCALGeometry * geom = (AliEMCALGetter::Instance())->EMCALGeometry();    
459   Int_t iDigit;
460
461   if (!fLocPos.X() || !fLocPos.Y() ) {
462     EvalLocalPosition(logWeight, digits);
463   }
464   
465   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
466     digit = (AliEMCALDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
467     Float_t etai = 0. ;
468     Float_t phii = 0. ;
469     geom->PosInAlice(digit->GetId(), etai, phii);
470     phii = phii * kDeg2Rad;
471   
472     Float_t distance = TMath::Sqrt((etai-fLocPos.X())*(etai-fLocPos.X())+(phii-fLocPos.Y())*(phii-fLocPos.Y())) ;
473     if(distance < fCoreRadius)
474       fCoreEnergy += fEnergyList[iDigit] ;
475   }
476   
477 }
478 //____________________________________________________________________________
479 void  AliEMCALRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
480 {
481   // Calculates the axis of the shower ellipsoid in eta and phi
482
483   Double_t wtot = 0. ;
484   Double_t x    = 0.;
485   Double_t z    = 0.;
486   Double_t dxx  = 0.;
487   Double_t dzz  = 0.;
488   Double_t dxz  = 0.;
489
490   const Float_t kDeg2Rad = TMath::DegToRad();
491   AliEMCALDigit * digit ;
492
493   AliEMCALGeometry * geom = (AliEMCALGetter::Instance())->EMCALGeometry();
494   TString gn(geom->GetName());
495
496   Int_t iDigit;
497   
498   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
499     digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit])  ;
500     Float_t etai = 0. ;
501     Float_t phii = 0. ; 
502     if(gn.Contains("SHISH")) {
503     //copied for shish-kebab geometry, ieta,iphi is cast as float as eta,phi conversion
504     // for this geometry does not exist
505       int nSupMod=0, nTower=0, nIphi=0, nIeta=0;
506       int iphi=0, ieta=0;
507       geom->GetCellIndex(digit->GetId(), nSupMod,nTower,nIphi,nIeta);
508       geom->GetCellPhiEtaIndexInSModule(nSupMod,nTower,nIphi,nIeta, iphi,ieta);
509       etai=(Float_t)ieta;
510       phii=(Float_t)iphi;
511     } else {
512       geom->EtaPhiFromIndex(digit->GetId(), etai, phii);
513       phii = phii * kDeg2Rad;
514     }
515
516     Double_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
517
518     dxx  += w * etai * etai ;
519     x    += w * etai ;
520     dzz  += w * phii * phii ;
521     z    += w * phii ; 
522
523     dxz  += w * etai * phii ; 
524
525     wtot += w ;
526   }
527
528   if ( wtot > 0 ) { 
529     dxx /= wtot ;
530     x   /= wtot ;
531     dxx -= x * x ;
532     dzz /= wtot ;
533     z   /= wtot ;
534     dzz -= z * z ;
535     dxz /= wtot ;
536     dxz -= x * z ;
537
538     fLambda[0] =  0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
539     if(fLambda[0] > 0)
540       fLambda[0] = TMath::Sqrt(fLambda[0]) ;
541     else
542       fLambda[0] = 0;
543     
544     fLambda[1] =  0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
545
546     if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
547       fLambda[1] = TMath::Sqrt(fLambda[1]) ;
548     else
549       fLambda[1]= 0. ;
550   } else { 
551     fLambda[0]= 0. ;
552     fLambda[1]= 0. ;
553   }
554
555   //  printf("Evalaxis: lambdas  = %f,%f", fLambda[0],fLambda[1]) ; 
556
557 }
558
559 //______________________________________________________________________________
560 void  AliEMCALRecPoint::EvalPrimaries(TClonesArray * digits)
561 {
562   // Constructs the list of primary particles (tracks) which have contributed to this RecPoint
563   
564   AliEMCALDigit * digit ;
565   Int_t * tempo    = new Int_t[fMaxTrack] ;
566
567   Int_t index ;  
568   for ( index = 0 ; index < GetDigitsMultiplicity() ; index++ ) { // all digits
569     digit = dynamic_cast<AliEMCALDigit *>(digits->At( fDigitsList[index] )) ; 
570     Int_t nprimaries = digit->GetNprimary() ;
571     if ( nprimaries == 0 ) continue ;
572     Int_t * newprimaryarray = new Int_t[nprimaries] ;
573     Int_t ii ; 
574     for ( ii = 0 ; ii < nprimaries ; ii++)
575       newprimaryarray[ii] = digit->GetPrimary(ii+1) ; 
576
577     Int_t jndex ;
578     for ( jndex = 0 ; jndex < nprimaries ; jndex++ ) { // all primaries in digit
579       if ( fMulTrack > fMaxTrack ) {
580         fMulTrack = fMaxTrack ;
581         Error("GetNprimaries", "increase fMaxTrack ")  ;
582         break ;
583       }
584       Int_t newprimary = newprimaryarray[jndex] ;
585       Int_t kndex ;
586       Bool_t already = kFALSE ;
587       for ( kndex = 0 ; kndex < fMulTrack ; kndex++ ) { //check if not already stored
588         if ( newprimary == tempo[kndex] ){
589           already = kTRUE ;
590           break ;
591         }
592       } // end of check
593       if ( !already && (fMulTrack < fMaxTrack)) { // store it
594         tempo[fMulTrack] = newprimary ; 
595         fMulTrack++ ;
596       } // store it
597     } // all primaries in digit
598     delete [] newprimaryarray ; 
599   } // all digits
600
601   
602   fTracksList = new Int_t[fMulTrack] ;
603   for(index = 0; index < fMulTrack; index++)
604    fTracksList[index] = tempo[index] ;
605  
606   delete [] tempo ;
607
608 }
609
610 //______________________________________________________________________________
611 void  AliEMCALRecPoint::EvalParents(TClonesArray * digits)
612 {
613   // Constructs the list of parent particles (tracks) which have contributed to this RecPoint
614  
615   AliEMCALDigit * digit ;
616   Int_t * tempo    = new Int_t[fMaxParent] ;
617
618   Int_t index ;  
619   for ( index = 0 ; index < GetDigitsMultiplicity() ; index++ ) { // all digits
620     digit = dynamic_cast<AliEMCALDigit *>(digits->At( fDigitsList[index] )) ; 
621     Int_t nparents = digit->GetNiparent() ;
622     if ( nparents == 0 ) continue ;
623     Int_t * newparentarray = new Int_t[nparents] ;
624     Int_t ii ; 
625     for ( ii = 0 ; ii < nparents ; ii++)
626       newparentarray[ii] = digit->GetIparent(ii+1) ; 
627
628     Int_t jndex ;
629     for ( jndex = 0 ; jndex < nparents ; jndex++ ) { // all primaries in digit
630       if ( fMulParent > fMaxParent ) {
631         fMulTrack = - 1 ;
632         Error("GetNiparent", "increase fMaxParent")  ;
633         break ;
634       }
635       Int_t newparent = newparentarray[jndex] ;
636       Int_t kndex ;
637       Bool_t already = kFALSE ;
638       for ( kndex = 0 ; kndex < fMulParent ; kndex++ ) { //check if not already stored
639         if ( newparent == tempo[kndex] ){
640           already = kTRUE ;
641           break ;
642         }
643       } // end of check
644       if ( !already && (fMulTrack < fMaxTrack)) { // store it
645         tempo[fMulParent] = newparent ; 
646         fMulParent++ ;
647       } // store it
648     } // all parents in digit
649     delete [] newparentarray ; 
650   } // all digits
651
652   if (fMulParent>0) {
653     fParentsList = new Int_t[fMulParent] ;
654     for(index = 0; index < fMulParent; index++)
655       fParentsList[index] = tempo[index] ;
656   }
657  
658   delete [] tempo ;
659
660 }
661
662 //____________________________________________________________________________
663 void AliEMCALRecPoint::GetLocalPosition(TVector3 & lpos) const
664 {
665   // returns the position of the cluster in the local reference system of ALICE
666   // X = eta, Y = phi, Z = r (a constant for the EMCAL)
667   
668   lpos.SetX(fLocPos.X()) ;
669   lpos.SetY(fLocPos.Y()) ;
670   lpos.SetZ(fLocPos.Z()) ;
671 }
672
673 //____________________________________________________________________________
674 void AliEMCALRecPoint::GetGlobalPosition(TVector3 & gpos) const
675 {
676   // returns the position of the cluster in the global reference system of ALICE
677   // These are now the Cartesian X, Y and Z
678
679   AliEMCALGeometry * geom = (AliEMCALGetter::Instance())->EMCALGeometry();
680   Int_t absid = geom->TowerIndexFromEtaPhi(fLocPos.X(), TMath::RadToDeg()*fLocPos.Y());
681   geom->XYZFromIndex(absid, gpos);
682 }
683
684 //____________________________________________________________________________
685 Float_t AliEMCALRecPoint::GetMaximalEnergy(void) const
686 {
687   // Finds the maximum energy in the cluster
688   
689   Float_t menergy = 0. ;
690
691   Int_t iDigit;
692
693   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
694  
695     if(fEnergyList[iDigit] > menergy) 
696       menergy = fEnergyList[iDigit] ;
697   }
698   return menergy ;
699 }
700
701 //____________________________________________________________________________
702 Int_t AliEMCALRecPoint::GetMultiplicityAtLevel(Float_t H) const
703 {
704   // Calculates the multiplicity of digits with energy larger than H*energy 
705   
706   Int_t multipl   = 0 ;
707   Int_t iDigit ;
708   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
709
710     if(fEnergyList[iDigit] > H * fAmp) 
711       multipl++ ;
712   }
713   return multipl ;
714 }
715
716 //____________________________________________________________________________
717 Int_t  AliEMCALRecPoint::GetNumberOfLocalMax(AliEMCALDigit **  maxAt, Float_t * maxAtEnergy,
718                                                Float_t locMaxCut,TClonesArray * digits) const
719
720   // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
721   // energy difference between two local maxima
722
723   AliEMCALDigit * digit ;
724   AliEMCALDigit * digitN ;
725   
726   Int_t iDigitN ;
727   Int_t iDigit ;
728
729   for(iDigit = 0; iDigit < fMulDigit; iDigit++)
730     maxAt[iDigit] = (AliEMCALDigit*) digits->At(fDigitsList[iDigit])  ;
731   
732   for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {   
733     if(maxAt[iDigit]) {
734       digit = maxAt[iDigit] ;
735           
736       for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {        
737         digitN = (AliEMCALDigit *) digits->At(fDigitsList[iDigitN]) ; 
738         
739         if ( AreNeighbours(digit, digitN) ) {
740           if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {    
741             maxAt[iDigitN] = 0 ;
742             // but may be digit too is not local max ?
743             if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut) 
744               maxAt[iDigit] = 0 ;
745           }
746           else {
747             maxAt[iDigit] = 0 ;
748             // but may be digitN too is not local max ?
749             if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut) 
750               maxAt[iDigitN] = 0 ; 
751           } 
752         } // if Areneighbours
753       } // while digitN
754     } // slot not empty
755   } // while digit
756   
757   iDigitN = 0 ;
758   for(iDigit = 0; iDigit < fMulDigit; iDigit++) { 
759     if(maxAt[iDigit] ){
760       maxAt[iDigitN] = maxAt[iDigit] ;
761       maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
762       iDigitN++ ; 
763     }
764   }
765   return iDigitN ;
766 }
767 //____________________________________________________________________________
768 void AliEMCALRecPoint::EvalTime(TClonesArray * digits){
769   // time is set to the time of the digit with the maximum energy
770
771   Float_t maxE = 0;
772   Int_t maxAt = 0;
773   for(Int_t idig=0; idig < fMulDigit; idig++){
774     if(fEnergyList[idig] > maxE){
775       maxE = fEnergyList[idig] ;
776       maxAt = idig;
777     }
778   }
779   fTime = ((AliEMCALDigit*) digits->At(fDigitsList[maxAt]))->GetTime() ;
780   
781 }
782
783 //______________________________________________________________________________
784 void AliEMCALRecPoint::Paint(Option_t *)
785 {
786   // Paint this ALiRecPoint as a TMarker  with its current attributes
787   
788   TVector3 pos(0.,0.,0.)  ;
789   GetLocalPosition(pos)   ;
790   Coord_t x = pos.X()     ;
791   Coord_t y = pos.Z()     ;
792   Color_t markercolor = 1 ;
793   Size_t  markersize = 1. ;
794   Style_t markerstyle = 5 ;
795   
796   if (!gPad->IsBatch()) {
797     gVirtualX->SetMarkerColor(markercolor) ;
798     gVirtualX->SetMarkerSize (markersize)  ;
799     gVirtualX->SetMarkerStyle(markerstyle) ;
800   }
801   gPad->SetAttMarkerPS(markercolor,markerstyle,markersize) ;
802   gPad->PaintPolyMarker(1,&x,&y,"") ;
803 }
804
805 //______________________________________________________________________________
806 Float_t AliEMCALRecPoint::EtaToTheta(Float_t arg) const
807 {
808   //Converts Theta (Radians) to Eta(Radians)
809   return (2.*TMath::ATan(TMath::Exp(-arg)));
810 }
811
812 //______________________________________________________________________________
813 Float_t AliEMCALRecPoint::ThetaToEta(Float_t arg) const
814 {
815   //Converts Eta (Radians) to Theta(Radians)
816   return (-1 * TMath::Log(TMath::Tan(0.5 * arg)));
817 }
818
819 //____________________________________________________________________________
820 void AliEMCALRecPoint::Print(Option_t *) const
821 {
822   // Print the list of digits belonging to the cluster
823   
824   TString message ; 
825   message  = "AliPHOSEmcRecPoint:\n" ;
826   message +=  " digits # = " ; 
827   Info("Print", message.Data()) ; 
828
829   Int_t iDigit;
830   for(iDigit=0; iDigit<fMulDigit; iDigit++)
831     printf(" %d ", fDigitsList[iDigit] ) ;  
832   
833   Info("Print", " Energies = ") ;
834   for(iDigit=0; iDigit<fMulDigit; iDigit++) 
835     printf(" %f ", fEnergyList[iDigit] ) ;
836   printf("\n") ; 
837    Info("Print", " Primaries  ") ;
838   for(iDigit = 0;iDigit < fMulTrack; iDigit++)
839     printf(" %d ", fTracksList[iDigit]) ;
840   printf("\n") ;        
841   message  = "       Multiplicity    = %d" ;
842   message += "       Cluster Energy  = %f" ; 
843   message += "       Core energy     = %f" ; 
844   message += "       Core radius     = %f" ; 
845   message += "       Number of primaries %d" ; 
846   message += "       Stored at position %d" ; 
847  
848   Info("Print", message.Data(), fMulDigit, fAmp, fCoreEnergy, fCoreRadius, fMulTrack, GetIndexInList() ) ;  
849 }