]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALRecPoint.cxx
Bug fix
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALRecPoint.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15 /* $Id$ */
16 //_________________________________________________________________________
17 //  Reconstructed Points for the EMCAL
18 //  A RecPoint is a cluster of digits
19 //*-- Author: Yves Schutz (SUBATECH)
20 //*-- Author: Dmitri Peressounko (RRC KI & SUBATECH)
21 //*-- Author: Heather Gray (LBL) merged AliEMCALRecPoint and AliEMCALTowerRecPoint 02/04
22
23 // --- ROOT system ---
24 #include "TPad.h"
25 #include "TGraph.h"
26 #include "TPaveText.h"
27 #include "TClonesArray.h"
28 #include "TMath.h"
29 #include "TGeoMatrix.h"
30 #include "TGeoManager.h"
31 #include "TGeoPhysicalNode.h"
32
33 // --- Standard library ---
34 #include <Riostream.h>
35
36 // --- AliRoot header files ---
37 //#include "AliGenerator.h"
38 class AliGenerator;
39 class AliEMCAL;
40 #include "AliLog.h"
41 #include "AliGeomManager.h"
42 #include "AliEMCALGeometry.h"
43 #include "AliEMCALHit.h"
44 #include "AliEMCALDigit.h"
45 #include "AliEMCALRecPoint.h"
46
47 ClassImp(AliEMCALRecPoint)
48
49 //____________________________________________________________________________
50 AliEMCALRecPoint::AliEMCALRecPoint()
51   : AliCluster(), fGeomPtr(0),
52     fAmp(0), fIndexInList(-1), //to be set when the point is already stored
53     fLocPos(0,0,0), fLocPosM(0),
54     fMaxDigit(100), fMulDigit(0), fMaxTrack(200),
55     fMulTrack(0), fDigitsList(0), fTracksList(0),
56     fClusterType(-1), fCoreEnergy(0), fDispersion(0),
57     fEnergyList(0), fTimeList(0), fAbsIdList(0),
58     fTime(0.), fCoreRadius(10),  //HG check this
59     fDETracksList(0), fMulParent(0), fMaxParent(0),
60     fParentsList(0), fDEParentsList(0), fSuperModuleNumber(0),
61     fDigitIndMax(-1)
62 {
63   // ctor
64   fGeomPtr = AliEMCALGeometry::GetInstance();
65   
66   fLambda[0] = 0;
67   fLambda[1] = 0;
68
69 }
70
71 //____________________________________________________________________________
72 AliEMCALRecPoint::AliEMCALRecPoint(const char *) 
73   : AliCluster(), fGeomPtr(0),
74     fAmp(0), fIndexInList(-1), //to be set when the point is already stored
75     fLocPos(0,0,0), fLocPosM(new TMatrixF(3,3)),
76     fMaxDigit(100), fMulDigit(0), fMaxTrack(1000), fMulTrack(0),
77     fDigitsList(new Int_t[fMaxDigit]), fTracksList(new Int_t[fMaxTrack]),
78     fClusterType(-1), fCoreEnergy(0), fDispersion(0),
79     fEnergyList(new Float_t[fMaxDigit]), fTimeList(new Float_t[fMaxDigit]), 
80     fAbsIdList(new Int_t[fMaxDigit]), fTime(-1.), fCoreRadius(10),
81     fDETracksList(new Float_t[fMaxTrack]), fMulParent(0), fMaxParent(1000),
82     fParentsList(new Int_t[fMaxParent]), fDEParentsList(new Float_t[fMaxParent]),
83     fSuperModuleNumber(0), fDigitIndMax(-1)
84 {
85   // ctor
86   for (Int_t i = 0; i < fMaxTrack; i++)
87     fDETracksList[i] = 0;
88   for (Int_t i = 0; i < fMaxParent; i++) {
89     fParentsList[i] = -1;
90     fDEParentsList[i] = 0;
91   }
92
93   fGeomPtr = AliEMCALGeometry::GetInstance();
94   fLambda[0] = 0;
95   fLambda[1] = 0;
96 }
97
98 //____________________________________________________________________________
99 AliEMCALRecPoint::AliEMCALRecPoint(const AliEMCALRecPoint & rp) 
100   : AliCluster(rp), fGeomPtr(rp.fGeomPtr),
101     fAmp(rp.fAmp), fIndexInList(rp.fIndexInList),
102     fLocPos(rp.fLocPos), fLocPosM(rp.fLocPosM),
103     fMaxDigit(rp.fMaxDigit), fMulDigit(rp.fMulDigit),
104     fMaxTrack(rp.fMaxTrack), fMulTrack(rp.fMaxTrack),
105     fDigitsList(new Int_t[rp.fMaxDigit]), fTracksList(new Int_t[rp.fMaxTrack]),
106     fClusterType(rp.fClusterType), fCoreEnergy(rp.fCoreEnergy), 
107     fDispersion(rp.fDispersion),
108     fEnergyList(new Float_t[rp.fMaxDigit]), fTimeList(new Float_t[rp.fMaxDigit]), 
109     fAbsIdList(new Int_t[rp.fMaxDigit]), fTime(rp.fTime), fCoreRadius(rp.fCoreRadius),
110     fDETracksList(new Float_t[rp.fMaxTrack]), fMulParent(rp.fMulParent), 
111     fMaxParent(rp.fMaxParent), fParentsList(new Int_t[rp.fMaxParent]), 
112     fDEParentsList(new Float_t[rp.fMaxParent]),
113     fSuperModuleNumber(rp.fSuperModuleNumber), fDigitIndMax(rp.fDigitIndMax)
114 {
115   //copy ctor
116   fLambda[0] = rp.fLambda[0];
117   fLambda[1] = rp.fLambda[1];
118
119   for(Int_t i = 0; i < rp.fMulDigit; i++) {
120     fEnergyList[i] = rp.fEnergyList[i];
121     fTimeList[i] = rp.fTimeList[i];
122     fAbsIdList[i] = rp.fAbsIdList[i];
123   }
124
125   for(Int_t i = 0; i < rp.fMulTrack; i++) fDETracksList[i] = rp.fDETracksList[i];
126
127   for(Int_t i = 0; i < rp.fMulParent; i++) {
128     fParentsList[i] = rp.fParentsList[i];
129     fDEParentsList[i] = rp.fDEParentsList[i];
130   }
131
132 }
133 //____________________________________________________________________________
134 AliEMCALRecPoint::~AliEMCALRecPoint()
135 {
136   // dtor
137   if ( fEnergyList )
138     delete[] fEnergyList ; 
139   if ( fTimeList )
140     delete[] fTimeList ; 
141   if ( fAbsIdList )
142     delete[] fAbsIdList ; 
143    if ( fDETracksList)
144     delete[] fDETracksList;
145    if ( fParentsList)
146     delete[] fParentsList;
147    if ( fDEParentsList)
148     delete[] fDEParentsList;
149
150    delete fLocPosM ;
151    delete [] fDigitsList ;
152    delete [] fTracksList ;
153 }
154
155 //____________________________________________________________________________
156 AliEMCALRecPoint& AliEMCALRecPoint::operator= (const AliEMCALRecPoint &rp)
157 {
158   // assignment operator
159
160   if(&rp == this) return *this;
161
162   fGeomPtr = rp.fGeomPtr;
163   fAmp = rp.fAmp;
164   fIndexInList = rp.fIndexInList;
165   fLocPos = rp.fLocPos;
166   fLocPosM = rp.fLocPosM;
167   fMaxDigit = rp.fMaxDigit;
168   fMulDigit = rp.fMulDigit;
169   fMaxTrack = rp.fMaxTrack;
170   fMulTrack = rp.fMaxTrack;
171   for(Int_t i = 0; i<fMaxDigit; i++) fDigitsList[i] = rp.fDigitsList[i];
172   for(Int_t i = 0; i<fMaxTrack; i++) fTracksList[i] = rp.fTracksList[i];
173   fClusterType = rp.fClusterType;
174   fCoreEnergy = rp.fCoreEnergy; 
175   fDispersion = rp.fDispersion;
176   for(Int_t i = 0; i<fMaxDigit; i++) {
177     fEnergyList[i] = rp.fEnergyList[i];
178     fTimeList[i] = rp.fTimeList[i]; 
179     fAbsIdList[i] = rp.fAbsIdList[i];
180   }
181   fTime = rp.fTime;
182   fCoreRadius = rp.fCoreRadius;
183   for(Int_t i = 0; i < fMaxTrack; i++) fDETracksList[i] = rp.fDETracksList[i];
184   fMulParent = rp.fMulParent;
185   fMaxParent = rp.fMaxParent;
186   for(Int_t i = 0; i < fMaxParent; i++) {
187     fParentsList[i] = rp.fParentsList[i]; 
188     fDEParentsList[i] = rp.fDEParentsList[i];
189   }
190   fSuperModuleNumber = rp.fSuperModuleNumber;
191   fDigitIndMax = rp.fDigitIndMax;
192
193   fLambda[0] = rp.fLambda[0];
194   fLambda[1] = rp.fLambda[1];
195
196   return *this;
197
198 }
199
200 //____________________________________________________________________________
201 void AliEMCALRecPoint::AddDigit(AliEMCALDigit & digit, Float_t Energy)
202 {
203   // Adds a digit to the RecPoint
204   // and accumulates the total amplitude and the multiplicity 
205   
206   if(fEnergyList == 0)
207     fEnergyList =  new Float_t[fMaxDigit]; 
208   if(fTimeList == 0)
209     fTimeList =  new Float_t[fMaxDigit]; 
210   if(fAbsIdList == 0) {
211     fAbsIdList =  new Int_t[fMaxDigit];
212   }
213
214   if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists 
215     fMaxDigit*=2 ; 
216     Int_t * tempo = new Int_t[fMaxDigit]; 
217     Float_t * tempoE =  new Float_t[fMaxDigit];
218     Float_t * tempoT =  new Float_t[fMaxDigit];
219     Int_t * tempoId = new Int_t[fMaxDigit]; 
220
221     Int_t index ;     
222     for ( index = 0 ; index < fMulDigit ; index++ ){
223       tempo[index]   = fDigitsList[index] ;
224       tempoE[index]  = fEnergyList[index] ; 
225       tempoT[index]  = fTimeList[index] ; 
226       tempoId[index] = fAbsIdList[index] ; 
227     }
228     
229     delete [] fDigitsList ;
230     delete [] fEnergyList ;
231     delete [] fTimeList ;
232     delete [] fAbsIdList ;
233
234     fDigitsList = tempo;
235     fEnergyList = tempoE; 
236     fTimeList = tempoT;
237     fAbsIdList = tempoId;
238   } // if
239   
240   fDigitsList[fMulDigit]   = digit.GetIndexInList()  ; 
241   fEnergyList[fMulDigit]   = Energy ;
242   fTimeList[fMulDigit]     = digit.GetTimeR() ;
243   fAbsIdList[fMulDigit]    = digit.GetId();
244   fMulDigit++ ; 
245   fAmp += Energy ; 
246
247   //JLK 10-Oct-2007 this hasn't been filled before because it was in
248   //the wrong place in previous versions.
249   //Now we evaluate it only if the supermodulenumber for this recpoint
250   //has not yet been set (or is the 0th one)
251   if(fSuperModuleNumber == 0)
252     fSuperModuleNumber = fGeomPtr->GetSuperModuleNumber(digit.GetId());
253
254 }
255 //____________________________________________________________________________
256 Bool_t AliEMCALRecPoint::AreNeighbours(AliEMCALDigit * digit1, AliEMCALDigit * digit2 ) const
257 {
258   // Tells if (true) or not (false) two digits are neighbours
259   // A neighbour is defined as being two digits which share a corner
260   
261   static Bool_t areNeighbours = kFALSE ;
262   static Int_t nSupMod=0, nModule=0, nIphi=0, nIeta=0;
263   static int nSupMod1=0, nModule1=0, nIphi1=0, nIeta1=0;
264   static Int_t relid1[2] , relid2[2] ; // ieta, iphi
265   static Int_t rowdiff=0, coldiff=0;
266
267   areNeighbours = kFALSE ;
268
269   fGeomPtr->GetCellIndex(digit1->GetId(), nSupMod,nModule,nIphi,nIeta);
270   fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nModule,nIphi,nIeta, relid1[0],relid1[1]);
271
272   fGeomPtr->GetCellIndex(digit2->GetId(), nSupMod1,nModule1,nIphi1,nIeta1);
273   fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod1,nModule1,nIphi1,nIeta1, relid2[0],relid2[1]);
274   
275   rowdiff = TMath::Abs( relid1[0] - relid2[0] ) ;  
276   coldiff = TMath::Abs( relid1[1] - relid2[1] ) ;  
277
278   if (( coldiff <= 1 )  && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0)) 
279   areNeighbours = kTRUE ;
280   
281   return areNeighbours;
282 }
283
284 //____________________________________________________________________________
285 Int_t AliEMCALRecPoint::Compare(const TObject * obj) const
286 {
287   // Compares two RecPoints according to their position in the EMCAL modules
288
289   Float_t delta = 1 ; //Width of "Sorting row". If you change this 
290                       //value (what is senseless) change as well delta in
291                       //AliEMCALTrackSegmentMakerv* and other RecPoints...
292   Int_t rv ; 
293
294   AliEMCALRecPoint * clu = (AliEMCALRecPoint *)obj ; 
295
296   TVector3 locpos1; 
297   GetLocalPosition(locpos1);
298   TVector3 locpos2;  
299   clu->GetLocalPosition(locpos2);  
300
301   Int_t rowdif = (Int_t)(TMath::Ceil(locpos1.X()/delta)-TMath::Ceil(locpos2.X()/delta)) ;
302   if (rowdif> 0) 
303     rv = 1 ;
304   else if(rowdif < 0) 
305     rv = -1 ;
306   else if(locpos1.Y()>locpos2.Y()) 
307     rv = -1 ;
308   else 
309     rv = 1 ; 
310
311   return rv ; 
312 }
313
314 //____________________________________________________________________________
315 Int_t AliEMCALRecPoint::DistancetoPrimitive(Int_t px, Int_t py)
316 {
317   // Compute distance from point px,py to  a AliEMCALRecPoint considered as a Tmarker
318   // Compute the closest distance of approach from point px,py to this marker.
319   // The distance is computed in pixels units.
320   // HG Still need to update -> Not sure what this should achieve
321
322   TVector3 pos(0.,0.,0.) ;
323   GetLocalPosition(pos) ;
324   Float_t x =  pos.X() ;
325   Float_t y =  pos.Y() ;
326   const Int_t kMaxDiff = 10;
327   Int_t pxm  = gPad->XtoAbsPixel(x);
328   Int_t pym  = gPad->YtoAbsPixel(y);
329   Int_t dist = (px-pxm)*(px-pxm) + (py-pym)*(py-pym);
330   
331   if (dist > kMaxDiff) return 9999;
332   return dist;
333 }
334
335 //___________________________________________________________________________
336  void AliEMCALRecPoint::Draw(Option_t *option)
337  {
338    // Draw this AliEMCALRecPoint with its current attributes
339    
340    AppendPad(option);
341  }
342
343 //______________________________________________________________________________
344 void AliEMCALRecPoint::ExecuteEvent(Int_t /*event*/, Int_t, Int_t)
345 {
346   // Execute action corresponding to one event
347   // This member function is called when a AliEMCALRecPoint is clicked with the locator
348   //
349   // If Left button is clicked on AliEMCALRecPoint, the digits are switched on    
350   // and switched off when the mouse button is released.
351
352   //  static Int_t pxold, pyold;
353
354   /*  static TGraph *  digitgraph = 0 ;
355   static TPaveText* clustertext = 0 ;
356   
357   if (!gPad->IsEditable()) return;
358   
359   switch (event) {
360     
361     
362   case kButton1Down:{
363     AliEMCALDigit * digit ;
364
365     Int_t iDigit;
366     Int_t relid[2] ;
367   
368     const Int_t kMulDigit=AliEMCALRecPoint::GetDigitsMultiplicity() ;
369     Float_t * xi = new Float_t [kMulDigit] ; 
370     Float_t * zi = new Float_t [kMulDigit] ;
371     
372     for(iDigit = 0; iDigit < kMulDigit; iDigit++) {
373       Fatal("AliEMCALRecPoint::ExecuteEvent", " -> Something wrong with the code"); 
374       digit = 0 ; //dynamic_cast<AliEMCALDigit *>((fDigitsList)[iDigit]);
375       fGeomPtr->AbsToRelNumbering(digit->GetId(), relid) ;
376       fGeomPtr->PosInAlice(relid, xi[iDigit], zi[iDigit]) ;
377     }
378     
379     if (!digitgraph) {
380       digitgraph = new TGraph(fMulDigit,xi,zi);
381       digitgraph-> SetMarkerStyle(5) ; 
382       digitgraph-> SetMarkerSize(1.) ;
383       digitgraph-> SetMarkerColor(1) ;
384       digitgraph-> Draw("P") ;
385     }
386     if (!clustertext) {
387       
388       TVector3 pos(0.,0.,0.) ;
389       GetLocalPosition(pos) ;
390       clustertext = new TPaveText(pos.X()-10,pos.Z()+10,pos.X()+50,pos.Z()+35,"") ;
391       Text_t  line1[40] ;
392       Text_t  line2[40] ;
393       sprintf(line1,"Energy=%1.2f GeV",GetEnergy()) ;
394       sprintf(line2,"%d Digits",GetDigitsMultiplicity()) ;
395       clustertext ->AddText(line1) ;
396       clustertext ->AddText(line2) ;
397       clustertext ->Draw("");
398     }
399     gPad->Update() ; 
400     Print("") ;
401     delete[] xi ; 
402     delete[] zi ; 
403    }
404   
405 break;
406   
407   case kButton1Up:
408     if (digitgraph) {
409       delete digitgraph  ;
410       digitgraph = 0 ;
411     }
412     if (clustertext) {
413       delete clustertext ;
414       clustertext = 0 ;
415     }
416     
417     break;
418     
419     }*/
420 }
421 //____________________________________________________________________________
422 void AliEMCALRecPoint::EvalAll(Float_t logWeight,TClonesArray * digits) 
423 {
424   // Evaluates all shower parameters
425   EvalLocalPosition(logWeight, digits) ;
426   EvalElipsAxis(logWeight, digits) ;
427   EvalDispersion(logWeight, digits) ;
428   //EvalCoreEnergy(logWeight, digits);
429   EvalTime(digits) ;
430   EvalPrimaries(digits) ;
431   EvalParents(digits);
432
433   //Called last because it sets the global position of the cluster?
434   EvalLocal2TrackingCSTransform();
435
436 }
437
438 //____________________________________________________________________________
439 void  AliEMCALRecPoint::EvalDispersion(Float_t logWeight, TClonesArray * digits)
440 {
441   // Calculates the dispersion of the shower at the origin of the RecPoint
442   // in cell units - Nov 16,2006
443
444   Double_t d = 0., wtot = 0., w = 0.;
445   Int_t iDigit=0, nstat=0;
446   AliEMCALDigit * digit ;
447   
448   // Calculates the dispersion in cell units 
449   Double_t etai, phii, etaMean=0.0, phiMean=0.0; 
450   int nSupMod=0, nModule=0, nIphi=0, nIeta=0;
451   int iphi=0, ieta=0;
452   // Calculate mean values
453   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
454     digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit])  ;
455
456     if (fAmp>0 && fEnergyList[iDigit]>0) {
457       fGeomPtr->GetCellIndex(digit->GetId(), nSupMod,nModule,nIphi,nIeta);
458       fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nModule,nIphi,nIeta, iphi,ieta);
459       etai=(Double_t)ieta;
460       phii=(Double_t)iphi;
461       w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
462
463       if(w>0.0) {
464         phiMean += phii*w;
465         etaMean += etai*w;
466         wtot    += w;
467       }
468     }
469   }
470   if (wtot>0) {
471     phiMean /= wtot ;
472     etaMean /= wtot ;
473   } else AliError(Form("Wrong weight %f\n", wtot));
474
475   // Calculate dispersion
476   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
477     digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit])  ;
478
479     if (fAmp>0 && fEnergyList[iDigit]>0) {
480       fGeomPtr->GetCellIndex(digit->GetId(), nSupMod,nModule,nIphi,nIeta);
481       fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nModule,nIphi,nIeta, iphi,ieta);
482       etai=(Double_t)ieta;
483       phii=(Double_t)iphi;
484       w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
485
486       if(w>0.0) {
487         nstat++;
488         d += w*((etai-etaMean)*(etai-etaMean)+(phii-phiMean)*(phii-phiMean));
489       }
490     }
491   }
492   
493   if ( wtot > 0 && nstat>1) d /= wtot ;
494   else                      d = 0. ; 
495
496   fDispersion = TMath::Sqrt(d) ;
497 }
498
499 //____________________________________________________________________________
500 void AliEMCALRecPoint::EvalLocalPosition(Float_t logWeight, TClonesArray * digits)
501 {
502   // Calculates the center of gravity in the local EMCAL-module coordinates 
503   //  Info("Print", " logWeight %f : cluster energy %f ", logWeight, fAmp); // for testing
504   
505   static Double_t dist;
506
507   AliEMCALDigit * digit;
508   Int_t i=0, nstat=0, idMax=-1;
509   Double_t clXYZ[3]={0.,0.,0.}, clRmsXYZ[3]={0.,0.,0.}, xyzi[3], wtot=0., w=0.;
510
511   //printf(" dist : %f  e : %f \n", dist, fAmp);
512   for(Int_t iDigit=0; iDigit<fMulDigit; iDigit++) {
513     digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
514     if(iDigit==0) {
515       idMax = digit->GetId(); // is it correct
516       dist  = TmaxInCm(Double_t(fAmp));
517     }
518     fGeomPtr->RelPosCellInSModule(digit->GetId(), idMax, dist, xyzi[0], xyzi[1], xyzi[2]);
519     //printf(" Id %i : dist %f Local x,y,z %f %f %f \n", digit->GetId(), dist, xyzi[0], xyzi[1], xyzi[2]);
520
521     //fGeomPtr->RelPosCellInSModule(digit->GetId(), xyzi[0], xyzi[1], xyzi[2]);
522     //printf(" Id %i : dist %f Local x,y,z %f %f %f \n", digit->GetId(), 0.0, xyzi[0], xyzi[1], xyzi[2]);
523     //    if(fAmp>102.) assert(0);
524
525     if(logWeight > 0.0)  w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ));
526     else  w = fEnergyList[iDigit]; // just energy
527
528     if(w>0.0) {
529       wtot += w ;
530       nstat++;
531       for(i=0; i<3; i++ ) {
532         clXYZ[i]    += (w*xyzi[i]);
533         clRmsXYZ[i] += (w*xyzi[i]*xyzi[i]);
534       }
535     }
536   }
537   //  cout << " wtot " << wtot << endl;
538   if ( wtot > 0 ) { 
539     //    xRMS   = TMath::Sqrt(x2m - xMean*xMean);
540     for(i=0; i<3; i++ ) {
541       clXYZ[i] /= wtot;
542       if(nstat>1) {
543         clRmsXYZ[i] /= (wtot*wtot);  
544         clRmsXYZ[i]  =  clRmsXYZ[i] - clXYZ[i]*clXYZ[i];
545         if(clRmsXYZ[i] > 0.0) {
546           clRmsXYZ[i] =  TMath::Sqrt(clRmsXYZ[i]);
547         } else      clRmsXYZ[i] = 0;
548       } else        clRmsXYZ[i] = 0;
549     }    
550   } else {
551     for(i=0; i<3; i++ ) {
552       clXYZ[i] = clRmsXYZ[i] = -1.;
553     }
554   }
555   // clRmsXYZ[i] ??
556   fLocPos.SetX(clXYZ[0]);
557   fLocPos.SetY(clXYZ[1]);
558   fLocPos.SetZ(clXYZ[2]);
559     
560 //  if (gDebug==2)
561 //    printf("EvalLocalPosition: eta,phi,r = %f,%f,%f", fLocPos.X(), fLocPos.Y(), fLocPos.Z()) ; 
562   fLocPosM = 0 ; // covariance matrix
563 }
564
565 //____________________________________________________________________________
566 void AliEMCALRecPoint::EvalLocalPositionFit(Double_t deff, Double_t logWeight, 
567 Double_t phiSlope, TClonesArray * digits)
568 {
569   // Aug 14-16, 2007 - for fit 
570   // Aug 31 - should be static ??
571   static Double_t dist, ycorr;
572   static AliEMCALDigit *digit;
573
574   Int_t i=0, nstat=0, idMax=-1;
575   Double_t clXYZ[3]={0.,0.,0.}, clRmsXYZ[3]={0.,0.,0.}, xyzi[3], wtot=0., w=0.; 
576
577   for(Int_t iDigit=0; iDigit<digits->GetEntries(); iDigit++) {
578     digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
579     if(iDigit==0) {
580       idMax = digit->GetId(); // is it correct
581       dist  = TmaxInCm(Double_t(fAmp));
582     }
583
584     dist = deff;
585     fGeomPtr->RelPosCellInSModule(digit->GetId(), idMax, dist, xyzi[0], xyzi[1], xyzi[2]);
586
587     if(logWeight > 0.0)  w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ));
588     else                 w = fEnergyList[iDigit]; // just energy
589
590     if(w>0.0) {
591       wtot += w ;
592       nstat++;
593       for(i=0; i<3; i++ ) {
594         clXYZ[i]    += (w*xyzi[i]);
595         clRmsXYZ[i] += (w*xyzi[i]*xyzi[i]);
596       }
597     }
598   }
599   //  cout << " wtot " << wtot << endl;
600   if ( wtot > 0 ) { 
601     //    xRMS   = TMath::Sqrt(x2m - xMean*xMean);
602     for(i=0; i<3; i++ ) {
603       clXYZ[i] /= wtot;
604       if(nstat>1) {
605         clRmsXYZ[i] /= (wtot*wtot);  
606         clRmsXYZ[i]  =  clRmsXYZ[i] - clXYZ[i]*clXYZ[i];
607         if(clRmsXYZ[i] > 0.0) {
608           clRmsXYZ[i] =  TMath::Sqrt(clRmsXYZ[i]);
609         } else      clRmsXYZ[i] = 0;
610       } else        clRmsXYZ[i] = 0;
611     }    
612   } else {
613     for(i=0; i<3; i++ ) {
614       clXYZ[i] = clRmsXYZ[i] = -1.;
615     }
616   }
617   // clRmsXYZ[i] ??
618   if(phiSlope != 0.0 && logWeight > 0.0 && wtot) { 
619     // Correction in phi direction (y - coords here); Aug 16;
620     // May be put to global level or seperate method
621     ycorr = clXYZ[1] * (1. + phiSlope);
622     //printf(" y %f : ycorr %f : slope %f \n", clXYZ[1], ycorr, phiSlope); 
623     clXYZ[1] = ycorr;
624   }
625   fLocPos.SetX(clXYZ[0]);
626   fLocPos.SetY(clXYZ[1]);
627   fLocPos.SetZ(clXYZ[2]);
628     
629 //  if (gDebug==2)
630 //    printf("EvalLocalPosition: eta,phi,r = %f,%f,%f", fLocPos.X(), fLocPos.Y(), fLocPos.Z()) ; 
631   fLocPosM = 0 ; // covariance matrix
632 }
633
634 //_____________________________________________________________________________
635 Bool_t AliEMCALRecPoint::EvalLocalPosition2(TClonesArray * digits, TArrayD &ed)
636 {
637   // Evaluated local position of rec.point using digits 
638   // and parametrisation of w0 and deff
639   //printf(" <I> AliEMCALRecPoint::EvalLocalPosition2() \n"); 
640   return AliEMCALRecPoint::EvalLocalPositionFromDigits(digits, ed, fLocPos);
641 }
642
643 //_____________________________________________________________________________
644 Bool_t AliEMCALRecPoint::EvalLocalPositionFromDigits(TClonesArray *digits, TArrayD &ed, TVector3 &locPos)
645 {
646   // Used when digits should be recalibrated
647   static Double_t deff, w0, esum;
648   static Int_t iDigit;
649   //  static AliEMCALDigit *digit;
650
651   if(ed.GetSize() && (digits->GetEntries()!=ed.GetSize())) return kFALSE;
652
653   // Calculate sum energy of digits
654   esum = 0.0;
655   for(iDigit=0; iDigit<ed.GetSize(); iDigit++) esum += ed[iDigit];
656
657   GetDeffW0(esum, deff, w0);
658   
659   return EvalLocalPositionFromDigits(esum, deff, w0, digits, ed, locPos); 
660 }
661
662 //_____________________________________________________________________________
663 Bool_t AliEMCALRecPoint::EvalLocalPositionFromDigits(const Double_t esum, const Double_t deff, const Double_t w0, TClonesArray *digits, TArrayD &ed, TVector3 &locPos)
664 {
665   static AliEMCALDigit *digit;
666
667   Int_t i=0, nstat=0, idMax=-1;
668   Double_t clXYZ[3]={0.,0.,0.}, xyzi[3], wtot=0., w=0.; 
669
670   // Get pointer to EMCAL geometry
671   // (can't use fGeomPtr in static method)
672   AliEMCALGeometry* geo = AliEMCALGeometry::GetInstance(); 
673
674   for(Int_t iDigit=0; iDigit<digits->GetEntries(); iDigit++) {
675     digit = dynamic_cast<AliEMCALDigit *>(digits->At(iDigit));
676
677     geo->RelPosCellInSModule(digit->GetId(), idMax, deff, xyzi[0], xyzi[1], xyzi[2]);
678
679     if(w0 > 0.0)  w = TMath::Max( 0., w0 + TMath::Log(ed[iDigit] / esum));
680     else          w = ed[iDigit]; // just energy
681
682     if(w>0.0) {
683       wtot += w ;
684       nstat++;
685       for(i=0; i<3; i++ ) {
686         clXYZ[i] += (w*xyzi[i]);
687       }
688     }
689   }
690   //  cout << " wtot " << wtot << endl;
691   if (wtot > 0) { 
692     for(i=0; i<3; i++ ) {
693       clXYZ[i] /= wtot;
694     }
695     locPos.SetX(clXYZ[0]);
696     locPos.SetY(clXYZ[1]);
697     locPos.SetZ(clXYZ[2]);
698     return kTRUE;
699   } else {
700     return kFALSE;
701   }
702
703 }
704
705 //_____________________________________________________________________________
706 void AliEMCALRecPoint::GetDeffW0(const Double_t esum , Double_t &deff,  Double_t &w0)
707 {
708   //
709   // Aug 31, 2001 
710   // Applied for simulation data with threshold 3 adc
711   // Calculate efective distance (deff) and weigh parameter (w0) 
712   // for coordinate calculation; 0.5 GeV < esum <100 GeV.
713   // Look to:  http://rhic.physics.wayne.edu/~pavlinov/ALICE/SHISHKEBAB/RES/CALIB/GEOMCORR/deffandW0VaEgamma_2.gif
714   //
715   static Double_t e=0.0;
716   const  Double_t kdp0=9.25147, kdp1=1.16700; // Hard coded now
717   const  Double_t kwp0=4.83713, kwp1=-2.77970e-01, kwp2 = 4.41116;
718
719   // No extrapolation here
720   e = esum<0.5?0.5:esum;
721   e = e>100.?100.:e;
722
723   deff = kdp0 + kdp1*TMath::Log(e);
724   w0   = kwp0 / (1. + TMath::Exp(kwp1*(e+kwp2)));
725   //printf("<I> AliEMCALRecPoint::GetDeffW0 esum %5.2f : deff %5.2f : w0 %5.2f \n", esum, deff, w0); 
726 }
727
728 //______________________________________________________________________________
729 void AliEMCALRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits)
730 {
731   // This function calculates energy in the core, 
732   // i.e. within a radius rad = fCoreEnergy around the center. Beyond this radius
733   // in accordance with shower profile the energy deposition 
734   // should be less than 2%
735   // Unfinished - Nov 15,2006
736   // Distance is calculate in (phi,eta) units
737
738   AliEMCALDigit * digit ;
739
740   Int_t iDigit;
741
742   if (!fLocPos.Mag()) {
743     EvalLocalPosition(logWeight, digits);
744   }
745   
746   Double_t phiPoint = fLocPos.Phi(), etaPoint = fLocPos.Eta();
747   Double_t eta, phi, distance;
748   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
749     digit = (AliEMCALDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
750     
751     eta = phi = 0.0;
752     fGeomPtr->EtaPhiFromIndex(digit->GetId(),eta, phi) ;
753     phi = phi * TMath::DegToRad();
754   
755     distance = TMath::Sqrt((eta-etaPoint)*(eta-etaPoint)+(phi-phiPoint)*(phi-phiPoint));
756     if(distance < fCoreRadius)
757       fCoreEnergy += fEnergyList[iDigit] ;
758   }
759   
760 }
761 //____________________________________________________________________________
762 void  AliEMCALRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
763 {
764   // Calculates the axis of the shower ellipsoid in eta and phi
765   // in cell units
766
767   static TString gn(fGeomPtr->GetName());
768
769   Double_t wtot = 0.;
770   Double_t x    = 0.;
771   Double_t z    = 0.;
772   Double_t dxx  = 0.;
773   Double_t dzz  = 0.;
774   Double_t dxz  = 0.;
775
776   AliEMCALDigit * digit = 0;
777
778   Double_t etai , phii, w; 
779   int nSupMod=0, nModule=0, nIphi=0, nIeta=0;
780   int iphi=0, ieta=0;
781   for(Int_t iDigit=0; iDigit<fMulDigit; iDigit++) {
782     digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit])  ;
783     etai = phii = 0.; 
784     if(gn.Contains("SHISH")) { 
785     // Nov 15,2006 - use cell numbers as coordinates
786     // Copied for shish-kebab geometry, ieta,iphi is cast as double as eta,phi
787     // We can use the eta,phi(or coordinates) of cell
788       nSupMod = nModule = nIphi = nIeta = iphi = ieta = 0;
789
790       fGeomPtr->GetCellIndex(digit->GetId(), nSupMod,nModule,nIphi,nIeta);
791       fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nModule,nIphi,nIeta, iphi,ieta);
792       etai=(Double_t)ieta;
793       phii=(Double_t)iphi;
794     } else { // 
795       fGeomPtr->EtaPhiFromIndex(digit->GetId(), etai, phii);
796       phii = phii * TMath::DegToRad();
797     }
798
799     w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
800     // fAmp summed amplitude of digits, i.e. energy of recpoint 
801     // Gives smaller value of lambda than log weight  
802     // w = fEnergyList[iDigit] / fAmp; // Nov 16, 2006 - try just energy
803
804     dxx  += w * etai * etai ;
805     x    += w * etai ;
806     dzz  += w * phii * phii ;
807     z    += w * phii ; 
808
809     dxz  += w * etai * phii ; 
810
811     wtot += w ;
812   }
813
814   if ( wtot > 0 ) { 
815     dxx /= wtot ;
816     x   /= wtot ;
817     dxx -= x * x ;
818     dzz /= wtot ;
819     z   /= wtot ;
820     dzz -= z * z ;
821     dxz /= wtot ;
822     dxz -= x * z ;
823
824     fLambda[0] =  0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
825     if(fLambda[0] > 0)
826       fLambda[0] = TMath::Sqrt(fLambda[0]) ;
827     else
828       fLambda[0] = 0;
829     
830     fLambda[1] =  0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
831
832     if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
833       fLambda[1] = TMath::Sqrt(fLambda[1]) ;
834     else
835       fLambda[1]= 0. ;
836   } else { 
837     fLambda[0]= 0. ;
838     fLambda[1]= 0. ;
839   }
840
841   //  printf("Evalaxis: lambdas  = %f,%f", fLambda[0],fLambda[1]) ; 
842
843 }
844
845 //______________________________________________________________________________
846 void  AliEMCALRecPoint::EvalPrimaries(TClonesArray * digits)
847 {
848   // Constructs the list of primary particles (tracks) which 
849   // have contributed to this RecPoint and calculate deposited energy 
850   // for each track
851   
852   AliEMCALDigit * digit ;
853   Int_t * primArray = new Int_t[fMaxTrack] ;
854   Float_t * dEPrimArray = new Float_t[fMaxTrack] ;
855
856   Int_t index ;  
857   for ( index = 0 ; index < GetDigitsMultiplicity() ; index++ ) { // all digits
858     digit = dynamic_cast<AliEMCALDigit *>(digits->At( fDigitsList[index] )) ; 
859     Int_t nprimaries = digit->GetNprimary() ;
860     if ( nprimaries == 0 ) continue ;
861     Int_t jndex ;
862     for ( jndex = 0 ; jndex < nprimaries ; jndex++ ) { // all primaries in digit
863       if ( fMulTrack > fMaxTrack ) {
864         fMulTrack = fMaxTrack ;
865         Error("EvalPrimaries", "increase fMaxTrack ")  ;
866         break ;
867       }
868       Int_t newPrimary = digit->GetPrimary(jndex+1);
869       Float_t dEPrimary = digit->GetDEPrimary(jndex+1);
870       Int_t kndex ;
871       Bool_t already = kFALSE ;
872       for ( kndex = 0 ; kndex < fMulTrack ; kndex++ ) { //check if not already stored
873         if ( newPrimary == primArray[kndex] ){
874           already = kTRUE ;
875           dEPrimArray[kndex] += dEPrimary; 
876           break ;
877         }
878       } // end of check
879       if ( !already && (fMulTrack < fMaxTrack)) { // store it
880         primArray[fMulTrack] = newPrimary ; 
881         dEPrimArray[fMulTrack] = dEPrimary ; 
882         fMulTrack++ ;
883       } // store it
884     } // all primaries in digit
885   } // all digits
886
887   Int_t *sortIdx = new Int_t[fMulTrack];
888   TMath::Sort(fMulTrack,dEPrimArray,sortIdx); 
889   for(index = 0; index < fMulTrack; index++) {
890     fTracksList[index] = primArray[sortIdx[index]] ;    
891     fDETracksList[index] = dEPrimArray[sortIdx[index]] ;
892   }
893   delete [] sortIdx;
894   delete [] primArray ;
895   delete [] dEPrimArray ;
896
897 }
898
899 //______________________________________________________________________________
900 void  AliEMCALRecPoint::EvalParents(TClonesArray * digits)
901 {
902   // Constructs the list of parent particles (tracks) which have contributed to this RecPoint
903  
904   AliEMCALDigit * digit ;
905   Int_t * parentArray = new Int_t[fMaxTrack] ;
906   Float_t * dEParentArray = new Float_t[fMaxTrack] ;
907
908   Int_t index ;  
909   for ( index = 0 ; index < GetDigitsMultiplicity() ; index++ ) { // all digits
910     if (fDigitsList[index] >= digits->GetEntries() || fDigitsList[index] < 0)
911        AliError(Form("Trying to get invalid digit %d (idx in WriteRecPoint %d)",fDigitsList[index],index));
912     digit = dynamic_cast<AliEMCALDigit *>(digits->At( fDigitsList[index] )) ; 
913     Int_t nparents = digit->GetNiparent() ;
914     if ( nparents == 0 ) continue ;
915
916     Int_t jndex ;
917     for ( jndex = 0 ; jndex < nparents ; jndex++ ) { // all primaries in digit
918       if ( fMulParent > fMaxParent ) {
919         fMulTrack = - 1 ;
920         Error("EvalParents", "increase fMaxParent")  ;
921         break ;
922       }
923       Int_t newParent = digit->GetIparent(jndex+1) ;
924       Float_t newdEParent = digit->GetDEParent(jndex+1) ;
925       Int_t kndex ;
926       Bool_t already = kFALSE ;
927       for ( kndex = 0 ; kndex < fMulParent ; kndex++ ) { //check if not already stored
928         if ( newParent == parentArray[kndex] ){
929           dEParentArray[kndex] += newdEParent;
930           already = kTRUE ;
931           break ;
932         }
933       } // end of check
934       if ( !already && (fMulParent < fMaxParent)) { // store it
935         parentArray[fMulParent] = newParent ; 
936         dEParentArray[fMulParent] = newdEParent ; 
937         fMulParent++ ;
938       } // store it
939     } // all parents in digit
940   } // all digits
941
942   if (fMulParent>0) {
943     Int_t *sortIdx = new Int_t[fMulParent];
944     TMath::Sort(fMulParent,dEParentArray,sortIdx); 
945     for(index = 0; index < fMulParent; index++) {
946       fParentsList[index] = parentArray[sortIdx[index]] ;      
947       fDEParentsList[index] = dEParentArray[sortIdx[index]] ;
948     }
949     delete [] sortIdx;
950   }
951  
952   delete [] parentArray;
953   delete [] dEParentArray;
954 }
955
956 //____________________________________________________________________________
957 void AliEMCALRecPoint::GetLocalPosition(TVector3 & lpos) const
958 {
959   // returns the position of the cluster in the local reference system
960   // of the sub-detector
961   
962   lpos = fLocPos;
963 }
964
965 //____________________________________________________________________________
966 void AliEMCALRecPoint::GetGlobalPosition(TVector3 & gpos) const
967 {
968   // returns the position of the cluster in the global reference system of ALICE
969   // These are now the Cartesian X, Y and Z
970   //  cout<<" geom "<<geom<<endl;
971   fGeomPtr->GetGlobal(fLocPos, gpos, fSuperModuleNumber);
972
973 }
974
975 //____________________________________________________________________________
976 void AliEMCALRecPoint::GetGlobalPosition(TVector3 & gpos, TMatrixF & gmat) const
977 {
978   // returns the position of the cluster in the global reference system of ALICE
979   // These are now the Cartesian X, Y and Z
980   //  cout<<" geom "<<geom<<endl;
981
982   //To be implemented
983   fGeomPtr->GetGlobalEMCAL(this, gpos, gmat);
984
985 }
986
987 //_____________________________________________________________________________
988 void AliEMCALRecPoint::EvalLocal2TrackingCSTransform()
989 {
990   //Evaluates local to "tracking" c.s. transformation (B.P.).
991   //All evaluations should be completed before calling for this
992   //function.                           
993   //See ALICE PPR Chapter 5 p.18 for "tracking" c.s. definition,
994   //or just ask Jouri Belikov. :) 
995
996   SetVolumeId(AliGeomManager::LayerToVolUID(AliGeomManager::kEMCAL,GetSuperModuleNumber()));
997
998   const TGeoHMatrix* tr2loc = GetTracking2LocalMatrix();
999   if(!tr2loc) AliFatal(Form("No Tracking2LocalMatrix found."));
1000
1001   Double_t lxyz[3] = {fLocPos.X(),fLocPos.Y(),fLocPos.Z()};
1002   Double_t txyz[3] = {0,0,0};
1003
1004   tr2loc->MasterToLocal(lxyz,txyz);
1005   SetX(txyz[0]); SetY(txyz[1]); SetZ(txyz[2]);
1006
1007   if(AliLog::GetGlobalDebugLevel()>0) {
1008     TVector3 gpos; TMatrixF gmat;
1009     GetGlobalPosition(gpos,gmat);
1010     Float_t gxyz[3];
1011     GetGlobalXYZ(gxyz);
1012     AliInfo(Form("lCS-->(%.3f,%.3f,%.3f), tCS-->(%.3f,%.3f,%.3f), gCS-->(%.3f,%.3f,%.3f),  gCScalc-\
1013 ->(%.3f,%.3f,%.3f), supermodule %d",
1014                  fLocPos.X(),fLocPos.Y(),fLocPos.Z(),
1015                  GetX(),GetY(),GetZ(),
1016                  gpos.X(),gpos.Y(),gpos.Z(),
1017                  gxyz[0],gxyz[1],gxyz[2],GetSuperModuleNumber()));
1018   }
1019
1020 }
1021
1022 //____________________________________________________________________________
1023 Float_t AliEMCALRecPoint::GetMaximalEnergy(void) const
1024 {
1025   // Finds the maximum energy in the cluster
1026   
1027   Float_t menergy = 0. ;
1028
1029   Int_t iDigit;
1030
1031   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
1032  
1033     if(fEnergyList[iDigit] > menergy) 
1034       menergy = fEnergyList[iDigit] ;
1035   }
1036   return menergy ;
1037 }
1038
1039 //____________________________________________________________________________
1040 Int_t AliEMCALRecPoint::GetMultiplicityAtLevel(Float_t H) const
1041 {
1042   // Calculates the multiplicity of digits with energy larger than H*energy 
1043   
1044   Int_t multipl   = 0 ;
1045   Int_t iDigit ;
1046   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
1047
1048     if(fEnergyList[iDigit] > H * fAmp) 
1049       multipl++ ;
1050   }
1051   return multipl ;
1052 }
1053
1054 //____________________________________________________________________________
1055 Int_t  AliEMCALRecPoint::GetNumberOfLocalMax(AliEMCALDigit **  maxAt, Float_t * maxAtEnergy,
1056                                                Float_t locMaxCut,TClonesArray * digits) const
1057
1058   // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
1059   // energy difference between two local maxima
1060
1061   AliEMCALDigit * digit ;
1062   AliEMCALDigit * digitN ;
1063   
1064   Int_t iDigitN ;
1065   Int_t iDigit ;
1066
1067   for(iDigit = 0; iDigit < fMulDigit; iDigit++)
1068     maxAt[iDigit] = (AliEMCALDigit*) digits->At(fDigitsList[iDigit])  ;
1069   
1070   for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {   
1071     if(maxAt[iDigit]) {
1072       digit = maxAt[iDigit] ;
1073           
1074       for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {        
1075         digitN = (AliEMCALDigit *) digits->At(fDigitsList[iDigitN]) ; 
1076         
1077         if ( AreNeighbours(digit, digitN) ) {
1078           if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {    
1079             maxAt[iDigitN] = 0 ;
1080             // but may be digit too is not local max ?
1081             if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut) 
1082               maxAt[iDigit] = 0 ;
1083           }
1084           else {
1085             maxAt[iDigit] = 0 ;
1086             // but may be digitN too is not local max ?
1087             if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut) 
1088               maxAt[iDigitN] = 0 ; 
1089           } 
1090         } // if Areneighbours
1091       } // while digitN
1092     } // slot not empty
1093   } // while digit
1094   
1095   iDigitN = 0 ;
1096   for(iDigit = 0; iDigit < fMulDigit; iDigit++) { 
1097     if(maxAt[iDigit] ){
1098       maxAt[iDigitN] = maxAt[iDigit] ;
1099       maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
1100       iDigitN++ ; 
1101     }
1102   }
1103   return iDigitN ;
1104 }
1105
1106 //____________________________________________________________________________
1107 Int_t AliEMCALRecPoint::GetPrimaryIndex() const  
1108 {
1109   // Get the primary track index in TreeK which deposits the most energy 
1110   // in Digits which forms RecPoint. 
1111
1112   if (fMulTrack)
1113     return fTracksList[0];
1114   return -12345;
1115 }
1116
1117 //____________________________________________________________________________
1118 void AliEMCALRecPoint::EvalTime(TClonesArray * digits){
1119   // time is set to the time of the digit with the maximum energy
1120
1121   Float_t maxE = 0;
1122   Int_t maxAt = 0;
1123   for(Int_t idig=0; idig < fMulDigit; idig++){
1124     if(fEnergyList[idig] > maxE){
1125       maxE = fEnergyList[idig] ;
1126       maxAt = idig;
1127     }
1128   }
1129   fTime = ((AliEMCALDigit*) digits->At(fDigitsList[maxAt]))->GetTime() ;
1130   
1131 }
1132
1133 //______________________________________________________________________________
1134 void AliEMCALRecPoint::Paint(Option_t *)
1135 {
1136   // Paint this ALiRecPoint as a TMarker  with its current attributes
1137   
1138   TVector3 pos(0.,0.,0.)  ;
1139   GetLocalPosition(pos)   ;
1140   Coord_t x = pos.X()     ;
1141   Coord_t y = pos.Z()     ;
1142   Color_t markercolor = 1 ;
1143   Size_t  markersize = 1. ;
1144   Style_t markerstyle = 5 ;
1145   
1146   if (!gPad->IsBatch()) {
1147     gVirtualX->SetMarkerColor(markercolor) ;
1148     gVirtualX->SetMarkerSize (markersize)  ;
1149     gVirtualX->SetMarkerStyle(markerstyle) ;
1150   }
1151   gPad->SetAttMarkerPS(markercolor,markerstyle,markersize) ;
1152   gPad->PaintPolyMarker(1,&x,&y,"") ;
1153 }
1154
1155 //_____________________________________________________________________
1156 Double_t AliEMCALRecPoint::TmaxInCm(const Double_t e , const Int_t key)
1157
1158   // e energy in GeV)
1159   // key  =  0(gamma, default)
1160   //     !=  0(electron)
1161   static Double_t ca = 4.82;  // shower max parameter - first guess; ca=TMath::Log(1000./8.07)
1162   static Double_t x0 = 1.23;  // radiation lenght (cm)
1163   static Double_t tmax = 0.;   // position of electromagnetic shower max in cm
1164
1165   tmax = 0.0;
1166   if(e>0.1) {
1167     tmax = TMath::Log(e) + ca;
1168     if      (key==0) tmax += 0.5; 
1169     else             tmax -= 0.5;
1170     tmax *= x0; // convert to cm
1171   }
1172   return tmax;
1173 }
1174
1175 //______________________________________________________________________________
1176 Float_t AliEMCALRecPoint::EtaToTheta(Float_t arg) const
1177 {
1178   //Converts Theta (Radians) to Eta(Radians)
1179   return (2.*TMath::ATan(TMath::Exp(-arg)));
1180 }
1181
1182 //______________________________________________________________________________
1183 Float_t AliEMCALRecPoint::ThetaToEta(Float_t arg) const
1184 {
1185   //Converts Eta (Radians) to Theta(Radians)
1186   return (-1 * TMath::Log(TMath::Tan(0.5 * arg)));
1187 }
1188
1189 //____________________________________________________________________________
1190 void AliEMCALRecPoint::Print(Option_t *opt) const
1191 {
1192   // Print the list of digits belonging to the cluster
1193   if(strlen(opt)==0) return;
1194   TString message ; 
1195   message  = "AliEMCALRecPoint:\n" ;
1196   message +=  " digits # = " ; 
1197   Info("Print", message.Data()) ; 
1198
1199   Int_t iDigit;
1200   for(iDigit=0; iDigit<fMulDigit; iDigit++)
1201     printf(" %d ", fDigitsList[iDigit] ) ;  
1202   printf("\n");
1203
1204   Info("Print", " Energies = ") ;
1205   for(iDigit=0; iDigit<fMulDigit; iDigit++) 
1206     printf(" %f ", fEnergyList[iDigit] ) ;
1207   printf("\n");
1208
1209   Info("Print", "\n Abs Ids = ") ;
1210   for(iDigit=0; iDigit<fMulDigit; iDigit++) 
1211     printf(" %i ", fAbsIdList[iDigit] ) ;
1212   printf("\n");
1213
1214   Info("Print", " Primaries  ") ;
1215   for(iDigit = 0;iDigit < fMulTrack; iDigit++)
1216     printf(" %d ", fTracksList[iDigit]) ;
1217
1218   printf("\n Local x %6.2f y %7.2f z %7.1f \n", fLocPos[0], fLocPos[1], fLocPos[2]);
1219
1220   message  = "       ClusterType     = %d" ;
1221   message += "       Multiplicity    = %d" ;
1222   message += "       Cluster Energy  = %f" ; 
1223   message += "       Core energy     = %f" ; 
1224   message += "       Core radius     = %f" ; 
1225   message += "       Number of primaries %d" ; 
1226   message += "       Stored at position %d" ; 
1227   Info("Print", message.Data(), fClusterType, fMulDigit, fAmp, fCoreEnergy, fCoreRadius, fMulTrack, GetIndexInList() ) ;  
1228 }
1229
1230 //___________________________________________________________
1231 Double_t  AliEMCALRecPoint::GetPointEnergy() const
1232 {
1233   //Returns energy ....
1234   static double e;
1235   e=0.0;
1236   for(int ic=0; ic<GetMultiplicity(); ic++) e += double(fEnergyList[ic]);
1237   return e;
1238 }