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