]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALRecPoint.cxx
suppress off-line finders ( vertex etc.)
[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){ 
470     EvalLocal2TrackingCSTransform();
471   }
472
473 }
474
475 //____________________________________________________________________________
476 void  AliEMCALRecPoint::EvalDispersion(Float_t logWeight, TClonesArray * digits)
477 {
478   // Calculates the dispersion of the shower at the origin of the RecPoint
479   // in cell units - Nov 16,2006
480
481   Double_t d = 0., wtot = 0., w = 0.;
482   Int_t iDigit=0, nstat=0;
483   AliEMCALDigit * digit=0;
484         
485   // Calculates the dispersion in cell units 
486   Double_t etai, phii, etaMean=0.0, phiMean=0.0; 
487   int nSupMod=0, nModule=0, nIphi=0, nIeta=0;
488   int iphi=0, ieta=0;
489   // Calculate mean values
490   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
491     digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit])  ;
492
493     if (fAmp>0 && fEnergyList[iDigit]>0) {
494       fGeomPtr->GetCellIndex(digit->GetId(), nSupMod,nModule,nIphi,nIeta);
495       fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nModule,nIphi,nIeta, iphi,ieta);
496         
497       // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2
498       // C Side impair SM, nSupMod%2=1; A side pair SM nSupMod%2=0
499       if(fSharedCluster && nSupMod%2) ieta+=AliEMCALGeoParams::fgkEMCALCols;
500                 
501       etai=(Double_t)ieta;
502       phii=(Double_t)iphi;
503       w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
504
505       if(w>0.0) {
506         phiMean += phii*w;
507         etaMean += etai*w;
508         wtot    += w;
509       }
510     }
511   }
512   if (wtot>0) {
513     phiMean /= wtot ;
514     etaMean /= wtot ;
515   } else AliError(Form("Wrong weight %f\n", wtot));
516
517   // Calculate dispersion
518   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
519     digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit])  ;
520
521     if (fAmp>0 && fEnergyList[iDigit]>0) {
522       fGeomPtr->GetCellIndex(digit->GetId(), nSupMod,nModule,nIphi,nIeta);
523       fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nModule,nIphi,nIeta, iphi,ieta);
524                 
525       // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2
526       // C Side impair SM, nSupMod%2=1; A side pair SM, nSupMod%2=0
527       if(fSharedCluster && nSupMod%2) ieta+=AliEMCALGeoParams::fgkEMCALCols;
528       
529       etai=(Double_t)ieta;
530       phii=(Double_t)iphi;
531       w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
532
533       if(w>0.0) {
534         nstat++;
535         d += w*((etai-etaMean)*(etai-etaMean)+(phii-phiMean)*(phii-phiMean));
536       }
537     }
538   }
539   
540   if ( wtot > 0 && nstat>1) d /= wtot ;
541   else                      d = 0. ; 
542
543   fDispersion = TMath::Sqrt(d) ;
544   //printf("AliEMCALRecPoint::EvalDispersion() : Dispersion %f \n",fDispersion);
545 }
546
547 //____________________________________________________________________________
548 void AliEMCALRecPoint::EvalDistanceToBadChannels(AliCaloCalibPedestal* caloped)
549 {
550   //For each EMC rec. point set the distance to the nearest bad channel.
551   //AliInfo(Form("%d bad channel(s) found.\n", caloped->GetDeadTowerCount()));
552   //It is done in cell units and not in global or local position as before (Sept 2010)
553   
554   if(!caloped->GetDeadTowerCount()) return;
555   
556   //Get channels map of the supermodule where the cluster is.
557   TH2D* hMap  = caloped->GetDeadMap(fSuperModuleNumber);
558   
559   Int_t dRrow, dReta;   
560   Float_t  minDist = 10000.;
561   Float_t  dist    = 0.;
562   Int_t nSupMod, nModule;
563   Int_t nIphi, nIeta;
564   Int_t iphi, ieta;
565   fDigitIndMax  = GetMaximalEnergyIndex();
566   fGeomPtr->GetCellIndex(fAbsIdList[fDigitIndMax], nSupMod,nModule,nIphi,nIeta);
567   fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nModule,nIphi,nIeta, iphi,ieta);
568
569   //    TVector3 dR;    
570   //    TVector3 cellpos;
571   //    Float_t  minDist = 100000;
572   //    Float_t  dist    = 0;
573   //    Int_t    absId   = -1;  
574   
575   //Loop on tower status map 
576   for(Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++){
577     for(Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++){
578       //Check if tower is bad.
579       if(hMap->GetBinContent(icol,irow)==AliCaloCalibPedestal::kAlive) continue;
580       //printf("AliEMCALRecPoint::EvalDistanceToBadChannels() - Bad channel in SM %d, col %d, row %d\n",iSM,icol, irow);
581       
582       dRrow=TMath::Abs(irow-iphi);
583       dReta=TMath::Abs(icol-ieta);
584       dist=TMath::Sqrt(dRrow*dRrow+dReta*dReta);
585       if(dist < minDist) minDist = dist;
586       
587       //                        //Tower is bad, get the absId of the index.
588       //                        absId = fGeomPtr->GetAbsCellIdFromCellIndexes(fSuperModuleNumber, irow, icol); 
589       //                        
590       //                        //Get the position of this tower.
591       //                        
592       //                        //Calculate the distance in local coordinates
593       //                        //fGeomPtr->RelPosCellInSModule(absId,cellpos);
594       //                        //Calculate distance between this tower and cluster, set if is smaller than previous.
595       //                        //dR = cellpos-fLocPos;
596       //                        
597       //                        //Calculate the distance in global coordinates
598       //                        fGeomPtr->GetGlobal(absId,cellpos);
599       //                        //Calculate distance between this tower and cluster, set if it is smaller than previous.
600       //                        dR = cellpos-fGlobPos;
601       //                        
602       //                        dist = dR.Mag();
603       //                        if(dist < minDist) minDist = dist;
604       
605     }
606   }
607   
608   //In case the cluster is shared by 2 SuperModules, need to check the map of the second Super Module
609   if (fSharedCluster) {
610     TH2D* hMap2 = 0;
611     Int_t nSupMod2 = -1;
612     
613     //The only possible combinations are (0,1), (2,3) ... (10,11)
614     if(fSuperModuleNumber%2) nSupMod2 = fSuperModuleNumber-1;
615     else                     nSupMod2 = fSuperModuleNumber+1;
616     hMap2  = caloped->GetDeadMap(nSupMod2);
617     
618     //Loop on tower status map of second super module
619     for(Int_t irow = 0; irow < AliEMCALGeoParams::fgkEMCALRows; irow++){
620       for(Int_t icol = 0; icol < AliEMCALGeoParams::fgkEMCALCols; icol++){
621         //Check if tower is bad.
622         if(hMap2->GetBinContent(icol,irow)==AliCaloCalibPedestal::kAlive) continue;
623         //printf("AliEMCALRecPoint::EvalDistanceToBadChannels() - Bad channel in SM %d, col %d, row %d\n",iSM,icol, irow);
624         dRrow=TMath::Abs(irow-iphi);
625         
626         if(fSuperModuleNumber%2) {
627           dReta=TMath::Abs(icol-(AliEMCALGeoParams::fgkEMCALCols+ieta));
628        }
629        else {
630          dReta=TMath::Abs(AliEMCALGeoParams::fgkEMCALCols+icol-ieta);
631        }                    
632         
633        dist=TMath::Sqrt(dRrow*dRrow+dReta*dReta);
634         if(dist < minDist) minDist = dist;        
635         
636 //                              
637 //                              //Tower is bad, get the absId of the index.
638 //                              absId = fGeomPtr->GetAbsCellIdFromCellIndexes(nSupMod2, irow, icol); 
639 //                              
640 //                              //Get the position of this tower.
641 //                              
642 //                              //Calculate the distance in global coordinates
643 //                              fGeomPtr->GetGlobal(absId,cellpos);
644 //                              //Calculate distance between this tower and cluster, set if it is smaller than previous.
645 //                              dR = cellpos-fGlobPos;
646 //                              
647 //                              dist = dR.Mag();
648 //                              if(dist < minDist) minDist = dist;
649       }
650     }
651     
652   }// shared cluster in 2 SuperModules
653   
654   fDistToBadTower = minDist;
655   //printf("AliEMCALRecPoint::EvalDistanceToBadChannel() - Distance to Bad is %f cm, shared cluster? %d \n",fDistToBadTower,fSharedCluster);
656 }
657
658
659 //____________________________________________________________________________
660 void AliEMCALRecPoint::EvalLocalPosition(Float_t logWeight, TClonesArray * digits)
661 {
662   // Calculates the center of gravity in the local EMCAL-module coordinates 
663   //  Info("Print", " logWeight %f : cluster energy %f ", logWeight, fAmp); // for testing
664   
665   AliEMCALDigit * digit=0;
666   Int_t i=0, nstat=0;
667   
668   Double_t dist  = TmaxInCm(Double_t(fAmp));
669   //Int_t       idMax = GetAbsIdMaxDigit();// idMax is not used at all in RelPosCellInSModule, why use it?
670   
671   Double_t clXYZ[3]={0.,0.,0.}, clRmsXYZ[3]={0.,0.,0.}, xyzi[3], wtot=0., w=0.;
672   
673   //printf(" dist : %f  e : %f \n", dist, fAmp);
674   for(Int_t iDigit=0; iDigit<fMulDigit; iDigit++) {
675     digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
676     
677     if(!digit) {
678       AliError("No Digit!!");
679       continue;
680     }
681     
682     //fGeomPtr->RelPosCellInSModule(digit->GetId(), idMax, dist, xyzi[0], xyzi[1], xyzi[2]);
683     fGeomPtr->RelPosCellInSModule(digit->GetId(), dist, xyzi[0], xyzi[1], xyzi[2]);
684     
685     //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
686     if(fSharedCluster && fSuperModuleNumber != fGeomPtr->GetSuperModuleNumber(digit->GetId())) xyzi[1]*=-1;
687     
688     //printf("EvalLocalPosition Cell:  Id %i, SM %i : dist %f Local x,y,z %f %f %f \n", 
689     //          digit->GetId(), fGeomPtr->GetSuperModuleNumber(digit->GetId()), dist, xyzi[0], xyzi[1], xyzi[2]);
690     
691     if(logWeight > 0.0)  w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ));
692     else  w = fEnergyList[iDigit]; // just energy
693     
694     if(w>0.0) {
695       wtot += w ;
696       nstat++;
697       for(i=0; i<3; i++ ) {
698         clXYZ[i]    += (w*xyzi[i]);
699         clRmsXYZ[i] += (w*xyzi[i]*xyzi[i]);
700       }
701     }
702   }
703   //  cout << " wtot " << wtot << endl;
704   if ( wtot > 0 ) { 
705     //    xRMS   = TMath::Sqrt(x2m - xMean*xMean);
706     for(i=0; i<3; i++ ) {
707       clXYZ[i] /= wtot;
708       if(nstat>1) {
709         clRmsXYZ[i] /= (wtot*wtot);  
710         clRmsXYZ[i]  =  clRmsXYZ[i] - clXYZ[i]*clXYZ[i];
711         if(clRmsXYZ[i] > 0.0) {
712           clRmsXYZ[i] =  TMath::Sqrt(clRmsXYZ[i]);
713         } else      clRmsXYZ[i] = 0;
714       } else        clRmsXYZ[i] = 0;
715     }    
716   } else {
717     for(i=0; i<3; i++ ) {
718       clXYZ[i] = clRmsXYZ[i] = -1.;
719     }
720   }
721   // clRmsXYZ[i] ??
722   
723   //    // Cluster of one single digit, smear the position to avoid discrete position
724   //    // smear x and z with +- 3 cm to uniform (avoid discrete effects). Tower size is approx 6 cm.
725   //    // Rndm generates a number in ]0,1]
726   //    if (fMulDigit==1) { 
727   //      clXYZ[0] += fGeomPtr->GetPhiTileSize()*(0.5 - gRandom->Rndm()); 
728   //      clXYZ[2] += fGeomPtr->GetEtaTileSize()*(0.5 - gRandom->Rndm()); 
729   //    }
730   
731   //Set position in local vector
732   fLocPos.SetX(clXYZ[0]);
733   fLocPos.SetY(clXYZ[1]);
734   fLocPos.SetZ(clXYZ[2]);
735   
736   if (gDebug==2)
737     printf("EvalLocalPosition Cluster: Local (x,y,z) = (%f,%f,%f) \n", fLocPos.X(), fLocPos.Y(), fLocPos.Z()) ; 
738   
739 }
740
741
742 //____________________________________________________________________________
743 void AliEMCALRecPoint::EvalGlobalPosition(Float_t logWeight, TClonesArray * digits)
744 {
745   // Calculates the center of gravity in the global ALICE coordinates 
746   //  Info("Print", " logWeight %f : cluster energy %f ", logWeight, fAmp); // for testing
747   
748   AliEMCALDigit * digit=0;
749   Int_t i=0, nstat=0;
750         
751   Double_t dist  = TmaxInCm(Double_t(fAmp));
752   //Int_t       idMax = GetAbsIdMaxDigit();// idMax is not used at all in RelPosCellInSModule, why use it?
753         
754   Double_t clXYZ[3]={0.,0.,0.}, clRmsXYZ[3]={0.,0.,0.}, lxyzi[3], xyzi[3], wtot=0., w=0.;
755
756   //printf(" dist : %f  e : %f \n", dist, fAmp);
757   for(Int_t iDigit=0; iDigit<fMulDigit; iDigit++) {
758     digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
759
760     if(!digit) {
761       AliError("No Digit!!");
762       continue;
763     }    
764     
765     //Get the local coordinates of the cell
766     //fGeomPtr->RelPosCellInSModule(digit->GetId(), idMax, dist, lxyzi[0], lxyzi[1], lxyzi[2]);
767     fGeomPtr->RelPosCellInSModule(digit->GetId(), dist, lxyzi[0], lxyzi[1], lxyzi[2]);
768     
769     //Now get the global coordinate
770     fGeomPtr->GetGlobal(lxyzi,xyzi, fGeomPtr->GetSuperModuleNumber(digit->GetId()));
771     //TVector3 pos(xyzi[0], xyzi[1], xyzi[2]);
772     //printf("EvalGlobalPosition Cell:  Id %i, SM %i : dist %f Local (x,y,z) = (%f %f %f), eta %f, phi%f \n", 
773     //     digit->GetId(), fGeomPtr->GetSuperModuleNumber(digit->GetId()),dist, xyzi[0], xyzi[1], xyzi[2],pos.Eta(),pos.Phi()*TMath::RadToDeg());
774           
775     if(logWeight > 0.0)  w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ));
776     else  w = fEnergyList[iDigit]; // just energy
777
778     if(w>0.0) {
779       wtot += w ;
780       nstat++;
781       for(i=0; i<3; i++ ) {
782         clXYZ[i]    += (w*xyzi[i]);
783         clRmsXYZ[i] += (w*xyzi[i]*xyzi[i]);
784       }
785     }
786   }
787   //  cout << " wtot " << wtot << endl;
788   if ( wtot > 0 ) { 
789     //    xRMS   = TMath::Sqrt(x2m - xMean*xMean);
790     for(i=0; i<3; i++ ) {
791       clXYZ[i] /= wtot;
792       if(nstat>1) {
793         clRmsXYZ[i] /= (wtot*wtot);  
794         clRmsXYZ[i]  =  clRmsXYZ[i] - clXYZ[i]*clXYZ[i];
795         if(clRmsXYZ[i] > 0.0) {
796           clRmsXYZ[i] =  TMath::Sqrt(clRmsXYZ[i]);
797         } else      clRmsXYZ[i] = 0;
798       } else        clRmsXYZ[i] = 0;
799     }    
800   } else {
801     for(i=0; i<3; i++ ) {
802       clXYZ[i] = clRmsXYZ[i] = -1.;
803     }
804   }
805   // clRmsXYZ[i] ??
806
807 //  // Cluster of one single digit, smear the position to avoid discrete position
808 //  // smear x and z with +- 3 cm to uniform (avoid discrete effects). Tower size is approx 6 cm.
809 //  // Rndm generates a number in ]0,1]
810 //  if (fMulDigit==1) { 
811 //    clXYZ[0] += fGeomPtr->GetPhiTileSize()*(0.5 - gRandom->Rndm()); 
812 //    clXYZ[2] += fGeomPtr->GetEtaTileSize()*(0.5 - gRandom->Rndm());   
813 //  }
814         
815   //Set position in global vector
816   fGlobPos.SetX(clXYZ[0]);
817   fGlobPos.SetY(clXYZ[1]);
818   fGlobPos.SetZ(clXYZ[2]);
819                 
820   if (gDebug==2)
821         printf("EvalGlobalPosition Cluster: (x ,y ,z) = (%f,%f,%f), eta %f,phi %f\n", 
822                    fGlobPos.X(), fGlobPos.Y(), fGlobPos.Z(),fGlobPos.Eta(),fGlobPos.Phi()*TMath::RadToDeg()) ; 
823 }
824
825 //____________________________________________________________________________
826 void AliEMCALRecPoint::EvalLocalPositionFit(Double_t deff, Double_t logWeight, 
827 Double_t phiSlope, TClonesArray * digits)
828 {
829   // Evaluates local position of clusters in SM
830   
831   Double_t ycorr=0;
832   AliEMCALDigit *digit=0;
833   Int_t i=0, nstat=0;
834   Double_t clXYZ[3]={0.,0.,0.}, clRmsXYZ[3]={0.,0.,0.}, xyzi[3], wtot=0., w=0.; 
835
836   Double_t dist  = TmaxInCm(Double_t(fAmp));
837   //Int_t       idMax = GetAbsIdMaxDigit();// idMax is not used at all in RelPosCellInSModule, why use it?
838         
839   for(Int_t iDigit=0; iDigit<digits->GetEntries(); iDigit++) {
840     digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
841     if(digit){
842       dist = deff;
843       //fGeomPtr->RelPosCellInSModule(digit->GetId(), idMax, dist, xyzi[0], xyzi[1], xyzi[2]);
844       fGeomPtr->RelPosCellInSModule(digit->GetId(), dist, xyzi[0], xyzi[1], xyzi[2]);
845       
846       if(logWeight > 0.0)  w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ));
847       else                 w = fEnergyList[iDigit]; // just energy
848       
849       if(w>0.0) {
850         wtot += w ;
851         nstat++;
852         for(i=0; i<3; i++ ) {
853           clXYZ[i]    += (w*xyzi[i]);
854           clRmsXYZ[i] += (w*xyzi[i]*xyzi[i]);
855         }
856       }
857     }else AliError("Digit null");
858   }//loop
859   //  cout << " wtot " << wtot << endl;
860   if ( wtot > 0 ) { 
861     //    xRMS   = TMath::Sqrt(x2m - xMean*xMean);
862     for(i=0; i<3; i++ ) {
863       clXYZ[i] /= wtot;
864       if(nstat>1) {
865         clRmsXYZ[i] /= (wtot*wtot);  
866         clRmsXYZ[i]  =  clRmsXYZ[i] - clXYZ[i]*clXYZ[i];
867         if(clRmsXYZ[i] > 0.0) {
868           clRmsXYZ[i] =  TMath::Sqrt(clRmsXYZ[i]);
869         } else      clRmsXYZ[i] = 0;
870       } else        clRmsXYZ[i] = 0;
871     }    
872   } else {
873     for(i=0; i<3; i++ ) {
874       clXYZ[i] = clRmsXYZ[i] = -1.;
875     }
876   }
877   // clRmsXYZ[i] ??
878   if(phiSlope != 0.0 && logWeight > 0.0 && wtot) { 
879     // Correction in phi direction (y - coords here); Aug 16;
880     // May be put to global level or seperate method
881     ycorr = clXYZ[1] * (1. + phiSlope);
882     //printf(" y %f : ycorr %f : slope %f \n", clXYZ[1], ycorr, phiSlope); 
883     clXYZ[1] = ycorr;
884   }
885         
886   fLocPos.SetX(clXYZ[0]);
887   fLocPos.SetY(clXYZ[1]);
888   fLocPos.SetZ(clXYZ[2]);
889     
890 //  if (gDebug==2)
891 //    printf("EvalLocalPosition: eta,phi,r = %f,%f,%f", fLocPos.X(), fLocPos.Y(), fLocPos.Z()) ; 
892 }
893
894 //_____________________________________________________________________________
895 Bool_t AliEMCALRecPoint::EvalLocalPosition2(TClonesArray * digits, TArrayD &ed)
896 {
897   // Evaluated local position of rec.point using digits 
898   // and parametrisation of w0 and deff
899   //printf(" <I> AliEMCALRecPoint::EvalLocalPosition2() \n"); 
900   return AliEMCALRecPoint::EvalLocalPositionFromDigits(digits, ed, fLocPos);
901 }
902
903 //_____________________________________________________________________________
904 Bool_t AliEMCALRecPoint::EvalLocalPositionFromDigits(TClonesArray *digits, TArrayD &ed, TVector3 &locPos)
905 {
906   // Used when digits should be recalibrated
907   Double_t deff=0, w0=0, esum=0;
908   Int_t iDigit=0;
909   //  AliEMCALDigit *digit;
910
911   if(ed.GetSize() && (digits->GetEntries()!=ed.GetSize())) return kFALSE;
912
913   // Calculate sum energy of digits
914   esum = 0.0;
915   for(iDigit=0; iDigit<ed.GetSize(); iDigit++) esum += ed[iDigit];
916
917   GetDeffW0(esum, deff, w0);
918   
919   return EvalLocalPositionFromDigits(esum, deff, w0, digits, ed, locPos); 
920 }
921
922 //_____________________________________________________________________________
923 Bool_t AliEMCALRecPoint::EvalLocalPositionFromDigits(const Double_t esum, const Double_t deff, const Double_t w0, TClonesArray *digits, TArrayD &ed, TVector3 &locPos)
924 {
925   //Evaluate position of digits in supermodule.
926   AliEMCALDigit *digit=0;
927
928   Int_t i=0, nstat=0;
929   Double_t clXYZ[3]={0.,0.,0.}, xyzi[3], wtot=0., w=0.; 
930   //Int_t       idMax = GetAbsIdMaxDigit();// idMax is not used at all in RelPosCellInSModule, why use it?
931         
932   // Get pointer to EMCAL geometry
933   // (can't use fGeomPtr in static method)
934   AliEMCALGeometry* geo = AliEMCALGeometry::GetInstance(); 
935
936   for(Int_t iDigit=0; iDigit<digits->GetEntries(); iDigit++) {
937     digit = dynamic_cast<AliEMCALDigit *>(digits->At(iDigit));
938     if(digit){
939       //geo->RelPosCellInSModule(digit->GetId(), idMax, deff, xyzi[0], xyzi[1], xyzi[2]);
940       geo->RelPosCellInSModule(digit->GetId(), deff, xyzi[0], xyzi[1], xyzi[2]);
941       
942       if(w0 > 0.0)  w = TMath::Max( 0., w0 + TMath::Log(ed[iDigit] / esum));
943       else          w = ed[iDigit]; // just energy
944       
945       if(w>0.0) {
946         wtot += w ;
947         nstat++;
948         for(i=0; i<3; i++ ) {
949           clXYZ[i] += (w*xyzi[i]);
950         }
951       }
952     }else AliError("Digit null");
953   }//loop
954   //  cout << " wtot " << wtot << endl;
955   if (wtot > 0) { 
956     for(i=0; i<3; i++ ) {
957       clXYZ[i] /= wtot;
958     }
959     locPos.SetX(clXYZ[0]);
960     locPos.SetY(clXYZ[1]);
961     locPos.SetZ(clXYZ[2]);
962     return kTRUE;
963   } else {
964     return kFALSE;
965   }
966
967 }
968
969 //_____________________________________________________________________________
970 void AliEMCALRecPoint::GetDeffW0(const Double_t esum , Double_t &deff,  Double_t &w0)
971 {
972   //
973   // Aug 31, 2001 
974   // Applied for simulation data with threshold 3 adc
975   // Calculate efective distance (deff) and weigh parameter (w0) 
976   // for coordinate calculation; 0.5 GeV < esum <100 GeV.
977   // Look to:  http://rhic.physics.wayne.edu/~pavlinov/ALICE/SHISHKEBAB/RES/CALIB/GEOMCORR/deffandW0VaEgamma_2.gif
978   //
979   Double_t e=0.0;
980   const  Double_t kdp0=9.25147, kdp1=1.16700; // Hard coded now
981   const  Double_t kwp0=4.83713, kwp1=-2.77970e-01, kwp2 = 4.41116;
982
983   // No extrapolation here
984   e = esum<0.5?0.5:esum;
985   e = e>100.?100.:e;
986
987   deff = kdp0 + kdp1*TMath::Log(e);
988   w0   = kwp0 / (1. + TMath::Exp(kwp1*(e+kwp2)));
989   //printf("<I> AliEMCALRecPoint::GetDeffW0 esum %5.2f : deff %5.2f : w0 %5.2f \n", esum, deff, w0); 
990 }
991
992 //______________________________________________________________________________
993 void AliEMCALRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits)
994 {
995   // This function calculates energy in the core, 
996   // i.e. within a radius rad = fCoreEnergy around the center. Beyond this radius
997   // in accordance with shower profile the energy deposition 
998   // should be less than 2%
999   // Unfinished - Nov 15,2006
1000   // Distance is calculate in (phi,eta) units
1001
1002   AliEMCALDigit * digit = 0 ;
1003
1004   Int_t iDigit=0;
1005
1006   if (!fLocPos.Mag()) {
1007     EvalLocalPosition(logWeight, digits);
1008   }
1009   
1010   Double_t phiPoint = fLocPos.Phi(), etaPoint = fLocPos.Eta();
1011   Double_t eta, phi, distance;
1012   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
1013     digit = (AliEMCALDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
1014     
1015     eta = phi = 0.0;
1016     fGeomPtr->EtaPhiFromIndex(digit->GetId(),eta, phi) ;
1017     phi = phi * TMath::DegToRad();
1018   
1019     distance = TMath::Sqrt((eta-etaPoint)*(eta-etaPoint)+(phi-phiPoint)*(phi-phiPoint));
1020     if(distance < fCoreRadius)
1021       fCoreEnergy += fEnergyList[iDigit] ;
1022   }
1023   
1024 }
1025 //____________________________________________________________________________
1026 void  AliEMCALRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
1027 {
1028   // Calculates the axis of the shower ellipsoid in eta and phi
1029   // in cell units
1030
1031   TString gn(fGeomPtr->GetName());
1032
1033   Double_t wtot = 0.;
1034   Double_t x    = 0.;
1035   Double_t z    = 0.;
1036   Double_t dxx  = 0.;
1037   Double_t dzz  = 0.;
1038   Double_t dxz  = 0.;
1039
1040   AliEMCALDigit * digit = 0;
1041         
1042   Double_t etai =0, phii=0, w=0; 
1043   int nSupMod=0, nModule=0, nIphi=0, nIeta=0;
1044   int iphi=0, ieta=0;
1045   for(Int_t iDigit=0; iDigit<fMulDigit; iDigit++) {
1046     digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit])  ;
1047     etai = phii = 0.; 
1048     // Nov 15,2006 - use cell numbers as coordinates
1049     // Copied for shish-kebab geometry, ieta,iphi is cast as double as eta,phi
1050     // We can use the eta,phi(or coordinates) of cell
1051     nSupMod = nModule = nIphi = nIeta = iphi = ieta = 0;
1052
1053     fGeomPtr->GetCellIndex(digit->GetId(), nSupMod,nModule,nIphi,nIeta);
1054     fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nModule,nIphi,nIeta, iphi,ieta);
1055           
1056     // In case of a shared cluster, index of SM in C side, columns start at 48 and ends at 48*2
1057     // C Side impair SM, nSupMod%2=1; A side pair SM, nSupMod%2=0
1058     if(fSharedCluster && nSupMod%2) ieta+=AliEMCALGeoParams::fgkEMCALCols;
1059     
1060     etai=(Double_t)ieta;
1061     phii=(Double_t)iphi;
1062           
1063     w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
1064     // fAmp summed amplitude of digits, i.e. energy of recpoint 
1065     // Gives smaller value of lambda than log weight  
1066     // w = fEnergyList[iDigit] / fAmp; // Nov 16, 2006 - try just energy
1067
1068     dxx  += w * etai * etai ;
1069     x    += w * etai ;
1070     dzz  += w * phii * phii ;
1071     z    += w * phii ; 
1072
1073     dxz  += w * etai * phii ; 
1074
1075     wtot += w ;
1076   }
1077
1078   if ( wtot > 0 ) { 
1079     dxx /= wtot ;
1080     x   /= wtot ;
1081     dxx -= x * x ;
1082     dzz /= wtot ;
1083     z   /= wtot ;
1084     dzz -= z * z ;
1085     dxz /= wtot ;
1086     dxz -= x * z ;
1087
1088     fLambda[0] =  0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
1089     if(fLambda[0] > 0)
1090       fLambda[0] = TMath::Sqrt(fLambda[0]) ;
1091     else
1092       fLambda[0] = 0;
1093     
1094     fLambda[1] =  0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
1095
1096     if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
1097       fLambda[1] = TMath::Sqrt(fLambda[1]) ;
1098     else
1099       fLambda[1]= 0. ;
1100   } else { 
1101     fLambda[0]= 0. ;
1102     fLambda[1]= 0. ;
1103   }
1104
1105   //printf("AliEMCALRecPoint::EvalElipsAxis() lambdas  = %f,%f \n", fLambda[0],fLambda[1]) ; 
1106
1107 }
1108
1109 //______________________________________________________________________________
1110 void  AliEMCALRecPoint::EvalPrimaries(TClonesArray * digits)
1111 {
1112   // Constructs the list of primary particles (tracks) which 
1113   // have contributed to this RecPoint and calculate deposited energy 
1114   // for each track
1115     
1116   AliEMCALDigit * digit =0;
1117   Int_t * primArray = new Int_t[fMaxTrack] ;
1118   memset(primArray,-1,sizeof(Int_t)*fMaxTrack);
1119   Float_t * dEPrimArray = new Float_t[fMaxTrack] ;
1120   memset(dEPrimArray,-1,sizeof(Int_t)*fMaxTrack);
1121   
1122   Int_t index ;  
1123   for ( index = 0 ; index < GetDigitsMultiplicity() ; index++ ) { // all digits
1124     digit = dynamic_cast<AliEMCALDigit *>(digits->At( fDigitsList[index] )) ; 
1125     if(!digit) {
1126       AliError("No Digit!!");
1127       continue;
1128     }
1129     
1130     Int_t nprimaries = digit->GetNprimary() ;
1131     if ( nprimaries == 0 ) continue ;
1132     Int_t jndex ;
1133     for ( jndex = 0 ; jndex < nprimaries ; jndex++ ) { // all primaries in digit
1134       if ( fMulTrack > fMaxTrack ) {
1135         fMulTrack = fMaxTrack ;
1136         Error("EvalPrimaries", "increase fMaxTrack ")  ;
1137         break ;
1138       }
1139       Int_t newPrimary = digit->GetPrimary(jndex+1);
1140       Float_t dEPrimary = digit->GetDEPrimary(jndex+1);
1141       Int_t kndex ;
1142       Bool_t already = kFALSE ;
1143       for ( kndex = 0 ; kndex < fMulTrack ; kndex++ ) { //check if not already stored
1144         if ( newPrimary == primArray[kndex] ){
1145           already = kTRUE ;
1146           dEPrimArray[kndex] += dEPrimary; 
1147           break ;
1148         }
1149       } // end of check
1150       if ( !already && (fMulTrack < fMaxTrack)) { // store it
1151         primArray[fMulTrack] = newPrimary ; 
1152         dEPrimArray[fMulTrack] = dEPrimary ; 
1153         fMulTrack++ ;
1154       } // store it
1155     } // all primaries in digit
1156   } // all digits
1157   
1158   Int_t *sortIdx = new Int_t[fMulTrack];
1159   TMath::Sort(fMulTrack,dEPrimArray,sortIdx); 
1160   for(index = 0; index < fMulTrack; index++) {
1161     fTracksList[index] = primArray[sortIdx[index]] ;    
1162     fDETracksList[index] = dEPrimArray[sortIdx[index]] ;
1163   }
1164   delete [] sortIdx;
1165   delete [] primArray ;
1166   delete [] dEPrimArray ;
1167   
1168 }
1169
1170 //______________________________________________________________________________
1171 void  AliEMCALRecPoint::EvalParents(TClonesArray * digits)
1172 {
1173   // Constructs the list of parent particles (tracks) which have contributed to this RecPoint
1174
1175   AliEMCALDigit * digit=0 ;
1176   Int_t * parentArray = new Int_t[fMaxTrack] ;
1177   memset(parentArray,-1,sizeof(Int_t)*fMaxTrack);
1178   Float_t * dEParentArray = new Float_t[fMaxTrack] ;
1179   memset(dEParentArray,-1,sizeof(Int_t)*fMaxTrack);
1180   
1181   Int_t index ;  
1182   for ( index = 0 ; index < GetDigitsMultiplicity() ; index++ ) { // all digits
1183     if (fDigitsList[index] >= digits->GetEntries() || fDigitsList[index] < 0)
1184       AliError(Form("Trying to get invalid digit %d (idx in WriteRecPoint %d)",fDigitsList[index],index));
1185     digit = dynamic_cast<AliEMCALDigit *>(digits->At( fDigitsList[index] )) ; 
1186     if(!digit) {
1187       AliError("No Digit!!");
1188       continue;
1189     }
1190     
1191     Int_t nparents = digit->GetNiparent() ;
1192     if ( nparents == 0 ) continue ;
1193     
1194     Int_t jndex ;
1195     for ( jndex = 0 ; jndex < nparents ; jndex++ ) { // all primaries in digit
1196       if ( fMulParent > fMaxParent ) {
1197         fMulTrack = - 1 ;
1198         Error("EvalParents", "increase fMaxParent")  ;
1199         break ;
1200       }
1201       Int_t newParent = digit->GetIparent(jndex+1) ;
1202       Float_t newdEParent = digit->GetDEParent(jndex+1) ;
1203       Int_t kndex ;
1204       Bool_t already = kFALSE ;
1205       for ( kndex = 0 ; kndex < fMulParent ; kndex++ ) { //check if not already stored
1206         if ( newParent == parentArray[kndex] ){
1207           dEParentArray[kndex] += newdEParent;
1208           already = kTRUE ;
1209           break ;
1210         }
1211       } // end of check
1212       if ( !already && (fMulParent < fMaxParent)) { // store it
1213         parentArray[fMulParent] = newParent ; 
1214         dEParentArray[fMulParent] = newdEParent ; 
1215         fMulParent++ ;
1216       } // store it
1217     } // all parents in digit
1218   } // all digits
1219   
1220   if (fMulParent>0) {
1221     Int_t *sortIdx = new Int_t[fMulParent];
1222     TMath::Sort(fMulParent,dEParentArray,sortIdx); 
1223     for(index = 0; index < fMulParent; index++) {
1224       fParentsList[index] = parentArray[sortIdx[index]] ;      
1225       fDEParentsList[index] = dEParentArray[sortIdx[index]] ;
1226     }
1227     delete [] sortIdx;
1228   }
1229   
1230   delete [] parentArray;
1231   delete [] dEParentArray;
1232 }
1233
1234 //____________________________________________________________________________
1235 void AliEMCALRecPoint::GetLocalPosition(TVector3 & lpos) const
1236 {
1237   // returns the position of the cluster in the local reference system
1238   // of the sub-detector
1239   
1240   lpos = fLocPos;
1241 }
1242
1243 //____________________________________________________________________________
1244 void AliEMCALRecPoint::GetGlobalPosition(TVector3 & gpos) const
1245 {
1246   // returns the position of the cluster in the global reference system of ALICE
1247   // These are now the Cartesian X, Y and Z
1248   //  cout<<" geom "<<geom<<endl;
1249   // fGeomPtr->GetGlobal(fLocPos, gpos, fSuperModuleNumber);
1250   gpos = fGlobPos;
1251         
1252 }
1253
1254 //____________________________________________________________________________
1255 //void AliEMCALRecPoint::GetGlobalPosition(TVector3 & gpos, TMatrixF & gmat) const
1256 //{
1257 //  // returns the position of the cluster in the global reference system of ALICE
1258 //  // These are now the Cartesian X, Y and Z
1259 //  //  cout<<" geom "<<geom<<endl;
1260 //
1261 //  //To be implemented
1262 //  fGeomPtr->GetGlobalEMCAL(this, gpos, gmat);
1263 //
1264 //}
1265
1266 //_____________________________________________________________________________
1267 void AliEMCALRecPoint::EvalLocal2TrackingCSTransform()
1268 {
1269   //Evaluates local to "tracking" c.s. transformation (B.P.).
1270   //All evaluations should be completed before calling for this
1271   //function.                           
1272   //See ALICE PPR Chapter 5 p.18 for "tracking" c.s. definition,
1273   //or just ask Jouri Belikov. :) 
1274
1275   SetVolumeId(AliGeomManager::LayerToVolUID(AliGeomManager::kEMCAL,GetSuperModuleNumber()));
1276
1277   const TGeoHMatrix* tr2loc = GetTracking2LocalMatrix();
1278   if(!tr2loc) AliFatal(Form("No Tracking2LocalMatrix found."));
1279
1280   Double_t lxyz[3] = {fLocPos.X(),fLocPos.Y(),fLocPos.Z()};
1281   Double_t txyz[3] = {0,0,0};
1282
1283   tr2loc->MasterToLocal(lxyz,txyz);
1284   SetX(txyz[0]); SetY(txyz[1]); SetZ(txyz[2]);
1285
1286   if(AliLog::GetGlobalDebugLevel()>0) {
1287     TVector3 gpos; //TMatrixF gmat;
1288     //GetGlobalPosition(gpos,gmat); //Not doing anythin special, replace by next line.  
1289     fGeomPtr->GetGlobal(fLocPos, gpos, GetSuperModuleNumber());
1290     
1291     Float_t gxyz[3];
1292     GetGlobalXYZ(gxyz);
1293     AliInfo(Form("lCS-->(%.3f,%.3f,%.3f), tCS-->(%.3f,%.3f,%.3f), gCS-->(%.3f,%.3f,%.3f),  gCScalc-\
1294 ->(%.3f,%.3f,%.3f), supermodule %d",
1295                  fLocPos.X(),fLocPos.Y(),fLocPos.Z(),
1296                  GetX(),GetY(),GetZ(),
1297                  gpos.X(),gpos.Y(),gpos.Z(),
1298                  gxyz[0],gxyz[1],gxyz[2],GetSuperModuleNumber()));
1299   }
1300
1301 }
1302
1303 //____________________________________________________________________________
1304 Float_t AliEMCALRecPoint::GetMaximalEnergy(void) const
1305 {
1306   // Finds the maximum energy in the cluster
1307   
1308   Float_t menergy = 0. ;
1309
1310   Int_t iDigit;
1311   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
1312  
1313     if(fEnergyList[iDigit] > menergy) 
1314       menergy = fEnergyList[iDigit] ;
1315   }
1316   return menergy ;
1317 }
1318
1319 //____________________________________________________________________________
1320 Int_t AliEMCALRecPoint::GetMaximalEnergyIndex(void) const
1321 {
1322   // Finds the maximum energy in the cluster
1323   
1324   Float_t menergy = 0. ;
1325   Int_t mid       = 0  ;
1326   Int_t iDigit;
1327   
1328   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
1329     
1330     if(fEnergyList[iDigit] > menergy){ 
1331       menergy = fEnergyList[iDigit] ;
1332       mid     = iDigit ;
1333     }
1334   }//loop on cluster digits
1335   
1336   return mid ;
1337 }
1338
1339
1340 //____________________________________________________________________________
1341 Int_t AliEMCALRecPoint::GetMultiplicityAtLevel(Float_t H) const
1342 {
1343   // Calculates the multiplicity of digits with energy larger than H*energy 
1344   
1345   Int_t multipl   = 0 ;
1346   Int_t iDigit ;
1347   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
1348
1349     if(fEnergyList[iDigit] > H * fAmp) 
1350       multipl++ ;
1351   }
1352   return multipl ;
1353 }
1354
1355 //____________________________________________________________________________
1356 Int_t  AliEMCALRecPoint::GetNumberOfLocalMax(AliEMCALDigit **  maxAt, Float_t * maxAtEnergy,
1357                                                Float_t locMaxCut,TClonesArray * digits) const
1358
1359   // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
1360   // energy difference between two local maxima
1361
1362   AliEMCALDigit * digit  = 0;
1363   AliEMCALDigit * digitN = 0;
1364   
1365   Int_t iDigitN = 0 ;
1366   Int_t iDigit  = 0 ;
1367
1368   for(iDigit = 0; iDigit < fMulDigit; iDigit++)
1369     maxAt[iDigit] = (AliEMCALDigit*) digits->At(fDigitsList[iDigit])  ;
1370   
1371   for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {   
1372     if(maxAt[iDigit]) {
1373       digit = maxAt[iDigit] ;
1374           
1375       for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {        
1376         digitN = (AliEMCALDigit *) digits->At(fDigitsList[iDigitN]) ; 
1377         
1378         if ( AreNeighbours(digit, digitN) ) {
1379           if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {    
1380             maxAt[iDigitN] = 0 ;
1381             // but may be digit too is not local max ?
1382             if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut) 
1383               maxAt[iDigit] = 0 ;
1384           }
1385           else {
1386             maxAt[iDigit] = 0 ;
1387             // but may be digitN too is not local max ?
1388             if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut) 
1389               maxAt[iDigitN] = 0 ; 
1390           } 
1391         } // if Areneighbours
1392       } // while digitN
1393     } // slot not empty
1394   } // while digit
1395   
1396   iDigitN = 0 ;
1397   for(iDigit = 0; iDigit < fMulDigit; iDigit++) { 
1398     if(maxAt[iDigit] ){
1399       maxAt[iDigitN] = maxAt[iDigit] ;
1400       maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
1401       iDigitN++ ; 
1402     }
1403   }
1404   return iDigitN ;
1405 }
1406
1407 //____________________________________________________________________________
1408 Int_t AliEMCALRecPoint::GetPrimaryIndex() const  
1409 {
1410   // Get the primary track index in TreeK which deposits the most energy 
1411   // in Digits which forms RecPoint. 
1412
1413   if (fMulTrack)
1414     return fTracksList[0];
1415   return -12345;
1416 }
1417
1418 //____________________________________________________________________________
1419 void AliEMCALRecPoint::EvalTime(TClonesArray * digits){
1420   // time is set to the time of the digit with the maximum energy
1421
1422   Float_t maxE  = 0;
1423   Int_t   maxAt = 0;
1424   for(Int_t idig=0; idig < fMulDigit; idig++){
1425     if(fEnergyList[idig] > maxE){
1426       maxE  = fEnergyList[idig] ;
1427       maxAt = idig;
1428     }
1429   }
1430   fTime = ((AliEMCALDigit*) digits->At(fDigitsList[maxAt]))->GetTime() ;
1431   
1432 }
1433
1434 //______________________________________________________________________________
1435 void AliEMCALRecPoint::Paint(Option_t *)
1436 {
1437   // Paint this ALiRecPoint as a TMarker  with its current attributes
1438   
1439   TVector3 pos(0.,0.,0.)  ;
1440   GetLocalPosition(pos)   ;
1441   Coord_t x = pos.X()     ;
1442   Coord_t y = pos.Z()     ;
1443   Color_t markercolor = 1 ;
1444   Size_t  markersize  = 1.;
1445   Style_t markerstyle = 5 ;
1446   
1447   if (!gPad->IsBatch()) {
1448     gVirtualX->SetMarkerColor(markercolor) ;
1449     gVirtualX->SetMarkerSize (markersize)  ;
1450     gVirtualX->SetMarkerStyle(markerstyle) ;
1451   }
1452   gPad->SetAttMarkerPS(markercolor,markerstyle,markersize) ;
1453   gPad->PaintPolyMarker(1,&x,&y,"") ;
1454 }
1455
1456 //_____________________________________________________________________
1457 Double_t AliEMCALRecPoint::TmaxInCm(const Double_t e , const Int_t key)
1458
1459   // e energy in GeV)
1460   // key  =  0(gamma, default)
1461   //     !=  0(electron)
1462   const Double_t ca   = 4.82;  // shower max parameter - first guess; ca=TMath::Log(1000./8.07)
1463   Double_t tmax = 0.;    // position of electromagnetic shower max in cm
1464
1465   Double_t x0   = 1.31;  // radiation lenght (cm)
1466   //If old geometry in use
1467   if(!((fGeomPtr->GetEMCGeometry()->GetGeoName()).Contains("V1"))) x0 = 1.28;
1468
1469   if(e>0.1) {
1470     tmax = TMath::Log(e) + ca;
1471     if      (key==0) tmax += 0.5; 
1472     else             tmax -= 0.5;
1473     tmax *= x0; // convert to cm
1474   }
1475   return tmax;
1476 }
1477
1478 //______________________________________________________________________________
1479 Float_t AliEMCALRecPoint::EtaToTheta(Float_t arg) const
1480 {
1481   //Converts Theta (Radians) to Eta(Radians)
1482   return (2.*TMath::ATan(TMath::Exp(-arg)));
1483 }
1484
1485 //______________________________________________________________________________
1486 Float_t AliEMCALRecPoint::ThetaToEta(Float_t arg) const
1487 {
1488   //Converts Eta (Radians) to Theta(Radians)
1489   return (-1 * TMath::Log(TMath::Tan(0.5 * arg)));
1490 }
1491
1492 //____________________________________________________________________________
1493 void AliEMCALRecPoint::Print(Option_t *opt) const
1494 {
1495   // Print the list of digits belonging to the cluster
1496   if(strlen(opt)==0) return;
1497   TString message ; 
1498   message  = "AliEMCALRecPoint:\n" ;
1499   message +=  " digits # = " ; 
1500   AliInfo(message.Data()) ; 
1501
1502   Int_t iDigit;
1503   for(iDigit=0; iDigit<fMulDigit; iDigit++)
1504     printf(" %d ", fDigitsList[iDigit] ) ;  
1505   printf("\n");
1506
1507   AliInfo(" Energies = ") ;
1508   for(iDigit=0; iDigit<fMulDigit; iDigit++) 
1509     printf(" %f ", fEnergyList[iDigit] ) ;
1510   printf("\n");
1511
1512   AliInfo("\n Abs Ids = ") ;
1513   for(iDigit=0; iDigit<fMulDigit; iDigit++) 
1514     printf(" %i ", fAbsIdList[iDigit] ) ;
1515   printf("\n");
1516
1517   AliInfo(" Primaries  ") ;
1518   for(iDigit = 0;iDigit < fMulTrack; iDigit++)
1519     printf(" %d ", fTracksList[iDigit]) ;
1520
1521   printf("\n Local x %6.2f y %7.2f z %7.1f \n", fLocPos[0], fLocPos[1], fLocPos[2]);
1522
1523   message  = "       ClusterType     = %d" ;
1524   message += "       Multiplicity    = %d" ;
1525   message += "       Cluster Energy  = %f" ; 
1526   message += "       Core energy     = %f" ; 
1527   message += "       Core radius     = %f" ; 
1528   message += "       Number of primaries %d" ; 
1529   message += "       Stored at position %d" ; 
1530   AliInfo(Form(message.Data(), fClusterType, fMulDigit, fAmp, fCoreEnergy, fCoreRadius, fMulTrack, GetIndexInList()) ) ;  
1531 }
1532
1533 //___________________________________________________________
1534 Double_t  AliEMCALRecPoint::GetPointEnergy() const
1535 {
1536   //Returns energy ....
1537   Double_t e=0.0;
1538   for(int ic=0; ic<GetMultiplicity(); ic++) e += double(fEnergyList[ic]);
1539   return e;
1540 }