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