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