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