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