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