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