]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALRecPoint.cxx
updates for Effective C++ compiler flags
[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 //*-- Author: Yves Schutz (SUBATECH)
20 //*-- Author: Dmitri Peressounko (RRC KI & SUBATECH)
21 //*-- Author: Heather Gray (LBL) merged AliEMCALRecPoint and AliEMCALTowerRecPoint 02/04
22
23 // --- ROOT system ---
24 class Riostream;
25 #include <TPad.h>
26 class TGraph;
27 class TPaveText;
28 #include <TClonesArray.h>
29 #include <TMath.h>
30
31 // --- Standard library ---
32
33 // --- AliRoot header files ---
34 //#include "AliGenerator.h"
35 class AliGenerator;
36 #include "AliRunLoader.h"
37 #include "AliRun.h"
38 class AliEMCAL;
39 #include "AliEMCALLoader.h"
40 #include "AliEMCALGeometry.h"
41 #include "AliEMCALHit.h"
42 #include "AliEMCALDigit.h"
43 #include "AliEMCALRecPoint.h"
44
45 ClassImp(AliEMCALRecPoint)
46
47 //____________________________________________________________________________
48 AliEMCALRecPoint::AliEMCALRecPoint()
49   : AliRecPoint(),
50     fGeomPtr(0),
51     fClusterType(-1),
52     fCoreEnergy(0),
53     fDispersion(0),
54     fEnergyList(0),
55     fTimeList(0),
56     fAbsIdList(0),
57     fTime(0.),
58     fCoreRadius(10),  //HG check this
59     fMulParent(0),
60     fMaxParent(0),
61     fParentsList(0),
62     fSuperModuleNumber(0)
63 {
64   // ctor
65   fMaxTrack = 0 ;
66   fMulDigit = 0 ;
67   fAmp   = 0. ;   
68   //  fLocPos.SetX(1.e+6)  ;      //Local position should be evaluated
69
70   AliRunLoader *rl = AliRunLoader::GetRunLoader();
71   fGeomPtr = dynamic_cast<AliEMCAL*>(rl->GetAliRun()->GetDetector("EMCAL"))->GetGeometry();
72   //fGeomPtr = AliEMCALGeometry::GetInstance();
73   fGeomPtr->GetTransformationForSM(); // Global <-> Local
74 }
75
76 //____________________________________________________________________________
77 AliEMCALRecPoint::AliEMCALRecPoint(const char * opt) 
78   : AliRecPoint(opt),
79     fGeomPtr(0),
80     fClusterType(-1),
81     fCoreEnergy(0),
82     fDispersion(0),
83     fEnergyList(0),
84     fTimeList(0),
85     fAbsIdList(0),
86     fTime(-1.),
87     fCoreRadius(10),  //HG check this
88     fMulParent(0),
89     fMaxParent(1000),
90     fParentsList(0),
91     fSuperModuleNumber(0)
92 {
93   // ctor
94   fMaxTrack = 1000 ;
95   fMulDigit   = 0 ; 
96   fAmp   = 0. ;   
97   fParentsList = new Int_t[fMaxParent];
98
99   //fLocPos.SetX(1.e+6)  ;      //Local position should be evaluated
100   //fGeomPtr = AliEMCALGeometry::GetInstance();
101   AliRunLoader *rl = AliRunLoader::GetRunLoader();
102   fGeomPtr = dynamic_cast<AliEMCAL*>(rl->GetAliRun()->GetDetector("EMCAL"))->GetGeometry();
103   fGeomPtr->GetTransformationForSM(); // Global <-> Local
104 }
105
106 //____________________________________________________________________________
107 AliEMCALRecPoint::AliEMCALRecPoint(const AliEMCALRecPoint & rp) 
108   : AliRecPoint(rp),
109     fGeomPtr(rp.fGeomPtr),
110     fClusterType(rp.fClusterType),
111     fCoreEnergy(rp.fCoreEnergy),
112     fDispersion(rp.fDispersion),
113     fEnergyList(0),
114     fTimeList(0),
115     fAbsIdList(0),
116     fTime(rp.fTime),
117     fCoreRadius(rp.fCoreRadius),
118     fMulParent(rp.fMulParent),
119     fMaxParent(rp.fMaxParent),
120     fParentsList(0),
121     fSuperModuleNumber(rp.fSuperModuleNumber)
122 {
123   //copy ctor
124   fLambda[0] = rp.fLambda[0];
125   fLambda[1] = rp.fLambda[1];
126
127   fEnergyList = new Float_t[rp.fMulDigit];
128   fTimeList = new Float_t[rp.fMulDigit];
129   fAbsIdList = new Int_t[rp.fMulDigit];
130   for(Int_t i = 0; i < rp.fMulDigit; i++) {
131     fEnergyList[i] = rp.fEnergyList[i];
132     fTimeList[i] = rp.fTimeList[i];
133     fAbsIdList[i] = rp.fAbsIdList[i];
134   }
135   fParentsList = new Int_t[rp.fMulParent];
136   for(Int_t i = 0; i < rp.fMulParent; i++) fParentsList[i] = rp.fParentsList[i];
137
138 }
139 //____________________________________________________________________________
140 AliEMCALRecPoint::~AliEMCALRecPoint()
141 {
142   // dtor
143   if ( fEnergyList )
144     delete[] fEnergyList ; 
145   if ( fTimeList )
146     delete[] fTimeList ; 
147   if ( fAbsIdList )
148     delete[] fAbsIdList ; 
149    if ( fParentsList)
150     delete[] fParentsList;
151 }
152
153 //____________________________________________________________________________
154 void AliEMCALRecPoint::AddDigit(AliEMCALDigit & digit, Float_t Energy)
155 {
156   // Adds a digit to the RecPoint
157   // and accumulates the total amplitude and the multiplicity 
158   
159   if(fEnergyList == 0)
160     fEnergyList =  new Float_t[fMaxDigit]; 
161   if(fTimeList == 0)
162     fTimeList =  new Float_t[fMaxDigit]; 
163   if(fAbsIdList == 0) {
164     fAbsIdList =  new Int_t[fMaxDigit];
165     fSuperModuleNumber = fGeomPtr->GetSuperModuleNumber(digit.GetId());
166   }
167
168   if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists 
169     fMaxDigit*=2 ; 
170     Int_t * tempo = new Int_t[fMaxDigit]; 
171     Float_t * tempoE =  new Float_t[fMaxDigit];
172     Float_t * tempoT =  new Float_t[fMaxDigit];
173     Int_t * tempoId = new Int_t[fMaxDigit]; 
174
175     Int_t index ;     
176     for ( index = 0 ; index < fMulDigit ; index++ ){
177       tempo[index]   = fDigitsList[index] ;
178       tempoE[index]  = fEnergyList[index] ; 
179       tempoT[index]  = fTimeList[index] ; 
180       tempoId[index] = fAbsIdList[index] ; 
181     }
182     
183     delete [] fDigitsList ; 
184     fDigitsList =  new Int_t[fMaxDigit];
185  
186     delete [] fEnergyList ;
187     fEnergyList =  new Float_t[fMaxDigit];
188
189     delete [] fTimeList ;
190     fTimeList =  new Float_t[fMaxDigit];
191
192     delete [] fAbsIdList ;
193     fAbsIdList =  new Int_t[fMaxDigit];
194
195     for ( index = 0 ; index < fMulDigit ; index++ ){
196       fDigitsList[index] = tempo[index] ;
197       fEnergyList[index] = tempoE[index] ; 
198       fTimeList[index] = tempoT[index] ; 
199       fAbsIdList[index]  = tempoId[index] ; 
200     }
201  
202     delete [] tempo ;
203     delete [] tempoE ; 
204     delete [] tempoT ; 
205     delete [] tempoId ; 
206   } // if
207   
208   fDigitsList[fMulDigit]   = digit.GetIndexInList()  ; 
209   fEnergyList[fMulDigit]   = Energy ;
210   fTimeList[fMulDigit]     = digit.GetTime() ;
211   fAbsIdList[fMulDigit]    = digit.GetId();
212   fMulDigit++ ; 
213   fAmp += Energy ; 
214
215 }
216 //____________________________________________________________________________
217 Bool_t AliEMCALRecPoint::AreNeighbours(AliEMCALDigit * digit1, AliEMCALDigit * digit2 ) const
218 {
219   // Tells if (true) or not (false) two digits are neighbours
220   // A neighbour is defined as being two digits which share a corner
221   
222   static Bool_t areNeighbours = kFALSE ;
223   static Int_t nSupMod=0, nTower=0, nIphi=0, nIeta=0;
224   static int nSupMod1=0, nTower1=0, nIphi1=0, nIeta1=0;
225   static Int_t relid1[2] , relid2[2] ; // ieta, iphi
226   static Int_t rowdiff=0, coldiff=0;
227
228   areNeighbours = kFALSE ;
229
230   fGeomPtr->GetCellIndex(digit1->GetId(), nSupMod,nTower,nIphi,nIeta);
231   fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nTower,nIphi,nIeta, relid1[0],relid1[1]);
232
233   fGeomPtr->GetCellIndex(digit2->GetId(), nSupMod1,nTower1,nIphi1,nIeta1);
234   fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod1,nTower1,nIphi1,nIeta1, relid2[0],relid2[1]);
235   
236   rowdiff = TMath::Abs( relid1[0] - relid2[0] ) ;  
237   coldiff = TMath::Abs( relid1[1] - relid2[1] ) ;  
238
239   if (( coldiff <= 1 )  && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0)) 
240   areNeighbours = kTRUE ;
241   
242   return areNeighbours;
243 }
244
245 //____________________________________________________________________________
246 Int_t AliEMCALRecPoint::Compare(const TObject * obj) const
247 {
248   // Compares two RecPoints according to their position in the EMCAL modules
249
250   Float_t delta = 1 ; //Width of "Sorting row". If you change this 
251                       //value (what is senseless) change as well delta in
252                       //AliEMCALTrackSegmentMakerv* and other RecPoints...
253   Int_t rv ; 
254
255   AliEMCALRecPoint * clu = (AliEMCALRecPoint *)obj ; 
256
257   TVector3 locpos1; 
258   GetLocalPosition(locpos1);
259   TVector3 locpos2;  
260   clu->GetLocalPosition(locpos2);  
261
262   Int_t rowdif = (Int_t)(TMath::Ceil(locpos1.X()/delta)-TMath::Ceil(locpos2.X()/delta)) ;
263   if (rowdif> 0) 
264     rv = 1 ;
265   else if(rowdif < 0) 
266     rv = -1 ;
267   else if(locpos1.Y()>locpos2.Y()) 
268     rv = -1 ;
269   else 
270     rv = 1 ; 
271
272   return rv ; 
273 }
274
275 //____________________________________________________________________________
276 Int_t AliEMCALRecPoint::DistancetoPrimitive(Int_t px, Int_t py)
277 {
278   // Compute distance from point px,py to  a AliEMCALRecPoint considered as a Tmarker
279   // Compute the closest distance of approach from point px,py to this marker.
280   // The distance is computed in pixels units.
281   // HG Still need to update -> Not sure what this should achieve
282
283   TVector3 pos(0.,0.,0.) ;
284   GetLocalPosition(pos) ;
285   Float_t x =  pos.X() ;
286   Float_t y =  pos.Y() ;
287   const Int_t kMaxDiff = 10;
288   Int_t pxm  = gPad->XtoAbsPixel(x);
289   Int_t pym  = gPad->YtoAbsPixel(y);
290   Int_t dist = (px-pxm)*(px-pxm) + (py-pym)*(py-pym);
291   
292   if (dist > kMaxDiff) return 9999;
293   return dist;
294 }
295
296 //___________________________________________________________________________
297  void AliEMCALRecPoint::Draw(Option_t *option)
298  {
299    // Draw this AliEMCALRecPoint with its current attributes
300    
301    AppendPad(option);
302  }
303
304 //______________________________________________________________________________
305 void AliEMCALRecPoint::ExecuteEvent(Int_t /*event*/, Int_t, Int_t)
306 {
307   // Execute action corresponding to one event
308   // This member function is called when a AliEMCALRecPoint is clicked with the locator
309   //
310   // If Left button is clicked on AliEMCALRecPoint, the digits are switched on    
311   // and switched off when the mouse button is released.
312
313   //  static Int_t pxold, pyold;
314
315   /* static TGraph *  digitgraph = 0 ;
316   static TPaveText* clustertext = 0 ;
317   
318   if (!gPad->IsEditable()) return;
319   
320   switch (event) {
321     
322     
323   case kButton1Down:{
324     AliEMCALDigit * digit ;
325     AliEMCALGeometry * emcalgeom =  (AliEMCALGetter::Instance())->EMCALGeometry() ;
326
327     Int_t iDigit;
328     Int_t relid[2] ;
329   
330     const Int_t kMulDigit=AliEMCALRecPoint::GetDigitsMultiplicity() ;
331     Float_t * xi = new Float_t [kMulDigit] ; 
332     Float_t * zi = new Float_t [kMulDigit] ;
333     
334     for(iDigit = 0; iDigit < kMulDigit; iDigit++) {
335       Fatal("AliEMCALRecPoint::ExecuteEvent", " -> Something wrong with the code"); 
336       digit = 0 ; //dynamic_cast<AliEMCALDigit *>((fDigitsList)[iDigit]);
337       emcalgeom->AbsToRelNumbering(digit->GetId(), relid) ;
338       emcalgeom->PosInAlice(relid, xi[iDigit], zi[iDigit]) ;
339     }
340     
341     if (!digitgraph) {
342       digitgraph = new TGraph(fMulDigit,xi,zi);
343       digitgraph-> SetMarkerStyle(5) ; 
344       digitgraph-> SetMarkerSize(1.) ;
345       digitgraph-> SetMarkerColor(1) ;
346       digitgraph-> Draw("P") ;
347     }
348     if (!clustertext) {
349       
350       TVector3 pos(0.,0.,0.) ;
351       GetLocalPosition(pos) ;
352       clustertext = new TPaveText(pos.X()-10,pos.Z()+10,pos.X()+50,pos.Z()+35,"") ;
353       Text_t  line1[40] ;
354       Text_t  line2[40] ;
355       sprintf(line1,"Energy=%1.2f GeV",GetEnergy()) ;
356       sprintf(line2,"%d Digits",GetDigitsMultiplicity()) ;
357       clustertext ->AddText(line1) ;
358       clustertext ->AddText(line2) ;
359       clustertext ->Draw("");
360     }
361     gPad->Update() ; 
362     Print("") ;
363     delete[] xi ; 
364     delete[] zi ; 
365    }
366   
367 break;
368   
369   case kButton1Up:
370     if (digitgraph) {
371       delete digitgraph  ;
372       digitgraph = 0 ;
373     }
374     if (clustertext) {
375       delete clustertext ;
376       clustertext = 0 ;
377     }
378     
379     break;
380     
381     }*/
382 }
383 //____________________________________________________________________________
384 void AliEMCALRecPoint::EvalAll(Float_t logWeight,TClonesArray * digits) 
385 {
386   // Evaluates all shower parameters
387
388   EvalLocalPosition(logWeight, digits) ;
389   //  printf("eval position done\n");
390   EvalElipsAxis(logWeight, digits) ;
391   //  printf("eval axis done\n");
392   EvalDispersion(logWeight, digits) ;
393   //  printf("eval dispersion done\n");
394   //EvalCoreEnergy(logWeight, digits);
395  // printf("eval energy done\n");
396   EvalTime(digits) ;
397   //  printf("eval time done\n");
398
399   EvalPrimaries(digits) ;
400   //  printf("eval pri done\n");
401   EvalParents(digits);
402   //  printf("eval parent done\n");
403 }
404
405 //____________________________________________________________________________
406 void  AliEMCALRecPoint::EvalDispersion(Float_t logWeight, TClonesArray * digits)
407 {
408   // Calculates the dispersion of the shower at the origin of the RecPoint
409
410   Double_t d = 0., wtot = 0., w = 0., xyzi[3], diff=0.;
411   Int_t iDigit=0, nstat=0, i=0;
412   AliEMCALDigit * digit ;
413
414   // Calculates the centre of gravity in the local EMCAL-module coordinates   
415   if (!fLocPos.Mag()) 
416     EvalLocalPosition(logWeight, digits) ;
417   
418   // Calculates the dispersion in coordinates 
419   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
420     digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit])  ;
421
422     fGeomPtr->RelPosCellInSModule(digit->GetId(), xyzi[0], xyzi[1], xyzi[2]);
423     w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
424
425     if(w>0.0) {
426       wtot += w;
427       nstat++;
428       for(i=0; i<3; i++ ) {
429         diff = xyzi[i] - double(fLocPos[i]);
430         d   += w * diff*diff;
431       }
432     }
433   }
434   
435   if ( wtot > 0 && nstat>1) d /= wtot ;
436   else                      d = 0. ; 
437
438   fDispersion = TMath::Sqrt(d) ;
439 }
440
441 //____________________________________________________________________________
442 void AliEMCALRecPoint::EvalLocalPosition(Float_t logWeight, TClonesArray * digits)
443 {
444   // Calculates the center of gravity in the local EMCAL-module coordinates 
445   //  Info("Print", " logWeight %f : cluster energy %f ", logWeight, fAmp); // for testing
446   
447   AliEMCALDigit * digit;
448   Int_t i=0, nstat=0;
449   Double_t clXYZ[3]={0.,0.,0.}, clRmsXYZ[3]={0.,0.,0.}, xyzi[3], wtot=0., w=0.;
450
451   for(Int_t iDigit=0; iDigit<fMulDigit; iDigit++) {
452     digit = dynamic_cast<AliEMCALDigit *>(digits->At(fDigitsList[iDigit])) ;
453
454     fGeomPtr->RelPosCellInSModule(digit->GetId(), xyzi[0], xyzi[1], xyzi[2]);
455     // printf(" Id %i : Local x,y,z %f %f %f \n", digit->GetId(), xyzi[0], xyzi[1], xyzi[2]);
456
457     if(logWeight > 0.0) w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ));
458     else                w = fEnergyList[iDigit]; // just energy
459
460     if(w>0.0) {
461       wtot += w ;
462       nstat++;
463       for(i=0; i<3; i++ ) {
464         clXYZ[i]    += (w*xyzi[i]);
465         clRmsXYZ[i] += (w*xyzi[i]*xyzi[i]);
466       }
467     }
468   }
469   //  cout << " wtot " << wtot << endl;
470   if ( wtot > 0 ) { 
471     //    xRMS   = TMath::Sqrt(x2m - xMean*xMean);
472     for(i=0; i<3; i++ ) {
473       clXYZ[i] /= wtot;
474       if(nstat>1) {
475         clRmsXYZ[i] /= (wtot*wtot);  
476         clRmsXYZ[i]  =  clRmsXYZ[i] - clXYZ[i]*clXYZ[i];
477         if(clRmsXYZ[i] > 0.0) {
478           clRmsXYZ[i] =  TMath::Sqrt(clRmsXYZ[i]);
479         } else      clRmsXYZ[i] = 0;
480       } else        clRmsXYZ[i] = 0;
481     }    
482   } else {
483     for(i=0; i<3; i++ ) {
484       clXYZ[i] = clRmsXYZ[i] = -1.;
485     }
486   }
487   // clRmsXYZ[i] ??
488   fLocPos.SetX(clXYZ[0]);
489   fLocPos.SetY(clXYZ[1]);
490   fLocPos.SetZ(clXYZ[2]);
491     
492 //  if (gDebug==2)
493 //    printf("EvalLocalPosition: eta,phi,r = %f,%f,%f", fLocPos.X(), fLocPos.Y(), fLocPos.Z()) ; 
494   fLocPosM = 0 ; // covariance matrix
495 }
496
497 //void AliEMCALRecPoint::EvalLocalPositionSimple()
498 //{ // Weight is proportional of cell energy 
499 //}
500
501 //______________________________________________________________________________
502 void AliEMCALRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits)
503 {
504   // This function calculates energy in the core, 
505   // i.e. within a radius rad = fCoreEnergy around the center. Beyond this radius
506   // in accordance with shower profile the energy deposition 
507   // should be less than 2%
508
509   AliEMCALDigit * digit ;
510   const Float_t kDeg2Rad = TMath::DegToRad() ;
511
512   Int_t iDigit;
513
514   if (!fLocPos.Mag()) {
515     EvalLocalPosition(logWeight, digits);
516   }
517   
518   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
519     digit = (AliEMCALDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
520     
521     Float_t ieta = 0.0;
522     Float_t iphi = 0.0;
523     //fGeomPtr->PosInAlice(digit->GetId(), ieta, iphi);
524     fGeomPtr->EtaPhiFromIndex(digit->GetId(),ieta, iphi) ;
525     iphi = iphi * kDeg2Rad;
526   
527     Float_t distance = TMath::Sqrt((ieta-fLocPos.X())*(ieta-fLocPos.X())+(iphi-fLocPos.Y())*(iphi-fLocPos.Y())) ;
528     if(distance < fCoreRadius)
529       fCoreEnergy += fEnergyList[iDigit] ;
530   }
531   
532 }
533 //____________________________________________________________________________
534 void  AliEMCALRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
535 {
536   // Calculates the axis of the shower ellipsoid in eta and phi
537
538   Double_t wtot = 0. ;
539   Double_t x    = 0.;
540   Double_t z    = 0.;
541   Double_t dxx  = 0.;
542   Double_t dzz  = 0.;
543   Double_t dxz  = 0.;
544
545   const Float_t kDeg2Rad = TMath::DegToRad();
546   AliEMCALDigit * digit ;
547
548   TString gn(fGeomPtr->GetName());
549
550   Int_t iDigit;
551   
552   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
553     digit = (AliEMCALDigit *) digits->At(fDigitsList[iDigit])  ;
554     Float_t etai = 0. ;
555     Float_t phii = 0. ; 
556     if(gn.Contains("SHISH")) { // have to be change - Feb 28, 2006
557     //copied for shish-kebab geometry, ieta,iphi is cast as float as eta,phi conversion
558     // for this geometry does not exist
559       int nSupMod=0, nTower=0, nIphi=0, nIeta=0;
560       int iphi=0, ieta=0;
561       fGeomPtr->GetCellIndex(digit->GetId(), nSupMod,nTower,nIphi,nIeta);
562       fGeomPtr->GetCellPhiEtaIndexInSModule(nSupMod,nTower,nIphi,nIeta, iphi,ieta);
563       etai=(Float_t)ieta;
564       phii=(Float_t)iphi;
565     } else {
566       fGeomPtr->EtaPhiFromIndex(digit->GetId(), etai, phii);
567       phii = phii * kDeg2Rad;
568     }
569
570     Double_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
571
572     dxx  += w * etai * etai ;
573     x    += w * etai ;
574     dzz  += w * phii * phii ;
575     z    += w * phii ; 
576
577     dxz  += w * etai * phii ; 
578
579     wtot += w ;
580   }
581
582   if ( wtot > 0 ) { 
583     dxx /= wtot ;
584     x   /= wtot ;
585     dxx -= x * x ;
586     dzz /= wtot ;
587     z   /= wtot ;
588     dzz -= z * z ;
589     dxz /= wtot ;
590     dxz -= x * z ;
591
592     fLambda[0] =  0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
593     if(fLambda[0] > 0)
594       fLambda[0] = TMath::Sqrt(fLambda[0]) ;
595     else
596       fLambda[0] = 0;
597     
598     fLambda[1] =  0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
599
600     if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
601       fLambda[1] = TMath::Sqrt(fLambda[1]) ;
602     else
603       fLambda[1]= 0. ;
604   } else { 
605     fLambda[0]= 0. ;
606     fLambda[1]= 0. ;
607   }
608
609   //  printf("Evalaxis: lambdas  = %f,%f", fLambda[0],fLambda[1]) ; 
610
611 }
612
613 //______________________________________________________________________________
614 void  AliEMCALRecPoint::EvalPrimaries(TClonesArray * digits)
615 {
616   // Constructs the list of primary particles (tracks) which have contributed to this RecPoint
617   
618   AliEMCALDigit * digit ;
619   Int_t * tempo    = new Int_t[fMaxTrack] ;
620
621   Int_t index ;  
622   for ( index = 0 ; index < GetDigitsMultiplicity() ; index++ ) { // all digits
623     digit = dynamic_cast<AliEMCALDigit *>(digits->At( fDigitsList[index] )) ; 
624     Int_t nprimaries = digit->GetNprimary() ;
625     if ( nprimaries == 0 ) continue ;
626     Int_t * newprimaryarray = new Int_t[nprimaries] ;
627     Int_t ii ; 
628     for ( ii = 0 ; ii < nprimaries ; ii++)
629       newprimaryarray[ii] = digit->GetPrimary(ii+1) ; 
630
631     Int_t jndex ;
632     for ( jndex = 0 ; jndex < nprimaries ; jndex++ ) { // all primaries in digit
633       if ( fMulTrack > fMaxTrack ) {
634         fMulTrack = fMaxTrack ;
635         Error("GetNprimaries", "increase fMaxTrack ")  ;
636         break ;
637       }
638       Int_t newprimary = newprimaryarray[jndex] ;
639       Int_t kndex ;
640       Bool_t already = kFALSE ;
641       for ( kndex = 0 ; kndex < fMulTrack ; kndex++ ) { //check if not already stored
642         if ( newprimary == tempo[kndex] ){
643           already = kTRUE ;
644           break ;
645         }
646       } // end of check
647       if ( !already && (fMulTrack < fMaxTrack)) { // store it
648         tempo[fMulTrack] = newprimary ; 
649         fMulTrack++ ;
650       } // store it
651     } // all primaries in digit
652     delete [] newprimaryarray ; 
653   } // all digits
654
655   
656   fTracksList = new Int_t[fMulTrack] ;
657   for(index = 0; index < fMulTrack; index++)
658    fTracksList[index] = tempo[index] ;
659  
660   delete [] tempo ;
661
662 }
663
664 //______________________________________________________________________________
665 void  AliEMCALRecPoint::EvalParents(TClonesArray * digits)
666 {
667   // Constructs the list of parent particles (tracks) which have contributed to this RecPoint
668  
669   AliEMCALDigit * digit ;
670   Int_t * tempo    = new Int_t[fMaxParent] ;
671
672   Int_t index ;  
673   for ( index = 0 ; index < GetDigitsMultiplicity() ; index++ ) { // all digits
674     digit = dynamic_cast<AliEMCALDigit *>(digits->At( fDigitsList[index] )) ; 
675     Int_t nparents = digit->GetNiparent() ;
676     if ( nparents == 0 ) continue ;
677     Int_t * newparentarray = new Int_t[nparents] ;
678     Int_t ii ; 
679     for ( ii = 0 ; ii < nparents ; ii++)
680       newparentarray[ii] = digit->GetIparent(ii+1) ; 
681
682     Int_t jndex ;
683     for ( jndex = 0 ; jndex < nparents ; jndex++ ) { // all primaries in digit
684       if ( fMulParent > fMaxParent ) {
685         fMulTrack = - 1 ;
686         Error("GetNiparent", "increase fMaxParent")  ;
687         break ;
688       }
689       Int_t newparent = newparentarray[jndex] ;
690       Int_t kndex ;
691       Bool_t already = kFALSE ;
692       for ( kndex = 0 ; kndex < fMulParent ; kndex++ ) { //check if not already stored
693         if ( newparent == tempo[kndex] ){
694           already = kTRUE ;
695           break ;
696         }
697       } // end of check
698       if ( !already && (fMulTrack < fMaxTrack)) { // store it
699         tempo[fMulParent] = newparent ; 
700         fMulParent++ ;
701       } // store it
702     } // all parents in digit
703     delete [] newparentarray ; 
704   } // all digits
705
706   if (fMulParent>0) {
707     fParentsList = new Int_t[fMulParent] ;
708     for(index = 0; index < fMulParent; index++)
709       fParentsList[index] = tempo[index] ;
710   }
711  
712   delete [] tempo ;
713
714 }
715
716 //____________________________________________________________________________
717 void AliEMCALRecPoint::GetLocalPosition(TVector3 & lpos) const
718 {
719   // returns the position of the cluster in the local reference system of ALICE
720   
721   lpos.SetX(fLocPos.X()) ;
722   lpos.SetY(fLocPos.Y()) ;
723   lpos.SetZ(fLocPos.Z()) ;
724 }
725
726 //____________________________________________________________________________
727 void AliEMCALRecPoint::GetGlobalPosition(TVector3 & gpos) const
728 {
729   // returns the position of the cluster in the global reference system of ALICE
730   // These are now the Cartesian X, Y and Z
731   //  cout<<" geom "<<geom<<endl;
732   fGeomPtr->GetGlobal(fLocPos, gpos, fSuperModuleNumber);
733 }
734
735 //____________________________________________________________________________
736 Float_t AliEMCALRecPoint::GetMaximalEnergy(void) const
737 {
738   // Finds the maximum energy in the cluster
739   
740   Float_t menergy = 0. ;
741
742   Int_t iDigit;
743
744   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
745  
746     if(fEnergyList[iDigit] > menergy) 
747       menergy = fEnergyList[iDigit] ;
748   }
749   return menergy ;
750 }
751
752 //____________________________________________________________________________
753 Int_t AliEMCALRecPoint::GetMultiplicityAtLevel(Float_t H) const
754 {
755   // Calculates the multiplicity of digits with energy larger than H*energy 
756   
757   Int_t multipl   = 0 ;
758   Int_t iDigit ;
759   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
760
761     if(fEnergyList[iDigit] > H * fAmp) 
762       multipl++ ;
763   }
764   return multipl ;
765 }
766
767 //____________________________________________________________________________
768 Int_t  AliEMCALRecPoint::GetNumberOfLocalMax(AliEMCALDigit **  maxAt, Float_t * maxAtEnergy,
769                                                Float_t locMaxCut,TClonesArray * digits) const
770
771   // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
772   // energy difference between two local maxima
773
774   AliEMCALDigit * digit ;
775   AliEMCALDigit * digitN ;
776   
777   Int_t iDigitN ;
778   Int_t iDigit ;
779
780   for(iDigit = 0; iDigit < fMulDigit; iDigit++)
781     maxAt[iDigit] = (AliEMCALDigit*) digits->At(fDigitsList[iDigit])  ;
782   
783   for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {   
784     if(maxAt[iDigit]) {
785       digit = maxAt[iDigit] ;
786           
787       for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {        
788         digitN = (AliEMCALDigit *) digits->At(fDigitsList[iDigitN]) ; 
789         
790         if ( AreNeighbours(digit, digitN) ) {
791           if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {    
792             maxAt[iDigitN] = 0 ;
793             // but may be digit too is not local max ?
794             if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut) 
795               maxAt[iDigit] = 0 ;
796           }
797           else {
798             maxAt[iDigit] = 0 ;
799             // but may be digitN too is not local max ?
800             if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut) 
801               maxAt[iDigitN] = 0 ; 
802           } 
803         } // if Areneighbours
804       } // while digitN
805     } // slot not empty
806   } // while digit
807   
808   iDigitN = 0 ;
809   for(iDigit = 0; iDigit < fMulDigit; iDigit++) { 
810     if(maxAt[iDigit] ){
811       maxAt[iDigitN] = maxAt[iDigit] ;
812       maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
813       iDigitN++ ; 
814     }
815   }
816   return iDigitN ;
817 }
818
819 //____________________________________________________________________________
820 Int_t AliEMCALRecPoint::GetPrimaryIndex() const  
821 {
822   // Get the primary track index in TreeK which deposits the most energy 
823   // in Digits which forms RecPoint. Kinematics, Hits and Digits must be 
824   // loaded before the call of the method.
825
826   AliRunLoader *rl = AliRunLoader::GetRunLoader(); 
827   if (!rl) 
828     AliError(Form(" No Runloader ")) ; 
829  
830   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>
831     (rl->GetDetectorLoader("EMCAL"));
832
833   // Get the list of digits forming this RecPoint
834   Int_t  nDigits   = fMulDigit   ;
835   Int_t *digitList = fDigitsList ;
836   
837   // Find the digit with maximum amplitude
838   AliEMCALDigit *digit = 0;
839   TClonesArray *digits = emcalLoader->Digits();
840   Int_t maxAmp = 0;
841   Int_t bestDigitIndex = -1;
842   for (Int_t iDigit=0; iDigit<nDigits; iDigit++) {
843     digit = static_cast<AliEMCALDigit *>(digits->At(digitList[iDigit]));
844     if (digit->GetAmp() > maxAmp) {
845       maxAmp = digit->GetAmp();
846       bestDigitIndex = iDigit;
847     }
848   }
849
850   digit = static_cast<AliEMCALDigit *>(digits->At(digitList[bestDigitIndex]));  
851
852   // Get the list of hits producing this digit,
853   // find which hit has deposited more energy 
854   // and find the primary track.
855
856   AliEMCALHit *hit = 0;
857   TClonesArray *hits = emcalLoader->Hits();
858
859   Double_t maxedep  =  0;
860   Int_t    maxtrack = -1;
861   Int_t    nHits    = hits ->GetEntries();
862   Int_t    id       = digit->GetId();
863   for (Int_t iHit=0; iHit<nHits; iHit++) {
864     hit = static_cast<AliEMCALHit*> (hits->At(iHit)) ;
865     if(hit->GetId() == id){
866       Double_t edep  = hit->GetEnergy();
867       Int_t    track = hit->GetIparent();//Primary();
868       if(edep > maxedep){
869         maxedep  = edep;
870         maxtrack = track;
871       }
872     }
873   }
874   if (maxtrack != -1) return maxtrack; 
875   return -12345;                       // no track found :(
876 }
877
878 //____________________________________________________________________________
879 void AliEMCALRecPoint::EvalTime(TClonesArray * digits){
880   // time is set to the time of the digit with the maximum energy
881
882   Float_t maxE = 0;
883   Int_t maxAt = 0;
884   for(Int_t idig=0; idig < fMulDigit; idig++){
885     if(fEnergyList[idig] > maxE){
886       maxE = fEnergyList[idig] ;
887       maxAt = idig;
888     }
889   }
890   fTime = ((AliEMCALDigit*) digits->At(fDigitsList[maxAt]))->GetTime() ;
891   
892 }
893
894 //______________________________________________________________________________
895 void AliEMCALRecPoint::Paint(Option_t *)
896 {
897   // Paint this ALiRecPoint as a TMarker  with its current attributes
898   
899   TVector3 pos(0.,0.,0.)  ;
900   GetLocalPosition(pos)   ;
901   Coord_t x = pos.X()     ;
902   Coord_t y = pos.Z()     ;
903   Color_t markercolor = 1 ;
904   Size_t  markersize = 1. ;
905   Style_t markerstyle = 5 ;
906   
907   if (!gPad->IsBatch()) {
908     gVirtualX->SetMarkerColor(markercolor) ;
909     gVirtualX->SetMarkerSize (markersize)  ;
910     gVirtualX->SetMarkerStyle(markerstyle) ;
911   }
912   gPad->SetAttMarkerPS(markercolor,markerstyle,markersize) ;
913   gPad->PaintPolyMarker(1,&x,&y,"") ;
914 }
915
916 //______________________________________________________________________________
917 Float_t AliEMCALRecPoint::EtaToTheta(Float_t arg) const
918 {
919   //Converts Theta (Radians) to Eta(Radians)
920   return (2.*TMath::ATan(TMath::Exp(-arg)));
921 }
922
923 //______________________________________________________________________________
924 Float_t AliEMCALRecPoint::ThetaToEta(Float_t arg) const
925 {
926   //Converts Eta (Radians) to Theta(Radians)
927   return (-1 * TMath::Log(TMath::Tan(0.5 * arg)));
928 }
929
930 //____________________________________________________________________________
931 void AliEMCALRecPoint::Print(Option_t *) const
932 {
933   // Print the list of digits belonging to the cluster
934   return;
935   TString message ; 
936   message  = "AliEMCALRecPoint:\n" ;
937   message +=  " digits # = " ; 
938   Info("Print", message.Data()) ; 
939
940   Int_t iDigit;
941   for(iDigit=0; iDigit<fMulDigit; iDigit++)
942     printf(" %d ", fDigitsList[iDigit] ) ;  
943   printf("\n");
944
945   Info("Print", " Energies = ") ;
946   for(iDigit=0; iDigit<fMulDigit; iDigit++) 
947     printf(" %f ", fEnergyList[iDigit] ) ;
948   printf("\n");
949
950   Info("Print", "\n Abs Ids = ") ;
951   for(iDigit=0; iDigit<fMulDigit; iDigit++) 
952     printf(" %i ", fAbsIdList[iDigit] ) ;
953   printf("\n");
954
955   Info("Print", " Primaries  ") ;
956   for(iDigit = 0;iDigit < fMulTrack; iDigit++)
957     printf(" %d ", fTracksList[iDigit]) ;
958
959   printf("\n Local x %6.2f y %7.2f z %7.1f \n", fLocPos[0], fLocPos[1], fLocPos[2]);
960
961   message  = "       ClusterType     = %d" ;
962   message += "       Multiplicity    = %d" ;
963   message += "       Cluster Energy  = %f" ; 
964   message += "       Core energy     = %f" ; 
965   message += "       Core radius     = %f" ; 
966   message += "       Number of primaries %d" ; 
967   message += "       Stored at position %d" ; 
968   Info("Print", message.Data(), fClusterType, fMulDigit, fAmp, fCoreEnergy, fCoreRadius, fMulTrack, GetIndexInList() ) ;  
969 }