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