New Version from Alla
[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
16//_________________________________________________________________________
17// Rec Point in the PHOS EM calorimeter
18//*-- Author : Dmitri Peressounko RRC KI
19//////////////////////////////////////////////////////////////////////////////
20
21// --- ROOT system ---
9f616d61 22#include "TPad.h"
23#include "TH2.h"
d15a28e7 24#include "TMath.h"
9f616d61 25#include "TCanvas.h"
d15a28e7 26
27// --- Standard library ---
28
9f616d61 29#include <iostream>
d15a28e7 30
31// --- AliRoot header files ---
32
33#include "AliPHOSGeometry.h"
34#include "AliPHOSEmcRecPoint.h"
35#include "AliRun.h"
36
37ClassImp(AliPHOSEmcRecPoint)
38
39//____________________________________________________________________________
40AliPHOSEmcRecPoint::AliPHOSEmcRecPoint(Float_t W0, Float_t LocMaxCut)
41 : AliPHOSRecPoint()
42{
43 // ctor
44
45 fMulDigit = 0 ;
46 fAmp = 0. ;
47 fEnergyList = new Float_t[fMaxDigit];
48 AliPHOSGeometry * PHOSGeom = (AliPHOSGeometry *) fGeom ;
49 fDelta = PHOSGeom->GetCrystalSize(0) ;
50 fW0 = W0 ;
51 fLocMaxCut = LocMaxCut ;
52 fLocPos.SetX(1000000.) ; //Local position should be evaluated
53}
54
9f616d61 55//____________________________________________________________________________
56AliPHOSEmcRecPoint::~AliPHOSEmcRecPoint()
57{
58 // dtor
59}
d15a28e7 60
61//____________________________________________________________________________
62void AliPHOSEmcRecPoint::AddDigit(AliDigitNew & digit, Float_t Energy)
63{
64 // adds a digit to the digits list
65 // and accumulates the total amplitude and the multiplicity
66
67 if ( fMulDigit >= fMaxDigit ) { // increase the size of the lists
9f616d61 68 fMaxDigit*=2 ;
69 int * tempo = new ( int[fMaxDigit] ) ;
70 Float_t * tempoE = new ( Float_t[fMaxDigit] ) ;
71
72 Int_t index ;
d15a28e7 73 for ( index = 0 ; index < fMulDigit ; index++ ){
74 tempo[index] = fDigitsList[index] ;
75 tempoE[index] = fEnergyList[index] ;
76 }
77
9f616d61 78 delete [] fDigitsList ;
79 fDigitsList = new ( int[fMaxDigit] ) ;
80
81 delete [] fEnergyList ;
82 fEnergyList = new ( Float_t[fMaxDigit] ) ;
83
84 for ( index = 0 ; index < fMulDigit ; index++ ){
85 fDigitsList[index] = tempo[index] ;
86 fEnergyList[index] = tempoE[index] ;
87 }
88
89 delete [] tempo ;
90 delete [] tempoE ;
91 } // if
d15a28e7 92
9f616d61 93 fDigitsList[fMulDigit] = (int) &digit ;
94 fEnergyList[fMulDigit++] = Energy ;
d15a28e7 95 fAmp += Energy ;
96}
97
98//____________________________________________________________________________
99Bool_t AliPHOSEmcRecPoint::AreNeighbours(AliPHOSDigit * digit1, AliPHOSDigit * digit2 )
100{
101
102 Bool_t aren = kFALSE ;
103
104 AliPHOSGeometry * PHOSGeom = (AliPHOSGeometry *) fGeom ;
105 Int_t relid1[4] ;
106 PHOSGeom->AbsToRelNumbering(digit1->GetId(), relid1) ;
107
108 Int_t relid2[4] ;
109 PHOSGeom->AbsToRelNumbering(digit2->GetId(), relid2) ;
110
111 Int_t RowDiff = TMath::Abs( relid1[2] - relid2[2] ) ;
112 Int_t ColDiff = TMath::Abs( relid1[3] - relid2[3] ) ;
113
114 if (( ColDiff<=1 ) && ( RowDiff <= 1 ) && (ColDiff+RowDiff > 0))
115 aren = kTRUE ;
116
117 return aren ;
118}
119
120//____________________________________________________________________________
121Int_t AliPHOSEmcRecPoint::Compare(TObject * obj)
122{
123 Int_t rv ;
124
125 AliPHOSEmcRecPoint * clu = (AliPHOSEmcRecPoint *)obj ;
126
127
128 Int_t PHOSMod1 = this->GetPHOSMod() ;
129 Int_t PHOSMod2 = clu->GetPHOSMod() ;
130
131 TVector3 LocPos1;
132 this->GetLocalPosition(LocPos1) ;
133 TVector3 LocPos2;
134 clu->GetLocalPosition(LocPos2) ;
135
136 if(PHOSMod1 == PHOSMod2 ) {
137 Int_t rowdif = (Int_t)TMath::Ceil(LocPos1.X()/fDelta)-(Int_t)TMath::Ceil(LocPos2.X()/fDelta) ;
138 if (rowdif> 0)
139 rv = -1 ;
140 else if(rowdif < 0)
141 rv = 1 ;
142 else if(LocPos1.Z()>LocPos2.Z())
143 rv = -1 ;
144 else
145 rv = 1 ;
146 }
147
148 else {
149 if(PHOSMod1 < PHOSMod2 )
150 rv = -1 ;
151 else
152 rv = 1 ;
153 }
154
155 return rv ;
156}
157
9f616d61 158//______________________________________________________________________________
159void AliPHOSEmcRecPoint::ExecuteEvent(Int_t event, Int_t px, Int_t py)
160{
161//*-*-*-*-*-*-*-*-*-*-*Execute action corresponding to one event*-*-*-*
162//*-* =========================================
163// This member function is called when a AliPHOSRecPoint is clicked with the locator
164//
165// If Left button is clicked on AliPHOSRecPoint, the digits are switched on
166// and switched off when the mouse button is released.
167//
168
169 // static Int_t pxold, pyold;
170
171 static TGraph * DigitGraph = 0 ;
172
173 if (!gPad->IsEditable()) return;
174
175 TH2F * Histo = 0 ;
176 TCanvas * HistoCanvas ;
177
178 switch (event) {
179
180 case kButton1Down: {
181 AliPHOSDigit * digit ;
182 AliPHOSGeometry * PHOSGeom = (AliPHOSGeometry *) fGeom ;
183 Int_t iDigit;
184 Int_t relid[4] ;
185 Float_t xi[fMulDigit] ;
186 Float_t zi[fMulDigit] ;
187
188 // create the histogram for the single cluster
189 // 1. gets histogram boundaries
190 Float_t ximax = -999. ;
191 Float_t zimax = -999. ;
192 Float_t ximin = 999. ;
193 Float_t zimin = 999. ;
194
195 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
196 digit = (AliPHOSDigit *) fDigitsList[iDigit];
197 PHOSGeom->AbsToRelNumbering(digit->GetId(), relid) ;
198 PHOSGeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
199 if ( xi[iDigit] > ximax )
200 ximax = xi[iDigit] ;
201 if ( xi[iDigit] < ximin )
202 ximin = xi[iDigit] ;
203 if ( zi[iDigit] > zimax )
204 zimax = zi[iDigit] ;
205 if ( zi[iDigit] < zimin )
206 zimin = zi[iDigit] ;
207 }
208 ximax += PHOSGeom->GetCrystalSize(0) / 2. ;
209 zimax += PHOSGeom->GetCrystalSize(2) / 2. ;
210 ximin -= PHOSGeom->GetCrystalSize(0) / 2. ;
211 zimin -= PHOSGeom->GetCrystalSize(2) / 2. ;
212 Int_t xdim = (int)( (ximax - ximin ) / PHOSGeom->GetCrystalSize(0) + 0.5 ) ;
213 Int_t zdim = (int)( (zimax - zimin ) / PHOSGeom->GetCrystalSize(2) + 0.5 ) ;
214
215 // 2. gets the histogram title
216
217 Text_t title[100] ;
218 sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
219
220 if (!Histo) {
221 delete Histo ;
222 Histo = 0 ;
223 }
224 Histo = new TH2F("cluster3D", title, xdim, ximin, ximax, zdim, zimin, zimax) ;
225
226 Float_t x, z ;
227 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
228 digit = (AliPHOSDigit *) fDigitsList[iDigit];
229 PHOSGeom->AbsToRelNumbering(digit->GetId(), relid) ;
230 PHOSGeom->RelPosInModule(relid, x, z);
231 Histo->Fill(x, z, fEnergyList[iDigit] ) ;
232 }
233
234 if (!DigitGraph) {
235 DigitGraph = new TGraph(fMulDigit,xi,zi);
236 DigitGraph-> SetMarkerStyle(5) ;
237 DigitGraph-> SetMarkerSize(1.) ;
238 DigitGraph-> SetMarkerColor(1) ;
239 DigitGraph-> Paint("P") ;
240 }
241
242 Print() ;
243 HistoCanvas = new TCanvas("cluser", "a single cluster", 600, 500) ;
244 HistoCanvas->Draw() ;
245 Histo->Draw("lego1") ;
246
247 break;
248 }
249
250 case kButton1Up:
251 if (DigitGraph) {
252 delete DigitGraph ;
253 DigitGraph = 0 ;
254 }
255 break;
256
257 }
258}
259
d15a28e7 260//____________________________________________________________________________
261Float_t AliPHOSEmcRecPoint::GetDispersion()
262{
263 Float_t D = 0 ;
264 Float_t wtot = 0 ;
265
266 TVector3 LocPos;
267 GetLocalPosition(LocPos);
268 Float_t x = LocPos.X() ;
269 Float_t z = LocPos.Z() ;
270 // Int_t i = GetPHOSMod() ;
271
272 AliPHOSDigit * digit ;
273 AliPHOSGeometry * PHOSGeom = (AliPHOSGeometry *) fGeom ;
274
275 Int_t iDigit;
276 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
277 digit = (AliPHOSDigit *) fDigitsList[iDigit];
278 Int_t relid[4] ;
279 Float_t xi ;
280 Float_t zi ;
281 PHOSGeom->AbsToRelNumbering(digit->GetId(), relid) ;
282 PHOSGeom->RelPosInModule(relid, xi, zi);
283 Float_t w = TMath::Max(0.,fW0+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
284 D += w*((xi-x)*(xi-x) + (zi-z)*(zi-z) ) ;
285 wtot+=w ;
286 }
287
288 D /= wtot ;
289
290 return TMath::Sqrt(D) ;
291}
292
293//____________________________________________________________________________
294void AliPHOSEmcRecPoint::GetElipsAxis(Float_t * lambda)
295{
296 Float_t wtot = 0. ;
297 Float_t x = 0.;
298 Float_t z = 0.;
299 Float_t Dxx = 0.;
300 Float_t Dzz = 0.;
301 Float_t Dxz = 0.;
302
303 AliPHOSDigit * digit ;
304 AliPHOSGeometry * PHOSGeom = (AliPHOSGeometry *) fGeom ;
305 Int_t iDigit;
306
307 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
308 digit = (AliPHOSDigit *) fDigitsList[iDigit];
309 Int_t relid[4] ;
310 Float_t xi ;
311 Float_t zi ;
312 PHOSGeom->AbsToRelNumbering(digit->GetId(), relid) ;
313 PHOSGeom->RelPosInModule(relid, xi, zi);
314 Float_t w = TMath::Max(0.,fW0+TMath::Log(fEnergyList[iDigit]/fAmp ) ) ;
315 Dxx += w * xi * xi ;
316 x += w * xi ;
317 Dzz += w * zi * zi ;
318 z += w * zi ;
319 Dxz += w * xi * zi ;
320 wtot += w ;
321 }
322
323 Dxx /= wtot ;
324 x /= wtot ;
325 Dxx -= x * x ;
326 Dzz /= wtot ;
327 z /= wtot ;
328 Dzz -= z * z ;
329 Dxz /= wtot ;
330 Dxz -= x * z ;
331
332 lambda[0] = TMath::Sqrt( 0.5 * (Dxx + Dzz) + TMath::Sqrt( 0.25 * (Dxx - Dzz) * (Dxx - Dzz) + Dxz * Dxz ) ) ;
333 lambda[1] = TMath::Sqrt( 0.5 * (Dxx + Dzz) - TMath::Sqrt( 0.25 * (Dxx - Dzz) * (Dxx - Dzz) + Dxz * Dxz ) ) ;
334}
335
336//____________________________________________________________________________
337Float_t AliPHOSEmcRecPoint::GetMaximalEnergy(void)
338{
339 Float_t menergy = 0. ;
340
341 Int_t iDigit;
342
343 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
344
345 if(fEnergyList[iDigit] > menergy)
346 menergy = fEnergyList[iDigit] ;
347 }
348 return menergy ;
349}
350
351//____________________________________________________________________________
352Int_t AliPHOSEmcRecPoint::GetMultiplicityAtLevel(Float_t H)
353{
354 Int_t multipl = 0 ;
355 Int_t iDigit ;
356 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
357
358 if(fEnergyList[iDigit] > H * fAmp)
359 multipl++ ;
360 }
361 return multipl ;
362}
363
364//____________________________________________________________________________
365Int_t AliPHOSEmcRecPoint::GetNumberOfLocalMax(int * maxAt, Float_t * maxAtEnergy)
366{
367 AliPHOSDigit * digit ;
368 AliPHOSDigit * digitN ;
369
370
371 Int_t iDigitN ;
372 Int_t iDigit ;
373
9f616d61 374 for(iDigit=0; iDigit<fMulDigit; iDigit++){
d15a28e7 375 maxAt[iDigit] = fDigitsList[iDigit] ;
376 }
377
378 for(iDigit=0 ; iDigit<fMulDigit; iDigit++) {
9f616d61 379 if(maxAt[iDigit] != -1) {
380 digit = (AliPHOSDigit *) maxAt[iDigit] ;
d15a28e7 381
382 for(iDigitN=0; iDigitN<fMulDigit; iDigitN++) {
9f616d61 383 digitN = (AliPHOSDigit *) fDigitsList[iDigitN] ;
d15a28e7 384
9f616d61 385 if ( AreNeighbours(digit, digitN) ) {
d15a28e7 386 if (fEnergyList[iDigit] > fEnergyList[iDigitN] ) {
387 maxAt[iDigitN] = -1 ;
388 //but may be digit is not local max too ?
389 if(fEnergyList[iDigit]<fEnergyList[iDigitN]+fLocMaxCut)
390 maxAt[iDigit] = -1 ;
391 }
392 else {
393 maxAt[iDigit] = -1 ;
394 //but may be digitN is not local max too ?
395 if(fEnergyList[iDigit] >fEnergyList[iDigitN]-fLocMaxCut)
396 maxAt[iDigitN] = -1 ;
397 }
398 } // if Areneighbours
399 } // while digitN
400 } // slot not empty
401 } // while digit
402
403 iDigitN = 0 ;
404 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
405 if(maxAt[iDigit] != -1){
406 maxAt[iDigitN] = maxAt[iDigit] ;
9f616d61 407 maxAtEnergy[iDigitN] = fEnergyList[iDigit] ;
408 iDigitN++ ;
d15a28e7 409 }
410 }
411 return iDigitN ;
412}
413
414//____________________________________________________________________________
415void AliPHOSEmcRecPoint::GetLocalPosition(TVector3 &LPos)
416{
417 if( fLocPos.X() < 1000000.) { //allready evaluated
418 LPos = fLocPos ;
419 return ;
420 }
421
422 Float_t wtot = 0. ;
423
424 Int_t relid[4] ;
425
426 Float_t x = 0. ;
427 Float_t z = 0. ;
428
429 AliPHOSDigit * digit ;
430
431 AliPHOSGeometry * PHOSGeom = (AliPHOSGeometry *) fGeom ;
432
433 Int_t iDigit;
434
435
436 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
437 digit = (AliPHOSDigit *) fDigitsList[iDigit];
438
439 Float_t xi ;
440 Float_t zi ;
441 PHOSGeom->AbsToRelNumbering(digit->GetId(), relid) ;
442 PHOSGeom->RelPosInModule(relid, xi, zi);
443 Float_t w = TMath::Max( 0., fW0 + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
444 x += xi * w ;
445 z += zi * w ;
446 wtot += w ;
447
448 }
449
450 x /= wtot ;
451 z /= wtot ;
452 fLocPos.SetX(x) ;
453 fLocPos.SetY(0.) ;
454 fLocPos.SetZ(z) ;
455
456 LPos = fLocPos ;
457}
458
459// //____________________________________________________________________________
460// AliPHOSEmcRecPoint& AliPHOSEmcRecPoint::operator = (AliPHOSEmcRecPoint Clu)
461// {
462// int * DL = Clu.GetDigitsList() ;
463
464// if(fDigitsList)
465// delete fDigitsList ;
466
467// AliPHOSDigit * digit ;
468
469// Int_t iDigit;
470
471// for(iDigit=0; iDigit<fMulDigit; iDigit++) {
472// digit = (AliPHOSDigit *) DL[iDigit];
473// AddDigit(*digit) ;
474// }
475
476// fAmp = Clu.GetTotalEnergy() ;
477// fGeom = Clu.GetGeom() ;
478// TVector3 LocPos;
479// Clu.GetLocalPosition(LocPos) ;
480// fLocPos = LocPos;
481// fMulDigit = Clu.GetMultiplicity() ;
482// fMaxDigit = Clu.GetMaximumMultiplicity() ;
483// fPHOSMod = Clu.GetPHOSMod() ;
484// fW0 = Clu.GetLogWeightCut() ;
485// fDelta = Clu.GetDelta() ;
486// fLocMaxCut = Clu.GetLocMaxCut() ;
487
488// delete DL ;
489
490// return *this ;
491// }
492
493//____________________________________________________________________________
494void AliPHOSEmcRecPoint::Print(Option_t * option)
495{
496 cout << "AliPHOSEmcRecPoint: " << endl ;
497
498 AliPHOSDigit * digit ;
499 Int_t iDigit;
9f616d61 500 AliPHOSGeometry * PHOSGeom = (AliPHOSGeometry *) fGeom ;
d15a28e7 501
9f616d61 502 Float_t xi ;
503 Float_t zi ;
504 Int_t relid[4] ;
505
d15a28e7 506 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
507 digit = (AliPHOSDigit *) fDigitsList[iDigit];
9f616d61 508 PHOSGeom->AbsToRelNumbering(digit->GetId(), relid) ;
509 PHOSGeom->RelPosInModule(relid, xi, zi);
510 cout << " Id = " << digit->GetId() ;
511 cout << " module = " << relid[0] ;
512 cout << " x = " << xi ;
513 cout << " z = " << zi ;
514 cout << " Energy = " << fEnergyList[iDigit] << endl ;
d15a28e7 515 }
516 cout << " Multiplicity = " << fMulDigit << endl ;
517 cout << " Cluster Energy = " << fAmp << endl ;
518
519}
520
521//______________________________________________________________________________
522void AliPHOSEmcRecPoint::Streamer(TBuffer &R__b)
523{
524 // Stream an object of class AliPHOSEmcRecPoint.
525
526 if (R__b.IsReading()) {
527 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
528 AliPHOSRecPoint::Streamer(R__b);
529 R__b >> fDelta;
530 R__b >> fLocMaxCut;
531 R__b.ReadArray(fEnergyList);
532 R__b >> fW0;
533 } else {
534 R__b.WriteVersion(AliPHOSEmcRecPoint::IsA());
535 AliPHOSRecPoint::Streamer(R__b);
536 R__b << fDelta;
537 R__b << fLocMaxCut;
538 R__b.WriteArray(fEnergyList, GetMaximumDigitMultiplicity() );
539 R__b << fW0;
540 }
541}
542
543