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