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