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