]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSEmcRecPoint.cxx
clustering version one
[u/mrichter/AliRoot.git] / PHOS / AliPHOSEmcRecPoint.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
16 /* $Id$ */
17
18 /* History of cvs commits:
19  *
20  * $Log$
21  * Revision 1.52  2005/05/28 14:19:04  schutz
22  * Compilation warnings fixed by T.P.
23  *
24  */
25
26 //_________________________________________________________________________
27 //  RecPoint implementation for PHOS-EMC 
28 //  An EmcRecPoint is a cluster of digits   
29 //*--
30 //*-- Author: Dmitri Peressounko (RRC KI & SUBATECH)
31
32
33 // --- ROOT system ---
34 #include "TH2.h"
35 #include "TMath.h" 
36 #include "TCanvas.h" 
37 #include "TGraph.h"
38
39 // --- Standard library ---
40
41 // --- AliRoot header files ---
42 #include "AliLog.h"
43 #include "AliPHOSLoader.h"
44 #include "AliGenerator.h"
45 #include "AliPHOSGeometry.h"
46 #include "AliPHOSDigit.h"
47 #include "AliPHOSEmcRecPoint.h"
48  
49 ClassImp(AliPHOSEmcRecPoint)
50
51 //____________________________________________________________________________
52 AliPHOSEmcRecPoint::AliPHOSEmcRecPoint() : AliPHOSRecPoint()
53 {
54   // ctor
55
56   fMulDigit   = 0 ;  
57   fAmp   = 0. ;   
58   fCoreEnergy = 0 ; 
59   fEnergyList = 0 ;
60   fNExMax     = 0 ;   //Not unfolded yet
61   fTime = -1. ;
62   fLocPos.SetX(1000000.)  ;      //Local position should be evaluated
63   fDebug=0;
64    
65 }
66
67 //____________________________________________________________________________
68 AliPHOSEmcRecPoint::AliPHOSEmcRecPoint(const char * opt) : AliPHOSRecPoint(opt)
69 {
70   // ctor
71   
72   fMulDigit   = 0 ;  
73   fAmp   = 0. ;   
74   fNExMax     = 0 ;   //Not unfolded yet
75   fCoreEnergy = 0 ; 
76   fEnergyList = 0 ;
77   fTime = -1. ;
78   fLocPos.SetX(1000000.)  ;      //Local position should be evaluated
79   fDebug=0;
80   
81 }
82
83 //____________________________________________________________________________
84 AliPHOSEmcRecPoint::AliPHOSEmcRecPoint(const AliPHOSEmcRecPoint & rp) : AliPHOSRecPoint(rp)
85 {
86   // cpy ctor
87
88   fMulDigit   = rp.fMulDigit ;  
89   fAmp        = rp.fAmp ;   
90   fCoreEnergy = rp.fCoreEnergy ; 
91   fEnergyList = new Float_t[rp.fMulDigit] ;
92   Int_t index ; 
93   for(index = 0 ; index < fMulDigit ; index++) 
94     fEnergyList[index] = rp.fEnergyList[index] ; 
95   fNExMax     = rp.fNExMax ;  
96   fTime       = rp.fTime ;   
97 }
98
99 //____________________________________________________________________________
100 AliPHOSEmcRecPoint::~AliPHOSEmcRecPoint()
101 {
102   // dtor
103
104   if ( fEnergyList )
105     delete[] fEnergyList ; 
106 }
107
108 //____________________________________________________________________________
109 void AliPHOSEmcRecPoint::AddDigit(AliPHOSDigit & digit, Float_t Energy)
110 {
111   // Adds a digit to the RecPoint
112   // and accumulates the total amplitude and the multiplicity 
113   
114   if(fEnergyList == 0)
115     fEnergyList =  new Float_t[fMaxDigit]; 
116
117   if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists 
118     fMaxDigit*=2 ; 
119     Int_t * tempo = new Int_t[fMaxDigit]; 
120     Float_t * tempoE =  new Float_t[fMaxDigit];
121
122     Int_t index ;     
123     for ( index = 0 ; index < fMulDigit ; index++ ){
124       tempo[index]  = fDigitsList[index] ;
125       tempoE[index] = fEnergyList[index] ; 
126     }
127     
128     delete [] fDigitsList ; 
129     fDigitsList =  new Int_t[fMaxDigit];
130  
131     delete [] fEnergyList ;
132     fEnergyList =  new Float_t[fMaxDigit];
133
134     for ( index = 0 ; index < fMulDigit ; index++ ){
135       fDigitsList[index] = tempo[index] ;
136       fEnergyList[index] = tempoE[index] ; 
137     }
138  
139     delete [] tempo ;
140     delete [] tempoE ; 
141   } // if
142   
143   fDigitsList[fMulDigit]   = digit.GetIndexInList()  ; 
144   fEnergyList[fMulDigit]   = Energy ;
145   fMulDigit++ ; 
146   fAmp += Energy ; 
147
148   EvalPHOSMod(&digit) ;
149 }
150
151 //____________________________________________________________________________
152 Bool_t AliPHOSEmcRecPoint::AreNeighbours(AliPHOSDigit * digit1, AliPHOSDigit * digit2 ) const
153 {
154   // Tells if (true) or not (false) two digits are neighbors
155   
156   Bool_t aren = kFALSE ;
157   
158   AliPHOSGeometry * phosgeom =  AliPHOSGeometry::GetInstance() ;
159
160   Int_t relid1[4] ; 
161   phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ; 
162
163   Int_t relid2[4] ; 
164   phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ; 
165   
166   Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;  
167   Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;  
168
169   if (( coldiff <= 1 )  && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0)) 
170     aren = kTRUE ;
171   
172   return aren ;
173 }
174
175 //____________________________________________________________________________
176 Int_t AliPHOSEmcRecPoint::Compare(const TObject * obj) const
177 {
178   // Compares two RecPoints according to their position in the PHOS modules
179   
180   const Float_t delta = 1 ; //Width of "Sorting row". If you changibg this 
181                       //value (what is senseless) change as vell delta in
182                       //AliPHOSTrackSegmentMakerv* and other RecPoints...
183   Int_t rv ; 
184
185   AliPHOSEmcRecPoint * clu = (AliPHOSEmcRecPoint *)obj ; 
186
187  
188   Int_t phosmod1 = GetPHOSMod() ;
189   Int_t phosmod2 = clu->GetPHOSMod() ;
190
191   TVector3 locpos1; 
192   GetLocalPosition(locpos1) ;
193   TVector3 locpos2;  
194   clu->GetLocalPosition(locpos2) ;  
195
196   if(phosmod1 == phosmod2 ) {
197     Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
198     if (rowdif> 0) 
199       rv = 1 ;
200      else if(rowdif < 0) 
201        rv = -1 ;
202     else if(locpos1.Z()>locpos2.Z()) 
203       rv = -1 ;
204     else 
205       rv = 1 ; 
206      }
207
208   else {
209     if(phosmod1 < phosmod2 ) 
210       rv = -1 ;
211     else 
212       rv = 1 ;
213   }
214
215   return rv ; 
216 }
217 //______________________________________________________________________________
218 void AliPHOSEmcRecPoint::ExecuteEvent(Int_t event, Int_t, Int_t) /*const*/
219 {
220   
221   // Execute action corresponding to one event
222   //  This member function is called when a AliPHOSRecPoint is clicked with the locator
223   //
224   //  If Left button is clicked on AliPHOSRecPoint, the digits are switched on    
225   //  and switched off when the mouse button is released.
226   
227   
228   AliPHOSGeometry * phosgeom =  AliPHOSGeometry::GetInstance();
229   
230   static TGraph *  digitgraph = 0 ;
231   
232   if (!gPad->IsEditable()) return;
233   
234   TH2F * histo = 0 ;
235   TCanvas * histocanvas ; 
236
237   
238   //try to get run loader from default event folder
239   AliRunLoader* rn = AliRunLoader::GetRunLoader(AliConfig::GetDefaultEventFolderName());
240   if (rn == 0x0) 
241     {
242       AliError(Form("Cannot find Run Loader in Default Event Folder"));
243       return;
244     }
245   AliPHOSLoader* gime = dynamic_cast<AliPHOSLoader*>(rn->GetLoader("PHOSLoader"));
246   if (gime == 0x0) 
247     {
248       AliError(Form("Cannot find PHOS Loader from Run Loader"));
249       return;
250     }
251   
252   
253   const TClonesArray * digits = gime->Digits() ;
254   
255   switch (event) {
256     
257   case kButton1Down: {
258     AliPHOSDigit * digit ;
259     Int_t iDigit;
260     Int_t relid[4] ;
261     
262     const Int_t kMulDigit = AliPHOSEmcRecPoint::GetDigitsMultiplicity() ; 
263     Float_t * xi = new Float_t[kMulDigit] ; 
264     Float_t * zi = new Float_t[kMulDigit] ; 
265     
266     // create the histogram for the single cluster 
267     // 1. gets histogram boundaries
268     Float_t ximax = -999. ; 
269     Float_t zimax = -999. ; 
270     Float_t ximin = 999. ; 
271     Float_t zimin = 999. ;
272     
273     for(iDigit=0; iDigit<kMulDigit; iDigit++) {
274       digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit])  ;
275       phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
276       phosgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
277       if ( xi[iDigit] > ximax )
278         ximax = xi[iDigit] ; 
279       if ( xi[iDigit] < ximin )
280         ximin = xi[iDigit] ; 
281       if ( zi[iDigit] > zimax )
282         zimax = zi[iDigit] ; 
283       if ( zi[iDigit] < zimin )
284         zimin = zi[iDigit] ;     
285     }
286     ximax += phosgeom->GetCrystalSize(0) / 2. ;
287     zimax += phosgeom->GetCrystalSize(2) / 2. ;
288     ximin -= phosgeom->GetCrystalSize(0) / 2. ;
289     zimin -= phosgeom->GetCrystalSize(2) / 2. ;
290     Int_t xdim = (int)( (ximax - ximin ) / phosgeom->GetCrystalSize(0) + 0.5  ) ; 
291     Int_t zdim = (int)( (zimax - zimin ) / phosgeom->GetCrystalSize(2) + 0.5 ) ;
292     
293     // 2. gets the histogram title
294     
295     Text_t title[100] ; 
296     sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
297     
298     if (!histo) {
299       delete histo ; 
300       histo = 0 ; 
301     }
302     histo = new TH2F("cluster3D", title,  xdim, ximin, ximax, zdim, zimin, zimax)  ;
303     
304     Float_t x, z ; 
305     for(iDigit=0; iDigit<kMulDigit; iDigit++) {
306       digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit])  ;
307       phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
308       phosgeom->RelPosInModule(relid, x, z);
309       histo->Fill(x, z, fEnergyList[iDigit] ) ;
310     }
311     
312     if (!digitgraph) {
313       digitgraph = new TGraph(kMulDigit,xi,zi);
314       digitgraph-> SetMarkerStyle(5) ; 
315       digitgraph-> SetMarkerSize(1.) ;
316       digitgraph-> SetMarkerColor(1) ;
317       digitgraph-> Paint("P") ;
318     }
319     
320     //    Print() ;
321     histocanvas = new TCanvas("cluster", "a single cluster", 600, 500) ; 
322     histocanvas->Draw() ; 
323     histo->Draw("lego1") ; 
324     
325     delete[] xi ; 
326     delete[] zi ; 
327     
328     break;
329   }
330   
331   case kButton1Up: 
332     if (digitgraph) {
333       delete digitgraph  ;
334       digitgraph = 0 ;
335     }
336     break;
337   
338    }
339 }
340
341 //____________________________________________________________________________
342 void  AliPHOSEmcRecPoint::EvalDispersion(Float_t logWeight,TClonesArray * digits)
343 {
344   // Calculates the dispersion of the shower at the origine of the RecPoint
345
346   Float_t d    = 0. ;
347   Float_t wtot = 0. ;
348
349   Float_t x = 0.;
350   Float_t z = 0.;
351
352   AliPHOSDigit * digit ;
353  
354   AliPHOSGeometry * phosgeom =  AliPHOSGeometry::GetInstance();
355
356  // Calculates the center of gravity in the local PHOS-module coordinates 
357   
358   Int_t iDigit;
359
360   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
361     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
362     Int_t relid[4] ;
363     Float_t xi ;
364     Float_t zi ;
365     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
366     phosgeom->RelPosInModule(relid, xi, zi);
367     if (fAmp>0 && fEnergyList[iDigit]>0) {
368       Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
369       x    += xi * w ;
370       z    += zi * w ;
371       wtot += w ;
372     }
373     else
374       AliError(Form("Wrong energy %f and/or amplitude %f\n", fEnergyList[iDigit], fAmp));
375   }
376   if (wtot>0) {
377     x /= wtot ;
378     z /= wtot ;
379   }
380   else
381     AliError(Form("Wrong weight %f\n", wtot));
382
383
384 // Calculates the dispersion in coordinates 
385   wtot = 0.;
386   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
387     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit])  ;
388     Int_t relid[4] ;
389     Float_t xi ;
390     Float_t zi ;
391     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
392     phosgeom->RelPosInModule(relid, xi, zi);
393     if (fAmp>0 && fEnergyList[iDigit]>0) {
394       Float_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
395       d += w*((xi-x)*(xi-x) + (zi-z)*(zi-z) ) ; 
396       wtot+=w ;
397     }
398     else
399       AliError(Form("Wrong energy %f and/or amplitude %f\n", fEnergyList[iDigit], fAmp));
400   }
401   
402
403   if (wtot>0) {
404     d /= wtot ;
405   }
406   else
407     AliError(Form("Wrong weight %f\n", wtot));
408
409   fDispersion = 0;
410   if (d>=0)
411     fDispersion = TMath::Sqrt(d) ;
412  
413 }
414 //______________________________________________________________________________
415 void AliPHOSEmcRecPoint::EvalCoreEnergy(Float_t logWeight, TClonesArray * digits)
416 {
417   // This function calculates energy in the core, 
418   // i.e. within a radius rad = 3cm around the center. Beyond this radius
419   // in accordance with shower profile the energy deposition 
420   // should be less than 2%
421
422   Float_t coreRadius = 3 ;
423
424   Float_t x = 0 ;
425   Float_t z = 0 ;
426
427   AliPHOSDigit * digit ;
428
429   AliPHOSGeometry * phosgeom =  AliPHOSGeometry::GetInstance();
430
431   Int_t iDigit;
432
433 // Calculates the center of gravity in the local PHOS-module coordinates 
434   Float_t wtot = 0;
435   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
436     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
437     Int_t relid[4] ;
438     Float_t xi ;
439     Float_t zi ;
440     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
441     phosgeom->RelPosInModule(relid, xi, zi);
442     if (fAmp>0 && fEnergyList[iDigit]>0) {
443       Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
444       x    += xi * w ;
445       z    += zi * w ;
446       wtot += w ;
447     }
448     else
449       AliError(Form("Wrong energy %f and/or amplitude %f\n", fEnergyList[iDigit], fAmp));
450   }
451   if (wtot>0) {
452     x /= wtot ;
453     z /= wtot ;
454   }
455   else
456     AliError(Form("Wrong weight %f\n", wtot));
457
458
459   for(iDigit=0; iDigit < fMulDigit; iDigit++) {
460     digit = (AliPHOSDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
461     Int_t relid[4] ;
462     Float_t xi ;
463     Float_t zi ;
464     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
465     phosgeom->RelPosInModule(relid, xi, zi);    
466     Float_t distance = TMath::Sqrt((xi-x)*(xi-x)+(zi-z)*(zi-z)) ;
467     if(distance < coreRadius)
468       fCoreEnergy += fEnergyList[iDigit] ;
469   }
470
471 }
472
473 //____________________________________________________________________________
474 void  AliPHOSEmcRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits)
475 {
476   // Calculates the axis of the shower ellipsoid
477
478   Double_t wtot = 0. ;
479   Double_t x    = 0.;
480   Double_t z    = 0.;
481   Double_t dxx  = 0.;
482   Double_t dzz  = 0.;
483   Double_t dxz  = 0.;
484
485   AliPHOSDigit * digit ;
486
487   AliPHOSGeometry * phosgeom =  AliPHOSGeometry::GetInstance();
488
489   Int_t iDigit;
490
491
492   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
493     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit])  ;
494     Int_t relid[4] ;
495     Float_t xi ;
496     Float_t zi ;
497     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
498     phosgeom->RelPosInModule(relid, xi, zi);
499     if (fAmp>0 && fEnergyList[iDigit]>0) {
500       Double_t w = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
501       dxx  += w * xi * xi ;
502       x    += w * xi ;
503       dzz  += w * zi * zi ;
504       z    += w * zi ; 
505       dxz  += w * xi * zi ; 
506       wtot += w ;
507     }
508     else
509       AliError(Form("Wrong energy %f and/or amplitude %f\n", fEnergyList[iDigit], fAmp));
510   }
511   if (wtot>0) {
512     dxx /= wtot ;
513     x   /= wtot ;
514     dxx -= x * x ;
515     dzz /= wtot ;
516     z   /= wtot ;
517     dzz -= z * z ;
518     dxz /= wtot ;
519     dxz -= x * z ;
520
521 //   //Apply correction due to non-perpendicular incidence
522 //   Double_t CosX ;
523 //   Double_t CosZ ;
524 //   AliPHOSGetter * gime = AliPHOSGetter::Instance() ; 
525 //   AliPHOSGeometry * phosgeom =  (AliPHOSGeometry*)gime->PHOSGeometry();
526   //   Double_t DistanceToIP= (Double_t ) phosgeom->GetIPtoCrystalSurface() ;
527   
528 //   CosX = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+x*x) ;
529 //   CosZ = DistanceToIP/TMath::Sqrt(DistanceToIP*DistanceToIP+z*z) ;
530
531 //   dxx = dxx/(CosX*CosX) ;
532 //   dzz = dzz/(CosZ*CosZ) ;
533 //   dxz = dxz/(CosX*CosZ) ;
534
535
536     fLambda[0] =  0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
537     if(fLambda[0] > 0)
538       fLambda[0] = TMath::Sqrt(fLambda[0]) ;
539
540     fLambda[1] =  0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
541     if(fLambda[1] > 0) //To avoid exception if numerical errors lead to negative lambda.
542       fLambda[1] = TMath::Sqrt(fLambda[1]) ;
543     else
544       fLambda[1]= 0. ;
545   }
546   else {
547     AliError(Form("Wrong weight %f\n", wtot));
548     fLambda[0]=fLambda[1]=0.;
549   }
550 }
551
552 //____________________________________________________________________________
553 void  AliPHOSEmcRecPoint::EvalMoments(Float_t logWeight,TClonesArray * digits)
554 {
555   // Calculate the shower moments in the eigen reference system
556   // M2x, M2z, M3x, M4z
557   // Calculate the angle between the shower position vector and the eigen vector
558
559   Double_t wtot = 0. ;
560   Double_t x    = 0.;
561   Double_t z    = 0.;
562   Double_t dxx  = 0.;
563   Double_t dzz  = 0.;
564   Double_t dxz  = 0.;
565   Double_t lambda0=0, lambda1=0;
566
567   AliPHOSDigit * digit ;
568
569   AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance() ;
570
571   Int_t iDigit;
572
573   // 1) Find covariance matrix elements:
574   //    || dxx dxz ||
575   //    || dxz dzz ||
576
577   Float_t xi ;
578   Float_t zi ;
579   Int_t relid[4] ;
580   Double_t w;
581   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
582     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit])  ;
583     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
584     phosgeom->RelPosInModule(relid, xi, zi);
585     if (fAmp>0 && fEnergyList[iDigit]>0) {
586       w     = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
587       x    += w * xi ;
588       z    += w * zi ; 
589       dxx  += w * xi * xi ;
590       dzz  += w * zi * zi ;
591       dxz  += w * xi * zi ; 
592       wtot += w ;
593     }
594     else
595       AliError(Form("Wrong energy %f and/or amplitude %f\n", fEnergyList[iDigit], fAmp));
596     
597   }
598   if (wtot>0) {
599     x   /= wtot ;
600     z   /= wtot ;
601     dxx /= wtot ;
602     dzz /= wtot ;
603     dxz /= wtot ;
604     dxx -= x * x ;
605     dzz -= z * z ;
606     dxz -= x * z ;
607
608   // 2) Find covariance matrix eigen values lambda0 and lambda1
609
610     lambda0 =  0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
611     lambda1 =  0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
612   }
613   else {
614     AliError(Form("Wrong weight %f\n", wtot));
615     lambda0=lambda1=0.;
616   }
617
618   // 3) Find covariance matrix eigen vectors e0 and e1
619
620   TVector2 e0,e1;
621   if (dxz != 0)
622     e0.Set(1.,(lambda0-dxx)/dxz);
623   else
624     e0.Set(0.,1.);
625
626   e0 = e0.Unit();
627   e1.Set(-e0.Y(),e0.X());
628
629   // 4) Rotate cluster tensor from (x,z) to (e0,e1) system
630   //    and calculate moments M3x and M4z
631
632   Float_t cosPhi = e0.X();
633   Float_t sinPhi = e0.Y();
634
635   Float_t xiPHOS ;
636   Float_t ziPHOS ;
637   Double_t dx3, dz3, dz4;
638   wtot = 0.;
639   x    = 0.;
640   z    = 0.;
641   dxx  = 0.;
642   dzz  = 0.;
643   dxz  = 0.;
644   dx3  = 0.;
645   dz3  = 0.;
646   dz4  = 0.;
647   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
648     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit])  ;
649     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
650     phosgeom->RelPosInModule(relid, xiPHOS, ziPHOS);
651     xi    = xiPHOS*cosPhi + ziPHOS*sinPhi;
652     zi    = ziPHOS*cosPhi - xiPHOS*sinPhi;
653     if (fAmp>0 && fEnergyList[iDigit]>0) {
654       w     = TMath::Max(0.,logWeight+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
655       x    += w * xi ;
656       z    += w * zi ; 
657       dxx  += w * xi * xi ;
658       dzz  += w * zi * zi ;
659       dxz  += w * xi * zi ; 
660       dx3  += w * xi * xi * xi;
661       dz3  += w * zi * zi * zi ;
662       dz4  += w * zi * zi * zi * zi ;
663       wtot += w ;
664     }
665     else
666       AliError(Form("Wrong energy %f and/or amplitude %f\n", fEnergyList[iDigit], fAmp));
667   }
668   if (wtot>0) {
669     x   /= wtot ;
670     z   /= wtot ;
671     dxx /= wtot ;
672     dzz /= wtot ;
673     dxz /= wtot ;
674     dx3 /= wtot ;
675     dz3 /= wtot ;
676     dz4 /= wtot ;
677     dx3 += -3*dxx*x + 2*x*x*x;
678     dz4 += -4*dz3*z + 6*dzz*z*z -3*z*z*z*z;
679     dxx -= x * x ;
680     dzz -= z * z ;
681     dxz -= x * z ;
682   }
683   else
684     AliError(Form("Wrong weight %f\n", wtot));
685
686   // 5) Find an angle between cluster center vector and eigen vector e0
687
688   Float_t phi = 0;
689   if ( (x*x+z*z) > 0 ) 
690     phi = TMath::ACos ((x*e0.X() + z*e0.Y()) / sqrt(x*x + z*z));
691
692   fM2x   = lambda0;
693   fM2z   = lambda1;
694   fM3x   = dx3;
695   fM4z   = dz4;
696   fPhixe = phi;
697
698 }
699 //____________________________________________________________________________
700 void AliPHOSEmcRecPoint::EvalAll(Float_t logWeight, TClonesArray * digits )
701 {
702   // Evaluates all shower parameters
703   EvalLocalPosition(logWeight, digits) ;
704   EvalElipsAxis(logWeight, digits) ;
705   EvalMoments(logWeight, digits) ;
706   EvalDispersion(logWeight, digits) ;
707   EvalCoreEnergy(logWeight, digits);
708   EvalTime(digits) ;
709   AliPHOSRecPoint::EvalAll(digits) ;
710 }
711 //____________________________________________________________________________
712 void AliPHOSEmcRecPoint::EvalLocalPosition(Float_t logWeight, TClonesArray * digits)
713 {
714   // Calculates the center of gravity in the local PHOS-module coordinates 
715   Float_t wtot = 0. ;
716  
717   Int_t relid[4] ;
718
719   Float_t x = 0. ;
720   Float_t z = 0. ;
721   
722   AliPHOSDigit * digit ;
723
724   AliPHOSGeometry * phosgeom =  AliPHOSGeometry::GetInstance() ;
725
726   Int_t iDigit;
727
728   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
729     digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
730
731     Float_t xi ;
732     Float_t zi ;
733     phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
734     phosgeom->RelPosInModule(relid, xi, zi);
735     if (fAmp>0 && fEnergyList[iDigit]>0) {
736       Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
737       x    += xi * w ;
738       z    += zi * w ;
739       wtot += w ;
740     }
741     else
742       AliError(Form("Wrong energy %f and/or amplitude %f\n", fEnergyList[iDigit], fAmp));
743   }
744   if (wtot>0) {
745     x /= wtot ;
746     z /= wtot ;
747   }
748   else
749     AliError(Form("Wrong weight %f\n", wtot));
750
751   // Correction for the depth of the shower starting point (TDR p 127)  
752   Float_t para = 0.925 ; 
753   Float_t parb = 6.52 ; 
754
755   Float_t xo,yo,zo ; //Coordinates of the origin
756   //We should check all 3 possibilities to avoid seg.v.
757   if(gAlice && gAlice->GetMCApp() && gAlice->Generator())
758     gAlice->Generator()->GetOrigin(xo,yo,zo) ;
759   else{
760     xo=yo=zo=0.;
761   }
762   Float_t phi = phosgeom->GetPHOSAngle(relid[0]) ;
763
764   //Transform to the local ref.frame
765   Float_t xoL,yoL ;
766   xoL = xo*TMath::Cos(phi)-yo*TMath::Sin(phi) ;
767   yoL = xo*TMath::Sin(phi)+yo*TMath::Cos(phi) ;
768   
769   Float_t radius = phosgeom->GetIPtoCrystalSurface()-yoL;
770  
771   Float_t incidencephi = TMath::ATan((x-xoL ) / radius) ; 
772   Float_t incidencetheta = TMath::ATan((z-zo) / radius) ;
773  
774   Float_t depthx = 0.; 
775   Float_t depthz = 0.;
776   if (fAmp>0) {
777     depthx = ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencephi) ; 
778     depthz = ( para * TMath::Log(fAmp) + parb ) * TMath::Sin(incidencetheta) ; 
779   }
780   else 
781     AliError(Form("Wrong amplitude %f\n", fAmp));
782
783   fLocPos.SetX(x - depthx)  ;
784   fLocPos.SetY(0.) ;
785   fLocPos.SetZ(z - depthz)  ;
786
787   fLocPosM = 0 ;
788 }
789
790 //____________________________________________________________________________
791 Float_t AliPHOSEmcRecPoint::GetMaximalEnergy(void) const
792 {
793   // Finds the maximum energy in the cluster
794   
795   Float_t menergy = 0. ;
796
797   Int_t iDigit;
798
799   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
800  
801     if(fEnergyList[iDigit] > menergy) 
802       menergy = fEnergyList[iDigit] ;
803   }
804   return menergy ;
805 }
806
807 //____________________________________________________________________________
808 Int_t AliPHOSEmcRecPoint::GetMultiplicityAtLevel(Float_t H) const
809 {
810   // Calculates the multiplicity of digits with energy larger than H*energy 
811   
812   Int_t multipl   = 0 ;
813   Int_t iDigit ;
814   for(iDigit=0; iDigit<fMulDigit; iDigit++) {
815
816     if(fEnergyList[iDigit] > H * fAmp) 
817       multipl++ ;
818   }
819   return multipl ;
820 }
821
822 //____________________________________________________________________________
823 Int_t  AliPHOSEmcRecPoint::GetNumberOfLocalMax( AliPHOSDigit **  maxAt, Float_t * maxAtEnergy,
824                                                Float_t locMaxCut,TClonesArray * digits) const
825
826   // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
827   // energy difference between two local maxima
828
829   AliPHOSDigit * digit ;
830   AliPHOSDigit * digitN ;
831   
832
833   Int_t iDigitN ;
834   Int_t iDigit ;
835
836   for(iDigit = 0; iDigit < fMulDigit; iDigit++)
837     maxAt[iDigit] = (AliPHOSDigit*) digits->At(fDigitsList[iDigit])  ;
838
839   
840   for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {   
841     if(maxAt[iDigit]) {
842       digit = maxAt[iDigit] ;
843           
844       for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {        
845         if(iDigit == iDigitN)
846           continue ;
847         
848         digitN = (AliPHOSDigit *) digits->At(fDigitsList[iDigitN]) ; 
849         
850         if ( AreNeighbours(digit, digitN) ) {
851           if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {    
852             maxAt[iDigitN] = 0 ;
853             // but may be digit too is not local max ?
854             if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut) 
855               maxAt[iDigit] = 0 ;
856           }
857           else {
858             maxAt[iDigit] = 0 ;
859             // but may be digitN too is not local max ?
860             if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut) 
861               maxAt[iDigitN] = 0 ; 
862           } 
863         } // if Areneighbours
864       } // while digitN
865     } // slot not empty
866   } // while digit
867   
868   iDigitN = 0 ;
869   for(iDigit = 0; iDigit < fMulDigit; iDigit++) { 
870     if(maxAt[iDigit]){
871       maxAt[iDigitN] = maxAt[iDigit] ;
872       maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
873       iDigitN++ ; 
874     }
875   }
876   return iDigitN ;
877 }
878 //____________________________________________________________________________
879 void AliPHOSEmcRecPoint::EvalTime(TClonesArray * digits)
880 {
881   // Define a rec.point time as a time in the cell with the maximum energy
882
883   Float_t maxE = 0;
884   Int_t maxAt = 0;
885   for(Int_t idig=0; idig < fMulDigit; idig++){
886     if(fEnergyList[idig] > maxE){
887       maxE = fEnergyList[idig] ;
888       maxAt = idig;
889     }
890   }
891   fTime = ((AliPHOSDigit*) digits->At(fDigitsList[maxAt]))->GetTime() ;
892   
893 }
894 //____________________________________________________________________________
895 void AliPHOSEmcRecPoint::Purify(Float_t threshold){
896   //Removes digits below threshold
897
898   Int_t * tempo = new Int_t[fMaxDigit]; 
899   Float_t * tempoE =  new Float_t[fMaxDigit];
900
901   Int_t mult = 0 ;
902   for(Int_t iDigit=0;iDigit< fMulDigit ;iDigit++){
903     if(fEnergyList[iDigit] > threshold){
904       tempo[mult] = fDigitsList[iDigit] ;
905       tempoE[mult] = fEnergyList[iDigit] ;
906       mult++ ;
907     }
908   }
909   
910   fMulDigit = mult ;
911   delete [] fDigitsList ;
912   delete [] fEnergyList ;
913   fDigitsList = new Int_t[fMulDigit];
914   fEnergyList = new Float_t[fMulDigit];
915
916   fAmp = 0.; //Recalculate total energy
917   for(Int_t iDigit=0;iDigit< fMulDigit ;iDigit++){
918      fDigitsList[iDigit] = tempo[iDigit];
919      fEnergyList[iDigit] = tempoE[iDigit] ;
920      fAmp+=tempoE[iDigit];
921   }
922       
923   delete [] tempo ;
924   delete [] tempoE ;
925
926 }
927 //____________________________________________________________________________
928 void AliPHOSEmcRecPoint::Print(Option_t *) const
929 {
930   // Print the list of digits belonging to the cluster
931   
932   TString message ; 
933   message  = "AliPHOSEmcRecPoint:\n" ;
934   message +=  " digits # = " ; 
935   AliInfo(Form(message.Data())) ; 
936
937   Int_t iDigit;
938   for(iDigit=0; iDigit<fMulDigit; iDigit++)
939     printf(" %d ", fDigitsList[iDigit] ) ;  
940   
941   printf(" Energies = ") ;
942   for(iDigit=0; iDigit<fMulDigit; iDigit++) 
943     printf(" %f ", fEnergyList[iDigit] ) ;
944   printf("\n") ; 
945    printf(" Primaries  ") ;
946   for(iDigit = 0;iDigit < fMulTrack; iDigit++)
947     printf(" %d ", fTracksList[iDigit]) ;
948   printf("\n") ;        
949   message  = "       Multiplicity    = %d" ;
950   message += "       Cluster Energy  = %f" ; 
951   message += "       Number of primaries %d" ; 
952   message += "       Stored at position %d" ; 
953  
954   printf(message.Data(), fMulDigit, fAmp, fMulTrack,GetIndexInList() ) ;  
955 }
956  
957