]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PHOS/AliPHOSEmcRecPoint.cxx
Update task for ITS tracking check and new macro for plots
[u/mrichter/AliRoot.git] / PHOS / AliPHOSEmcRecPoint.cxx
CommitLineData
d15a28e7 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
b2a60966 16/* $Id$ */
17
702ab87e 18/* History of cvs commits:
19 *
20 * $Log$
25779c68 21 * Revision 1.59 2007/10/18 15:12:22 kharlov
22 * Moved MakePrimary to EMCRecPoint to rpduce correct order of primaries
23 *
c307c629 24 * Revision 1.58 2007/04/16 09:03:37 kharlov
25 * Incedent angle correction fixed
26 *
753b19cd 27 * Revision 1.57 2007/04/05 10:18:58 policheh
28 * Introduced distance to nearest bad crystal.
29 *
dce76bbb 30 * Revision 1.56 2007/03/06 06:47:28 kharlov
31 * DP:Possibility to use actual vertex position added
32 *
91daaf24 33 * Revision 1.55 2007/01/19 20:31:19 kharlov
34 * Improved formatting for Print()
702ab87e 35 */
36
d15a28e7 37//_________________________________________________________________________
b2a60966 38// RecPoint implementation for PHOS-EMC
39// An EmcRecPoint is a cluster of digits
6c8cd883 40//--
41//-- Author: Dmitri Peressounko (RRC KI & SUBATECH)
b2a60966 42
d15a28e7 43
44// --- ROOT system ---
9f616d61 45#include "TH2.h"
d15a28e7 46#include "TMath.h"
9f616d61 47#include "TCanvas.h"
55fe9d13 48#include "TGraph.h"
d15a28e7 49
50// --- Standard library ---
51
d15a28e7 52// --- AliRoot header files ---
351dd634 53#include "AliLog.h"
e957fea8 54#include "AliPHOSLoader.h"
88cb7938 55#include "AliGenerator.h"
d15a28e7 56#include "AliPHOSGeometry.h"
e957fea8 57#include "AliPHOSDigit.h"
d15a28e7 58#include "AliPHOSEmcRecPoint.h"
48c5db5b 59#include "AliPHOSReconstructor.h"
e957fea8 60
d15a28e7 61ClassImp(AliPHOSEmcRecPoint)
62
63//____________________________________________________________________________
3663622c 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.),
dce76bbb 69 fPhixe(0.), fDistToBadCrystal(-1),fDebug(0)
d15a28e7 70{
71 // ctor
d15a28e7 72 fMulDigit = 0 ;
73 fAmp = 0. ;
a6eedfad 74 fLocPos.SetX(1000000.) ; //Local position should be evaluated
83974468 75}
76
73a68ccb 77//____________________________________________________________________________
3663622c 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.),
dce76bbb 83 fPhixe(0.), fDistToBadCrystal(-1), fDebug(0)
73a68ccb 84{
85 // ctor
73a68ccb 86 fMulDigit = 0 ;
87 fAmp = 0. ;
73a68ccb 88 fLocPos.SetX(1000000.) ; //Local position should be evaluated
73a68ccb 89}
90
55fe9d13 91//____________________________________________________________________________
3663622c 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),
dce76bbb 97 fPhixe(rp.fPhixe), fDistToBadCrystal(rp.fDistToBadCrystal), fDebug(rp.fDebug)
55fe9d13 98{
99 // cpy ctor
55fe9d13 100 fMulDigit = rp.fMulDigit ;
101 fAmp = rp.fAmp ;
25779c68 102 if (rp.fMulDigit>0) fEnergyList = new Float_t[rp.fMulDigit] ;
103 for(Int_t index = 0 ; index < fMulDigit ; index++)
55fe9d13 104 fEnergyList[index] = rp.fEnergyList[index] ;
55fe9d13 105}
106
83974468 107//____________________________________________________________________________
108AliPHOSEmcRecPoint::~AliPHOSEmcRecPoint()
109{
88714635 110 // dtor
83974468 111 if ( fEnergyList )
112 delete[] fEnergyList ;
d15a28e7 113}
114
d15a28e7 115//____________________________________________________________________________
6f47f50d 116void AliPHOSEmcRecPoint::AddDigit(AliPHOSDigit & digit, Float_t Energy, Float_t time)
d15a28e7 117{
b2a60966 118 // Adds a digit to the RecPoint
a4e98857 119 // and accumulates the total amplitude and the multiplicity
d15a28e7 120
7932f811 121 if(fEnergyList == 0)
122 fEnergyList = new Float_t[fMaxDigit];
123
d15a28e7 124 if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists
9f616d61 125 fMaxDigit*=2 ;
0fbb8738 126 Int_t * tempo = new Int_t[fMaxDigit];
127 Float_t * tempoE = new Float_t[fMaxDigit];
9f616d61 128
129 Int_t index ;
d15a28e7 130 for ( index = 0 ; index < fMulDigit ; index++ ){
83974468 131 tempo[index] = fDigitsList[index] ;
d15a28e7 132 tempoE[index] = fEnergyList[index] ;
133 }
134
9f616d61 135 delete [] fDigitsList ;
0fbb8738 136 fDigitsList = new Int_t[fMaxDigit];
9f616d61 137
138 delete [] fEnergyList ;
0fbb8738 139 fEnergyList = new Float_t[fMaxDigit];
9f616d61 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
d15a28e7 149
6f47f50d 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
83974468 164 fDigitsList[fMulDigit] = digit.GetIndexInList() ;
165 fEnergyList[fMulDigit] = Energy ;
166 fMulDigit++ ;
d15a28e7 167 fAmp += Energy ;
7932f811 168 EvalPHOSMod(&digit) ;
d15a28e7 169}
170
171//____________________________________________________________________________
ad8cfaf4 172Bool_t AliPHOSEmcRecPoint::AreNeighbours(AliPHOSDigit * digit1, AliPHOSDigit * digit2 ) const
d15a28e7 173{
a4e98857 174 // Tells if (true) or not (false) two digits are neighbors
d15a28e7 175
176 Bool_t aren = kFALSE ;
177
686b771f 178 AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance() ;
7b7c1533 179
d15a28e7 180 Int_t relid1[4] ;
92862013 181 phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ;
d15a28e7 182
183 Int_t relid2[4] ;
92862013 184 phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ;
d15a28e7 185
92862013 186 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
187 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
d15a28e7 188
92862013 189 if (( coldiff <= 1 ) && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0))
d15a28e7 190 aren = kTRUE ;
191
192 return aren ;
193}
194
195//____________________________________________________________________________
2a941f4e 196Int_t AliPHOSEmcRecPoint::Compare(const TObject * obj) const
d15a28e7 197{
b2a60966 198 // Compares two RecPoints according to their position in the PHOS modules
686b771f 199
200 const Float_t delta = 1 ; //Width of "Sorting row". If you changibg this
7932f811 201 //value (what is senseless) change as vell delta in
202 //AliPHOSTrackSegmentMakerv* and other RecPoints...
d15a28e7 203 Int_t rv ;
204
205 AliPHOSEmcRecPoint * clu = (AliPHOSEmcRecPoint *)obj ;
206
207
ad8cfaf4 208 Int_t phosmod1 = GetPHOSMod() ;
92862013 209 Int_t phosmod2 = clu->GetPHOSMod() ;
d15a28e7 210
92862013 211 TVector3 locpos1;
7932f811 212 GetLocalPosition(locpos1) ;
92862013 213 TVector3 locpos2;
214 clu->GetLocalPosition(locpos2) ;
d15a28e7 215
92862013 216 if(phosmod1 == phosmod2 ) {
7932f811 217 Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
d15a28e7 218 if (rowdif> 0)
7932f811 219 rv = 1 ;
d15a28e7 220 else if(rowdif < 0)
7932f811 221 rv = -1 ;
92862013 222 else if(locpos1.Z()>locpos2.Z())
d15a28e7 223 rv = -1 ;
224 else
225 rv = 1 ;
226 }
227
228 else {
92862013 229 if(phosmod1 < phosmod2 )
d15a28e7 230 rv = -1 ;
231 else
232 rv = 1 ;
233 }
234
235 return rv ;
236}
9f616d61 237//______________________________________________________________________________
702ab87e 238void AliPHOSEmcRecPoint::ExecuteEvent(Int_t event, Int_t, Int_t) /*const*/
9f616d61 239{
9f616d61 240
9688c1dd 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
88cb7938 247
686b771f 248 AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance();
83974468 249
9688c1dd 250 static TGraph * digitgraph = 0 ;
83974468 251
9688c1dd 252 if (!gPad->IsEditable()) return;
83974468 253
9688c1dd 254 TH2F * histo = 0 ;
255 TCanvas * histocanvas ;
256
88cb7938 257
258 //try to get run loader from default event folder
e191bb57 259 AliRunLoader* rn = AliRunLoader::GetRunLoader(AliConfig::GetDefaultEventFolderName());
88cb7938 260 if (rn == 0x0)
261 {
351dd634 262 AliError(Form("Cannot find Run Loader in Default Event Folder"));
88cb7938 263 return;
264 }
7fb9892d 265 AliPHOSLoader* phosLoader = dynamic_cast<AliPHOSLoader*>(rn->GetLoader("PHOSLoader"));
266 if (phosLoader == 0x0)
88cb7938 267 {
351dd634 268 AliError(Form("Cannot find PHOS Loader from Run Loader"));
88cb7938 269 return;
270 }
271
272
7fb9892d 273 const TClonesArray * digits = phosLoader->Digits() ;
9688c1dd 274
275 switch (event) {
83974468 276
9688c1dd 277 case kButton1Down: {
278 AliPHOSDigit * digit ;
279 Int_t iDigit;
280 Int_t relid[4] ;
83974468 281
9688c1dd 282 const Int_t kMulDigit = AliPHOSEmcRecPoint::GetDigitsMultiplicity() ;
283 Float_t * xi = new Float_t[kMulDigit] ;
284 Float_t * zi = new Float_t[kMulDigit] ;
83974468 285
9688c1dd 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. ;
83974468 292
9688c1dd 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 ) ;
83974468 312
9688c1dd 313 // 2. gets the histogram title
83974468 314
9688c1dd 315 Text_t title[100] ;
316 sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
83974468 317
9688c1dd 318 if (!histo) {
319 delete histo ;
320 histo = 0 ;
321 }
322 histo = new TH2F("cluster3D", title, xdim, ximin, ximax, zdim, zimin, zimax) ;
83974468 323
9688c1dd 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 }
83974468 331
9688c1dd 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 }
83974468 339
9688c1dd 340 // Print() ;
341 histocanvas = new TCanvas("cluster", "a single cluster", 600, 500) ;
342 histocanvas->Draw() ;
343 histo->Draw("lego1") ;
83974468 344
9688c1dd 345 delete[] xi ;
346 delete[] zi ;
83974468 347
9688c1dd 348 break;
349 }
83974468 350
9688c1dd 351 case kButton1Up:
352 if (digitgraph) {
353 delete digitgraph ;
354 digitgraph = 0 ;
355 }
356 break;
9f616d61 357
9688c1dd 358 }
9f616d61 359}
360
d15a28e7 361//____________________________________________________________________________
308fb942 362void AliPHOSEmcRecPoint::EvalDispersion(Float_t logWeight,TClonesArray * digits, TVector3 & /* vInc */)
d15a28e7 363{
b2a60966 364 // Calculates the dispersion of the shower at the origine of the RecPoint
91daaf24 365 //DP: should we correct dispersion for non-perpendicular hit????????
366
e5b16749 367 Float_t d = 0. ;
368 Float_t wtot = 0. ;
d15a28e7 369
d084d50d 370 Float_t x = 0.;
371 Float_t z = 0.;
d15a28e7 372
373 AliPHOSDigit * digit ;
7b7c1533 374
686b771f 375 AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance();
e5b16749 376
377 // Calculates the center of gravity in the local PHOS-module coordinates
378
d15a28e7 379 Int_t iDigit;
e5b16749 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);
f33a3764 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));
e5b16749 396 }
f33a3764 397 if (wtot>0) {
398 x /= wtot ;
399 z /= wtot ;
400 }
401 else
402 AliError(Form("Wrong weight %f\n", wtot));
e5b16749 403
404
405// Calculates the dispersion in coordinates
406 wtot = 0.;
88714635 407 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
7932f811 408 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
d15a28e7 409 Int_t relid[4] ;
410 Float_t xi ;
411 Float_t zi ;
92862013 412 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
413 phosgeom->RelPosInModule(relid, xi, zi);
f33a3764 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));
d15a28e7 421 }
e5b16749 422
7932f811 423
f33a3764 424 if (wtot>0) {
425 d /= wtot ;
426 }
427 else
428 AliError(Form("Wrong weight %f\n", wtot));
d15a28e7 429
f33a3764 430 fDispersion = 0;
431 if (d>=0)
432 fDispersion = TMath::Sqrt(d) ;
91daaf24 433
c63c49e9 434
d15a28e7 435}
fad3e5b9 436//______________________________________________________________________________
d55c9739 437void AliPHOSEmcRecPoint::EvalCoreEnergy(Float_t logWeight, Float_t coreRadius, TClonesArray * digits)
fad3e5b9 438{
a4e98857 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
fad3e5b9 442 // should be less than 2%
91daaf24 443//DP: non-perpendicular incidence??????????????
fad3e5b9 444
e5b16749 445 Float_t x = 0 ;
446 Float_t z = 0 ;
fad3e5b9 447
448 AliPHOSDigit * digit ;
7b7c1533 449
686b771f 450 AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance();
88cb7938 451
fad3e5b9 452 Int_t iDigit;
e5b16749 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);
f33a3764 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));
e5b16749 471 }
f33a3764 472 if (wtot>0) {
473 x /= wtot ;
474 z /= wtot ;
475 }
476 else
477 AliError(Form("Wrong weight %f\n", wtot));
e5b16749 478
479
fad3e5b9 480 for(iDigit=0; iDigit < fMulDigit; iDigit++) {
7932f811 481 digit = (AliPHOSDigit *) ( digits->At(fDigitsList[iDigit]) ) ;
fad3e5b9 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)
7932f811 489 fCoreEnergy += fEnergyList[iDigit] ;
fad3e5b9 490 }
491
91daaf24 492
fad3e5b9 493}
d15a28e7 494
495//____________________________________________________________________________
308fb942 496void AliPHOSEmcRecPoint::EvalElipsAxis(Float_t logWeight,TClonesArray * digits, TVector3 & /* vInc */)
d15a28e7 497{
b2a60966 498 // Calculates the axis of the shower ellipsoid
83974468 499
e8dbb96e 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.;
d15a28e7 506
507 AliPHOSDigit * digit ;
7b7c1533 508
686b771f 509 AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance();
7b7c1533 510
d15a28e7 511 Int_t iDigit;
512
e5b16749 513
d15a28e7 514 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
7932f811 515 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
d15a28e7 516 Int_t relid[4] ;
517 Float_t xi ;
518 Float_t zi ;
92862013 519 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
520 phosgeom->RelPosInModule(relid, xi, zi);
f33a3764 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));
d15a28e7 532 }
f33a3764 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 ;
d15a28e7 542
69183710 543// //Apply correction due to non-perpendicular incidence
544// Double_t CosX ;
545// Double_t CosZ ;
7fb9892d 546// AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance() ;
91daaf24 547// Double_t DistanceToIP= (Double_t ) phosgeom->GetIPtoCrystalSurface() ;
69183710 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
f33a3764 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]) ;
e8dbb96e 560
f33a3764 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 }
d15a28e7 571}
572
ce2a9a95 573//____________________________________________________________________________
308fb942 574void AliPHOSEmcRecPoint::EvalMoments(Float_t logWeight,TClonesArray * digits, TVector3 & /* vInc */)
ce2a9a95 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
686b771f 590 AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance() ;
ce2a9a95 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);
f33a3764 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
ce2a9a95 618 }
f33a3764 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 ;
ce2a9a95 628
629 // 2) Find covariance matrix eigen values lambda0 and lambda1
630
f33a3764 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 }
ce2a9a95 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;
f33a3764 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));
ce2a9a95 688 }
f33a3764 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));
ce2a9a95 706
707 // 5) Find an angle between cluster center vector and eigen vector e0
708
f33a3764 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));
ce2a9a95 712
713 fM2x = lambda0;
714 fM2z = lambda1;
715 fM3x = dx3;
716 fM4z = dz4;
717 fPhixe = phi;
718
719}
c307c629 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 }
7fb719de 742 for(fMulTrack=0; fMulTrack<nprimaries ; fMulTrack++){
743 tempo[fMulTrack] = digit->GetPrimary(fMulTrack+1) ;
c307c629 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
c307c629 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
b2a60966 783//____________________________________________________________________________
e347da7b 784void AliPHOSEmcRecPoint::EvalAll(TClonesArray * digits )
ad8cfaf4 785{
d55c9739 786// EvalCoreEnergy(logWeight, digits);
a6eedfad 787 EvalTime(digits) ;
c307c629 788 EvalPrimaries(digits) ;
e957fea8 789 AliPHOSRecPoint::EvalAll(digits) ;
ad8cfaf4 790}
791//____________________________________________________________________________
91daaf24 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)
b2a60966 803{
804 // Calculates the center of gravity in the local PHOS-module coordinates
b2a60966 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
686b771f 814 AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance() ;
b2a60966 815
816 Int_t iDigit;
817
b2a60966 818 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
7932f811 819 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
b2a60966 820
821 Float_t xi ;
822 Float_t zi ;
823 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
824 phosgeom->RelPosInModule(relid, xi, zi);
f33a3764 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));
9ce1a8d9 833 }
f33a3764 834 if (wtot>0) {
835 x /= wtot ;
836 z /= wtot ;
837 }
838 else
839 AliError(Form("Wrong weight %f\n", wtot));
ad8cfaf4 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
91daaf24 845 phosgeom->GetIncidentVector(vtx,GetPHOSMod(),x,z,vInc) ;
1b799736 846
f33a3764 847 Float_t depthx = 0.;
848 Float_t depthz = 0.;
91daaf24 849 if (fAmp>0&&vInc.Y()!=0.) {
753b19cd 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()) ;
f33a3764 852 }
853 else
854 AliError(Form("Wrong amplitude %f\n", fAmp));
ad8cfaf4 855
856 fLocPos.SetX(x - depthx) ;
b2a60966 857 fLocPos.SetY(0.) ;
ad8cfaf4 858 fLocPos.SetZ(z - depthz) ;
b2a60966 859
b2a60966 860}
861
d15a28e7 862//____________________________________________________________________________
ad8cfaf4 863Float_t AliPHOSEmcRecPoint::GetMaximalEnergy(void) const
d15a28e7 864{
b2a60966 865 // Finds the maximum energy in the cluster
866
d15a28e7 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//____________________________________________________________________________
fc7e2f43 880Int_t AliPHOSEmcRecPoint::GetMultiplicityAtLevel(Float_t H) const
d15a28e7 881{
b2a60966 882 // Calculates the multiplicity of digits with energy larger than H*energy
883
d15a28e7 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//____________________________________________________________________________
a0636361 895Int_t AliPHOSEmcRecPoint::GetNumberOfLocalMax( AliPHOSDigit ** maxAt, Float_t * maxAtEnergy,
7932f811 896 Float_t locMaxCut,TClonesArray * digits) const
d15a28e7 897{
b2a60966 898 // Calculates the number of local maxima in the cluster using fLocalMaxCut as the minimum
a4e98857 899 // energy difference between two local maxima
b2a60966 900
d15a28e7 901 AliPHOSDigit * digit ;
902 AliPHOSDigit * digitN ;
903
904
905 Int_t iDigitN ;
906 Int_t iDigit ;
907
7932f811 908 for(iDigit = 0; iDigit < fMulDigit; iDigit++)
a0636361 909 maxAt[iDigit] = (AliPHOSDigit*) digits->At(fDigitsList[iDigit]) ;
7932f811 910
d15a28e7 911
6ad0bfa0 912 for(iDigit = 0 ; iDigit < fMulDigit; iDigit++) {
a0636361 913 if(maxAt[iDigit]) {
914 digit = maxAt[iDigit] ;
83974468 915
6ad0bfa0 916 for(iDigitN = 0; iDigitN < fMulDigit; iDigitN++) {
686b771f 917 if(iDigit == iDigitN)
918 continue ;
919
7932f811 920 digitN = (AliPHOSDigit *) digits->At(fDigitsList[iDigitN]) ;
d15a28e7 921
9f616d61 922 if ( AreNeighbours(digit, digitN) ) {
d15a28e7 923 if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {
a0636361 924 maxAt[iDigitN] = 0 ;
6ad0bfa0 925 // but may be digit too is not local max ?
7932f811 926 if(fEnergyList[iDigit] < fEnergyList[iDigitN] + locMaxCut)
a0636361 927 maxAt[iDigit] = 0 ;
d15a28e7 928 }
929 else {
a0636361 930 maxAt[iDigit] = 0 ;
6ad0bfa0 931 // but may be digitN too is not local max ?
7932f811 932 if(fEnergyList[iDigit] > fEnergyList[iDigitN] - locMaxCut)
a0636361 933 maxAt[iDigitN] = 0 ;
d15a28e7 934 }
935 } // if Areneighbours
936 } // while digitN
937 } // slot not empty
938 } // while digit
939
940 iDigitN = 0 ;
6ad0bfa0 941 for(iDigit = 0; iDigit < fMulDigit; iDigit++) {
a0636361 942 if(maxAt[iDigit]){
d15a28e7 943 maxAt[iDigitN] = maxAt[iDigit] ;
9f616d61 944 maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
945 iDigitN++ ;
d15a28e7 946 }
947 }
948 return iDigitN ;
949}
9688c1dd 950//____________________________________________________________________________
6f47f50d 951void AliPHOSEmcRecPoint::EvalTime(TClonesArray * /*digits*/)
0bc3b8ed 952{
953 // Define a rec.point time as a time in the cell with the maximum energy
6f47f50d 954 //Time already evaluated during AddDigit()
0bc3b8ed 955
6f47f50d 956/*
9688c1dd 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() ;
6f47f50d 966*/
9688c1dd 967}
d15a28e7 968//____________________________________________________________________________
092b50ba 969void AliPHOSEmcRecPoint::Purify(Float_t threshold){
970 //Removes digits below threshold
971
0fbb8738 972 Int_t * tempo = new Int_t[fMaxDigit];
973 Float_t * tempoE = new Float_t[fMaxDigit];
092b50ba 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 ;
0fbb8738 987 fDigitsList = new Int_t[fMulDigit];
988 fEnergyList = new Float_t[fMulDigit];
092b50ba 989
0458aa58 990 fAmp = 0.; //Recalculate total energy
092b50ba 991 for(Int_t iDigit=0;iDigit< fMulDigit ;iDigit++){
0458aa58 992 fDigitsList[iDigit] = tempo[iDigit];
993 fEnergyList[iDigit] = tempoE[iDigit] ;
994 fAmp+=tempoE[iDigit];
092b50ba 995 }
996
997 delete [] tempo ;
998 delete [] tempoE ;
999
1000}
1001//____________________________________________________________________________
a8c47ab6 1002void AliPHOSEmcRecPoint::Print(Option_t *) const
d15a28e7 1003{
b2a60966 1004 // Print the list of digits belonging to the cluster
1005
21cd0c07 1006 TString message ;
1007 message = "AliPHOSEmcRecPoint:\n" ;
6f766ada 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() ) ;
d15a28e7 1013
d15a28e7 1014 Int_t iDigit;
6f766ada 1015 printf(" digits ids = ") ;
7932f811 1016 for(iDigit=0; iDigit<fMulDigit; iDigit++)
53ad38d8 1017 printf(" %d ", fDigitsList[iDigit] ) ;
7932f811 1018
6f766ada 1019 printf("\n digit energies = ") ;
7932f811 1020 for(iDigit=0; iDigit<fMulDigit; iDigit++)
53ad38d8 1021 printf(" %f ", fEnergyList[iDigit] ) ;
6f766ada 1022
1023 printf("\n digit primaries ") ;
1024 if (fMulTrack<1) printf("... no primaries");
bc68d12c 1025 for(iDigit = 0;iDigit < fMulTrack; iDigit++)
53ad38d8 1026 printf(" %d ", fTracksList[iDigit]) ;
1027 printf("\n") ;
6f766ada 1028
d15a28e7 1029}
88714635 1030
7932f811 1031