Using AliESDCaloCluster instead of AliESDtrack
[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 AliEMCALGeometry * geom = 0;
42
43 //____________________________________________________________________________
44 AliEMCALRecPoint::AliEMCALRecPoint()
45   : AliRecPoint()
46 {
47   // ctor
48   fClusterType = -1;
49   fMaxTrack = 0 ;
50   fMulDigit   = 0 ;  
51   fMaxParent = 0;
52   fMulParent = 0;
53   fAmp   = 0. ;   
54   fCoreEnergy = 0 ; 
55   fEnergyList = 0 ;
56   fTimeList = 0 ;
57   fAbsIdList  = 0;
58   fParentsList = 0;
59   fTime = 0. ;
60   //  fLocPos.SetX(1.e+6)  ;      //Local position should be evaluated
61   fCoreRadius = 10;        //HG Check this
62   geom = AliEMCALGeometry::GetInstance();
63   geom->GetTransformationForSM(); // Global <-> Local
64 }
65
66 //____________________________________________________________________________
67 AliEMCALRecPoint::AliEMCALRecPoint(const char * opt) : AliRecPoint(opt)
68 {
69   // ctor
70   fClusterType = -1;
71   fMaxTrack = 1000 ;
72   fMaxParent = 1000;
73   fMulDigit   = 0 ; 
74   fMulParent = 0; 
75   fAmp   = 0. ;   
76   fCoreEnergy = 0 ; 
77   fEnergyList = 0 ;
78   fTimeList = 0 ;
79   fAbsIdList  = 0;
80   fParentsList = new Int_t[fMaxParent];
81   fTime = -1. ;
82   //fLocPos.SetX(1.e+6)  ;      //Local position should be evaluated
83   fCoreRadius = 10;        //HG Check this
84   geom = AliEMCALGeometry::GetInstance();
85   geom->GetTransformationForSM(); // Global <-> Local
86 }
87 //____________________________________________________________________________
88 AliEMCALRecPoint::~AliEMCALRecPoint()
89 {
90   // dtor
91   if ( fEnergyList )
92     delete[] fEnergyList ; 
93   if ( fTimeList )
94     delete[] fTimeList ; 
95   if ( fAbsIdList )
96     delete[] fAbsIdList ; 
97    if ( fParentsList)
98     delete[] fParentsList;
99 }
100
101 //____________________________________________________________________________
102 void AliEMCALRecPoint::AddDigit(AliEMCALDigit & digit, Float_t Energy)
103 {
104   // Adds a digit to the RecPoint
105   // and accumulates the total amplitude and the multiplicity 
106   
107   if(fEnergyList == 0)
108     fEnergyList =  new Float_t[fMaxDigit]; 
109   if(fTimeList == 0)
110     fTimeList =  new Float_t[fMaxDigit]; 
111   if(fAbsIdList == 0) {
112     fAbsIdList =  new Int_t[fMaxDigit];
113     fSuperModuleNumber = geom->GetSuperModuleNumber(digit.GetId());
114   }
115
116   if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists 
117     fMaxDigit*=2 ; 
118     Int_t * tempo = new Int_t[fMaxDigit]; 
119     Float_t * tempoE =  new Float_t[fMaxDigit];
120     Float_t * tempoT =  new Float_t[fMaxDigit];
121     Int_t * tempoId = new Int_t[fMaxDigit]; 
122
123     Int_t index ;     
124     for ( index = 0 ; index < fMulDigit ; index++ ){
125       tempo[index]   = fDigitsList[index] ;
126       tempoE[index]  = fEnergyList[index] ; 
127       tempoT[index]  = fTimeList[index] ; 
128       tempoId[index] = fAbsIdList[index] ; 
129     }
130     
131     delete [] fDigitsList ; 
132     fDigitsList =  new Int_t[fMaxDigit];
133  
134     delete [] fEnergyList ;
135     fEnergyList =  new Float_t[fMaxDigit];
136
137     delete [] fTimeList ;
138     fTimeList =  new Float_t[fMaxDigit];
139
140     delete [] fAbsIdList ;
141     fAbsIdList =  new Int_t[fMaxDigit];
142
143     for ( index = 0 ; index < fMulDigit ; index++ ){
144       fDigitsList[index] = tempo[index] ;
145       fEnergyList[index] = tempoE[index] ; 
146       fTimeList[index] = tempoT[index] ; 
147       fAbsIdList[index]  = tempoId[index] ; 
148     }
149  
150     delete [] tempo ;
151     delete [] tempoE ; 
152     delete [] tempoT ; 
153     delete [] tempoId ; 
154   } // if
155   
156   fDigitsList[fMulDigit]   = digit.GetIndexInList()  ; 
157   fEnergyList[fMulDigit]   = Energy ;
158   fTimeList[fMulDigit]     = digit.GetTime() ;
159   fAbsIdList[fMulDigit]    = digit.GetId();
160   fMulDigit++ ; 
161   fAmp += Energy ; 
162
163 }
164 //____________________________________________________________________________
165 Bool_t AliEMCALRecPoint::AreNeighbours(AliEMCALDigit * digit1, AliEMCALDigit * digit2 ) const
166 {
167   // Tells if (true) or not (false) two digits are neighbours
168   // A neighbour is defined as being two digits which share a corner
169   
170   static Bool_t areNeighbours = kFALSE ;
171   static Int_t nSupMod=0, nTower=0, nIphi=0, nIeta=0;
172   static int nSupMod1=0, nTower1=0, nIphi1=0, nIeta1=0;
173   static Int_t relid1[2] , relid2[2] ; // ieta, iphi
174   static Int_t rowdiff=0, coldiff=0;
175
176   areNeighbours = kFALSE ;
177
178   geom->GetCellIndex(digit1->GetId(), nSupMod,nTower,nIphi,nIeta);
179   geom->GetCellPhiEtaIndexInSModule(nSupMod,nTower,nIphi,nIeta, relid1[0],relid1[1]);
180
181   geom->GetCellIndex(digit2->GetId(), nSupMod1,nTower1,nIphi1,nIeta1);
182   geom->GetCellPhiEtaIndexInSModule(nSupMod1,nTower1,nIphi1,nIeta1, relid2[0],relid2[1]);
183   
184   rowdiff = TMath::Abs( relid1[0] - relid2[0] ) ;  
185   coldiff = TMath::Abs( relid1[1] - relid2[1] ) ;  
186
187   if (( coldiff <= 1 )  && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0)) 
188   areNeighbours = kTRUE ;
189   
190   return areNeighbours;
191 }
192
193 //____________________________________________________________________________
194 Int_t AliEMCALRecPoint::Compare(const TObject * obj) const
195 {
196   // Compares two RecPoints according to their position in the EMCAL modules
197
198   Float_t delta = 1 ; //Width of "Sorting row". If you change this 
199                       //value (what is senseless) change as well delta in
200                       //AliEMCALTrackSegmentMakerv* and other RecPoints...
201   Int_t rv ; 
202
203   AliEMCALRecPoint * clu = (AliEMCALRecPoint *)obj ; 
204
205   TVector3 locpos1; 
206   GetLocalPosition(locpos1);
207   TVector3 locpos2;  
208   clu->GetLocalPosition(locpos2);  
209
210   Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
211   if (rowdif> 0) 
212     rv = 1 ;
213   else if(rowdif < 0) 
214     rv = -1 ;
215   else if(locpos1.Y()>locpos2.Y()) 
216     rv = -1 ;
217   else 
218     rv = 1 ; 
219
220   return rv ; 
221 }
222
223 //____________________________________________________________________________
224 Int_t AliEMCALRecPoint::DistancetoPrimitive(Int_t px, Int_t py)
225 {
226   // Compute distance from point px,py to  a AliEMCALRecPoint considered as a Tmarker
227   // Compute the closest distance of approach from point px,py to this marker.
228   // The distance is computed in pixels units.
229   // HG Still need to update -> Not sure what this should achieve
230
231   TVector3 pos(0.,0.,0.) ;
232   GetLocalPosition(pos) ;
233   Float_t x =  pos.X() ;
234   Float_t y =  pos.Y() ;
235   const Int_t kMaxDiff = 10;
236   Int_t pxm  = gPad->XtoAbsPixel(x);
237   Int_t pym  = gPad->YtoAbsPixel(y);
238   Int_t dist = (px-pxm)*(px-pxm) + (py-pym)*(py-pym);
239   
240   if (dist > kMaxDiff) return 9999;
241   return dist;
242 }
243
244 //___________________________________________________________________________
245  void AliEMCALRecPoint::Draw(Option_t *option)
246  {
247    // Draw this AliEMCALRecPoint with its current attributes
248    
249    AppendPad(option);
250  }
251
252 //______________________________________________________________________________
253 void AliEMCALRecPoint::ExecuteEvent(Int_t /*event*/, Int_t, Int_t)
254 {
255   // Execute action corresponding to one event
256   // This member function is called when a AliEMCALRecPoint is clicked with the locator
257   //
258   // If Left button is clicked on AliEMCALRecPoint, the digits are switched on    
259   // and switched off when the mouse button is released.
260
261   //  static Int_t pxold, pyold;
262
263   /* static TGraph *  digitgraph = 0 ;
264   static TPaveText* clustertext = 0 ;
265   
266   if (!gPad->IsEditable()) return;
267   
268   switch (event) {
269     
270     
271   case kButton1Down:{
272     AliEMCALDigit * digit ;
273     AliEMCALGeometry * emcalgeom =  (AliEMCALGetter::Instance())->EMCALGeometry() ;
274
275     Int_t iDigit;
276     Int_t relid[2] ;
277   
278     const Int_t kMulDigit=AliEMCALRecPoint::GetDigitsMultiplicity() ;
279     Float_t * xi = new Float_t [kMulDigit] ; 
280     Float_t * zi = new Float_t [kMulDigit] ;
281     
282     for(iDigit = 0; iDigit < kMulDigit; iDigit++) {
283       Fatal("AliEMCALRecPoint::ExecuteEvent", " -> Something wrong with the code"); 
284       digit = 0 ; //dynamic_cast<AliEMCALDigit *>((fDigitsList)[iDigit]);
285       emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
286       emcalgeom->PosInAlice(relid, xi[iDigit], zi[iDigit]) ;
287     }
288     
289     if (!digitgraph) {
290       digitgraph = new TGraph(fMulDigit,xi,zi);
291       digitgraph-> SetMarkerStyle(5) ; 
292       digitgraph-> SetMarkerSize(1.) ;
293       digitgraph-> SetMarkerColor(1) ;
294       digitgraph-> Draw("P") ;
295     }
296     if (!clustertext) {
297       
298       TVector3 pos(0.,0.,0.) ;
299       GetLocalPosition(pos) ;
300       clustertext = new TPaveText(pos.X()-10,pos.Z()+10,pos.X()+50,pos.Z()+35,"") ;
301       Text_t  line1[40] ;
302       Text_t  line2[40] ;
303       sprintf(line1,"Energy=%1.2f GeV",GetEnergy()) ;
304       sprintf(line2,"%d Digits",GetDigitsMultiplicity()) ;
305       clustertext ->AddText(line1) ;
306       clustertext ->AddText(line2) ;
307       clustertext ->Draw("");
308     }
309     gPad->Update() ; 
310     Print("") ;
311     delete[] xi ; 
312     delete[] zi ; 
313    }
314   
315 break;
316   
317   case kButton1Up:
318     if (digitgraph) {
319       delete digitgraph  ;
320       digitgraph = 0 ;
321     }
322     if (clustertext) {
323       delete clustertext ;
324       clustertext = 0 ;
325     }
326     
327     break;
328     
329     }*/
330 }
331 //____________________________________________________________________________
332 void AliEMCALRecPoint::EvalAll(Float_t logWeight,TClonesArray * digits) 
333 {
334   // Evaluates all shower parameters
335
336   EvalLocalPosition(logWeight, digits) ;
337   //  printf("eval position done\n");
338   EvalElipsAxis(logWeight, digits) ;
339   //  printf("eval axis done\n");
340   EvalDispersion(logWeight, digits) ;
341   //  printf("eval dispersion done\n");
342  // EvalCoreEnergy(logWeight, digits);
343  // printf("eval energy done\n");
344   EvalTime(digits) ;
345   //  printf("eval time done\n");
346
347   EvalPrimaries(digits) ;
348   //  printf("eval pri done\n");
349   EvalParents(digits);
350   //  printf("eval parent done\n");
351 }
352
353 //____________________________________________________________________________
354 void  AliEMCALRecPoint::EvalDispersion(Float_t logWeight, TClonesArray * digits)
355 {
356   // Calculates the dispersion of the shower at the origin of the RecPoint
357
358   Double_t d = 0., wtot = 0., w = 0., xyzi[3], diff=0.;
359   Int_t iDigit=0, nstat=0, i=0;
360   AliEMCALDigit * digit ;
361
362   // Calculates the centre of gravity in the local EMCAL-module coordinates   
363   if (!fLocPos.Mag()) 
364     EvalLocalPosition(logWeight, digits) ;
365   
366   // Calculates the dispersion in coordinates 
367   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
368     digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit])  ;
369
370     geom->RelPosCellInSModule(digit->GetId(), xyzi[0], xyzi[1], xyzi[2]);
371     w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
372
373     if(w>0.0) {
374       wtot += w;
375       nstat++;
376       for(i=0; i<3; i++ ) {
377         diff = xyzi[i] - double(fLocPos[i]);
378         d   += w * diff*diff;
379       }
380     }
381   }
382   
383   if ( wtot > 0 && nstat>1) d /= wtot ;
384   else                      d = 0. ; 
385
386   fDispersion = TMath::Sqrt(d) ;
387 }
388
389 //____________________________________________________________________________
390 void AliEMCALRecPoint::EvalLocalPosition(Float_t logWeight, TClonesArray * digits)
391 {
392   // Calculates the center of gravity in the local EMCAL-module coordinates 
393   //  Info("Print", " logWeight %f : cluster energy %f ", logWeight, fAmp); // for testing
394   
395   AliEMCALDigit * digit;
396   Int_t i=0, nstat=0;
397   Double_t clXYZ[3]={0.,0.,0.}, clRmsXYZ[3]={0.,0.,0.}, xyzi[3], wtot=0., w=0.;
398
399   for(Int_t iDigit=0; iDigit<fMulDigit; iDigit++) {
400     digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
401
402     geom->RelPosCellInSModule(digit->GetId(), xyzi[0], xyzi[1], xyzi[2]);
403     // printf(" Id %i : Local x,y,z %f %f %f \n", digit->GetId(), xyzi[0], xyzi[1], xyzi[2]);
404
405     if(logWeight > 0.0) w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ));
406     else                w = fEnergyList[iDigit]; // just energy
407
408     if(w>0.0) {
409       wtot += w ;
410       nstat++;
411       for(i=0; i<3; i++ ) {
412         clXYZ[i]    += (w*xyzi[i]);
413         clRmsXYZ[i] += (w*xyzi[i]*xyzi[i]);
414       }
415     }
416   }
417   //  cout << " wtot " << wtot << endl;
418   if ( wtot > 0 ) { 
419     //    xRMS   = TMath::Sqrt(x2m - xMean*xMean);
420     for(i=0; i<3; i++ ) {
421       clXYZ[i] /= wtot;
422       if(nstat>1) {
423         clRmsXYZ[i] /= (wtot*wtot);  
424         clRmsXYZ[i]  =  clRmsXYZ[i] - clXYZ[i]*clXYZ[i];
425         if(clRmsXYZ[i] > 0.0) {
426           clRmsXYZ[i] =  TMath::Sqrt(clRmsXYZ[i]);
427         } else      clRmsXYZ[i] = 0;
428       } else        clRmsXYZ[i] = 0;
429     }    
430   } else {
431     for(i=0; i<3; i++ ) {
432       clXYZ[i] = clRmsXYZ[i] = -1.;
433     }
434   }
435   // clRmsXYZ[i] ??
436   fLocPos.SetX(clXYZ[0]);
437   fLocPos.SetY(clXYZ[1]);
438   fLocPos.SetZ(clXYZ[2]);
439     
440 //  if (gDebug==2)
441 //    printf("EvalLocalPosition: eta,phi,r = %f,%f,%f", fLocPos.X(), fLocPos.Y(), fLocPos.Z()) ; 
442   fLocPosM = 0 ; // covariance matrix
443 }
444
445 //void AliEMCALRecPoint::EvalLocalPositionSimple()
446 //{ // Weight is proportional of cell energy 
447 //}
448
449 //______________________________________________________________________________
450 void AliEMCALRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits)
451 {
452   // This function calculates energy in the core, 
453   // i.e. within a radius rad = 3cm around the center. Beyond this radius
454   // in accordance with shower profile the energy deposition 
455   // should be less than 2%
456
457   AliEMCALDigit * digit ;
458   const Float_t kDeg2Rad = TMath::DegToRad() ;
459
460   Int_t iDigit;
461
462   if (!fLocPos.Mag()) {
463     EvalLocalPosition(logWeight, digits);
464   }
465   
466   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
467     digit = (AliEMCALDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
468     Float_t etai = 0. ;
469     Float_t phii = 0. ;
470     geom->PosInAlice(digit->GetId(), etai, phii);
471     phii = phii * kDeg2Rad;
472   
473     Float_t distance = TMath::Sqrt((etai-fLocPos.X())*(etai-fLocPos.X())+(phii-fLocPos.Y())*(phii-fLocPos.Y())) ;
474     if(distance < fCoreRadius)
475       fCoreEnergy += fEnergyList[iDigit] ;
476   }
477   
478 }
479 //____________________________________________________________________________
480 void  AliEMCALRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
481 {
482   // Calculates the axis of the shower ellipsoid in eta and phi
483
484   Double_t wtot = 0. ;
485   Double_t x    = 0.;
486   Double_t z    = 0.;
487   Double_t dxx  = 0.;
488   Double_t dzz  = 0.;
489   Double_t dxz  = 0.;
490
491   const Float_t kDeg2Rad = TMath::DegToRad();
492   AliEMCALDigit * digit ;
493
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")) { // have to be change - Feb 28, 2006
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   
667   lpos.SetX(fLocPos.X()) ;
668   lpos.SetY(fLocPos.Y()) ;
669   lpos.SetZ(fLocPos.Z()) ;
670 }
671
672 //____________________________________________________________________________
673 void AliEMCALRecPoint::GetGlobalPosition(TVector3 & gpos) const
674 {
675   // returns the position of the cluster in the global reference system of ALICE
676   // These are now the Cartesian X, Y and Z
677   geom = AliEMCALGeometry::GetInstance();
678   //  cout<<" geom "<<geom<<endl;
679   geom->GetGlobal(fLocPos, gpos, fSuperModuleNumber);
680 }
681
682 //____________________________________________________________________________
683 Float_t AliEMCALRecPoint::GetMaximalEnergy(void) const
684 {
685   // Finds the maximum energy in the cluster
686   
687   Float_t menergy = 0. ;
688
689   Int_t iDigit;
690
691   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
692  
693     if(fEnergyList[iDigit] > menergy) 
694       menergy = fEnergyList[iDigit] ;
695   }
696   return menergy ;
697 }
698
699 //____________________________________________________________________________
700 Int_t AliEMCALRecPoint::GetMultiplicityAtLevel(Float_t H) const
701 {
702   // Calculates the multiplicity of digits with energy larger than H*energy 
703   
704   Int_t multipl   = 0 ;
705   Int_t iDigit ;
706   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
707
708     if(fEnergyList[iDigit] > H * fAmp) 
709       multipl++ ;
710   }
711   return multipl ;
712 }
713
714 //____________________________________________________________________________
715 Int_t  AliEMCALRecPoint::GetNumberOfLocalMax(AliEMCALDigit **  maxAt, Float_t * maxAtEnergy,
716                                                Float_t locMaxCut,TClonesArray * digits) const
717
718   // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
719   // energy difference between two local maxima
720
721   AliEMCALDigit * digit ;
722   AliEMCALDigit * digitN ;
723   
724   Int_t iDigitN ;
725   Int_t iDigit ;
726
727   for(iDigit = 0; iDigit < fMulDigit; iDigit++)
728     maxAt[iDigit] = (AliEMCALDigit*) digits->At(fDigitsList[iDigit])  ;
729   
730   for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {   
731     if(maxAt[iDigit]) {
732       digit = maxAt[iDigit] ;
733           
734       for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {        
735         digitN = (AliEMCALDigit *) digits->At(fDigitsList[iDigitN]) ; 
736         
737         if ( AreNeighbours(digit, digitN) ) {
738           if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {    
739             maxAt[iDigitN] = 0 ;
740             // but may be digit too is not local max ?
741             if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut) 
742               maxAt[iDigit] = 0 ;
743           }
744           else {
745             maxAt[iDigit] = 0 ;
746             // but may be digitN too is not local max ?
747             if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut) 
748               maxAt[iDigitN] = 0 ; 
749           } 
750         } // if Areneighbours
751       } // while digitN
752     } // slot not empty
753   } // while digit
754   
755   iDigitN = 0 ;
756   for(iDigit = 0; iDigit < fMulDigit; iDigit++) { 
757     if(maxAt[iDigit] ){
758       maxAt[iDigitN] = maxAt[iDigit] ;
759       maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
760       iDigitN++ ; 
761     }
762   }
763   return iDigitN ;
764 }
765 //____________________________________________________________________________
766 void AliEMCALRecPoint::EvalTime(TClonesArray * digits){
767   // time is set to the time of the digit with the maximum energy
768
769   Float_t maxE = 0;
770   Int_t maxAt = 0;
771   for(Int_t idig=0; idig < fMulDigit; idig++){
772     if(fEnergyList[idig] > maxE){
773       maxE = fEnergyList[idig] ;
774       maxAt = idig;
775     }
776   }
777   fTime = ((AliEMCALDigit*) digits->At(fDigitsList[maxAt]))->GetTime() ;
778   
779 }
780
781 //______________________________________________________________________________
782 void AliEMCALRecPoint::Paint(Option_t *)
783 {
784   // Paint this ALiRecPoint as a TMarker  with its current attributes
785   
786   TVector3 pos(0.,0.,0.)  ;
787   GetLocalPosition(pos)   ;
788   Coord_t x = pos.X()     ;
789   Coord_t y = pos.Z()     ;
790   Color_t markercolor = 1 ;
791   Size_t  markersize = 1. ;
792   Style_t markerstyle = 5 ;
793   
794   if (!gPad->IsBatch()) {
795     gVirtualX->SetMarkerColor(markercolor) ;
796     gVirtualX->SetMarkerSize (markersize)  ;
797     gVirtualX->SetMarkerStyle(markerstyle) ;
798   }
799   gPad->SetAttMarkerPS(markercolor,markerstyle,markersize) ;
800   gPad->PaintPolyMarker(1,&x,&y,"") ;
801 }
802
803 //______________________________________________________________________________
804 Float_t AliEMCALRecPoint::EtaToTheta(Float_t arg) const
805 {
806   //Converts Theta (Radians) to Eta(Radians)
807   return (2.*TMath::ATan(TMath::Exp(-arg)));
808 }
809
810 //______________________________________________________________________________
811 Float_t AliEMCALRecPoint::ThetaToEta(Float_t arg) const
812 {
813   //Converts Eta (Radians) to Theta(Radians)
814   return (-1 * TMath::Log(TMath::Tan(0.5 * arg)));
815 }
816
817 //____________________________________________________________________________
818 void AliEMCALRecPoint::Print(Option_t *) const
819 {
820   // Print the list of digits belonging to the cluster
821   return;
822   TString message ; 
823   message  = "AliEMCALRecPoint:\n" ;
824   message +=  " digits # = " ; 
825   Info("Print", message.Data()) ; 
826
827   Int_t iDigit;
828   for(iDigit=0; iDigit<fMulDigit; iDigit++)
829     printf(" %d ", fDigitsList[iDigit] ) ;  
830   printf("\n");
831
832   Info("Print", " Energies = ") ;
833   for(iDigit=0; iDigit<fMulDigit; iDigit++) 
834     printf(" %f ", fEnergyList[iDigit] ) ;
835   printf("\n");
836
837   Info("Print", "\n Abs Ids = ") ;
838   for(iDigit=0; iDigit<fMulDigit; iDigit++) 
839     printf(" %i ", fAbsIdList[iDigit] ) ;
840   printf("\n");
841
842   Info("Print", " Primaries  ") ;
843   for(iDigit = 0;iDigit < fMulTrack; iDigit++)
844     printf(" %d ", fTracksList[iDigit]) ;
845
846   printf("\n Local x %6.2f y %7.2f z %7.1f \n", fLocPos[0], fLocPos[1], fLocPos[2]);
847
848   message  = "       ClusterType     = %d" ;
849   message += "       Multiplicity    = %d" ;
850   message += "       Cluster Energy  = %f" ; 
851   message += "       Core energy     = %f" ; 
852   message += "       Core radius     = %f" ; 
853   message += "       Number of primaries %d" ; 
854   message += "       Stored at position %d" ; 
855   Info("Print", message.Data(), fClusterType, fMulDigit, fAmp, fCoreEnergy, fCoreRadius, fMulTrack, GetIndexInList() ) ;  
856 }