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