1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 /* History of cvs commits:
21 * Revision 1.22 2005/05/28 14:19:04 schutz
22 * Compilation warnings fixed by T.P.
26 //_________________________________________________________________________
27 // RecPoint implementation for PHOS-CPV
28 // An CpvRecPoint is a cluster of digits
29 //*-- Author: Yuri Kharlov
30 // (after Dmitri Peressounko (RRC KI & SUBATECH))
33 // --- ROOT system ---
36 #include "TClonesArray.h"
38 // --- Standard library ---
40 // --- AliRoot header files ---
42 #include "AliPHOSGeometry.h"
43 #include "AliPHOSDigit.h"
44 #include "AliPHOSCpvRecPoint.h"
45 #include "AliPHOSLoader.h"
47 ClassImp(AliPHOSCpvRecPoint)
49 //____________________________________________________________________________
50 AliPHOSCpvRecPoint::AliPHOSCpvRecPoint() : AliPHOSEmcRecPoint()
58 //____________________________________________________________________________
59 AliPHOSCpvRecPoint::AliPHOSCpvRecPoint(const char * opt) : AliPHOSEmcRecPoint(opt)
67 //____________________________________________________________________________
68 AliPHOSCpvRecPoint::~AliPHOSCpvRecPoint()
73 //____________________________________________________________________________
74 Bool_t AliPHOSCpvRecPoint::AreNeighbours(AliPHOSDigit * digit1, AliPHOSDigit * digit2 ) const
76 // Tells if (true) or not (false) two digits are neighbors)
78 Bool_t aren = kFALSE ;
80 AliPHOSGeometry * phosgeom = AliPHOSLoader::GetPHOSGeometry();
83 phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ;
86 phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ;
88 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
89 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
91 if (( coldiff <= 1 ) && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0))
97 //____________________________________________________________________________
98 Int_t AliPHOSCpvRecPoint::Compare(const TObject * obj) const
100 // Compares two RecPoints according to their position in the PHOS modules
102 Float_t delta = 1 ; //Width of "Sorting row". If you changibg this
103 //value (what is senseless) change as vell delta in
104 //AliPHOSTrackSegmentMakerv* and other RecPoints...
108 AliPHOSCpvRecPoint * clu = (AliPHOSCpvRecPoint *) obj ;
110 Int_t phosmod1 = GetPHOSMod() ;
111 Int_t phosmod2 = clu->GetPHOSMod() ;
114 GetLocalPosition(locpos1) ;
116 clu->GetLocalPosition(locpos2) ;
118 if(phosmod1 == phosmod2 ) {
119 Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
124 else if(locpos1.Z()>locpos2.Z())
131 if(phosmod1 < phosmod2 )
141 //______________________________________________________________________________
142 void AliPHOSCpvRecPoint::ExecuteEvent(Int_t, Int_t, Int_t ) /*const*/
144 // // Execute action corresponding to one event
145 // // This member function is called when a AliPHOSRecPoint is clicked with the locator
147 // // If Left button is clicked on AliPHOSRecPoint, the digits are switched on
148 // // and switched off when the mouse button is released.
151 // // static Int_t pxold, pyold;
153 // AliPHOSLoader * gime = AliPHOSLoader::GetInstance() ;
155 // static TGraph * digitgraph = 0 ;
157 // if (!gPad->IsEditable()) return;
159 // TH2F * histo = 0 ;
160 // TCanvas * histocanvas ;
164 // case kButton1Down: {
165 // AliPHOSDigit * digit ;
166 // AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
170 // const Int_t kMulDigit = AliPHOSCpvRecPoint::GetDigitsMultiplicity() ;
171 // Float_t * xi = new Float_t[kMulDigit] ;
172 // Float_t * zi = new Float_t[kMulDigit] ;
174 // // create the histogram for the single cluster
175 // // 1. gets histogram boundaries
176 // Float_t ximax = -999. ;
177 // Float_t zimax = -999. ;
178 // Float_t ximin = 999. ;
179 // Float_t zimin = 999. ;
181 // for(iDigit=0; iDigit<kMulDigit; iDigit++) {
182 // digit = (AliPHOSDigit *) ( gime->Digit(fDigitsList[iDigit]) ) ;
183 // phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
184 // phosgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
185 // if ( xi[iDigit] > ximax )
186 // ximax = xi[iDigit] ;
187 // if ( xi[iDigit] < ximin )
188 // ximin = xi[iDigit] ;
189 // if ( zi[iDigit] > zimax )
190 // zimax = zi[iDigit] ;
191 // if ( zi[iDigit] < zimin )
192 // zimin = zi[iDigit] ;
194 // ximax += phosgeom->GetCrystalSize(0) / 2. ;
195 // zimax += phosgeom->GetCrystalSize(2) / 2. ;
196 // ximin -= phosgeom->GetCrystalSize(0) / 2. ;
197 // zimin -= phosgeom->GetCrystalSize(2) / 2. ;
198 // Int_t xdim = (int)( (ximax - ximin ) / phosgeom->GetCrystalSize(0) + 0.5 ) ;
199 // Int_t zdim = (int)( (zimax - zimin ) / phosgeom->GetCrystalSize(2) + 0.5 ) ;
201 // // 2. gets the histogram title
203 // Text_t title[100] ;
204 // sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
210 // histo = new TH2F("cluster3D", title, xdim, ximin, ximax, zdim, zimin, zimax) ;
213 // for(iDigit=0; iDigit<kMulDigit; iDigit++) {
214 // digit = (AliPHOSDigit *) ( gime->Digit(fDigitsList[iDigit]) ) ;
215 // phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
216 // phosgeom->RelPosInModule(relid, x, z);
217 // histo->Fill(x, z, fEnergyList[iDigit] ) ;
220 // if (!digitgraph) {
221 // digitgraph = new TGraph(kMulDigit,xi,zi);
222 // digitgraph-> SetMarkerStyle(5) ;
223 // digitgraph-> SetMarkerSize(1.) ;
224 // digitgraph-> SetMarkerColor(1) ;
225 // digitgraph-> Paint("P") ;
229 // histocanvas = new TCanvas("cluser", "a single cluster", 600, 500) ;
230 // histocanvas->Draw() ;
231 // histo->Draw("lego1") ;
241 // delete digitgraph ;
249 //____________________________________________________________________________
250 void AliPHOSCpvRecPoint::EvalAll(Float_t logWeight,TClonesArray * digits)
252 // wraps other methods
253 AliPHOSEmcRecPoint::EvalAll(logWeight,digits) ;
254 EvalClusterLengths(digits) ;
256 //____________________________________________________________________________
257 void AliPHOSCpvRecPoint::EvalLocalPosition(Float_t logWeight,TClonesArray * digits)
259 // Calculates the center of gravity in the local PHOS-module coordinates
268 AliPHOSDigit * digit ;
270 AliPHOSGeometry * phosgeom = AliPHOSLoader::GetPHOSGeometry();
274 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
275 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]);
279 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
280 phosgeom->RelPosInModule(relid, xi, zi);
281 if (fAmp>0 && fEnergyList[iDigit]>0) {
282 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
288 AliError(Form("Wrong energy %f and/or amplitude %f\n", fEnergyList[iDigit], fAmp));
298 AliWarning(Form("Too low log weight factor to evaluate cluster's center" )) ;
307 //____________________________________________________________________________
308 void AliPHOSCpvRecPoint::EvalClusterLengths(TClonesArray * digits)
310 //Modified 15.03.2001 by Dmitri Peressounko
312 // Calculates the cluster lengths along X and Z axes
313 // These characteristics are needed for CPV to tune
314 // digitization+reconstruction to experimental data
315 // Yuri Kharlov. 24 October 2000
319 AliPHOSDigit * digit ;
321 AliPHOSGeometry * phosgeom = AliPHOSLoader::GetPHOSGeometry();
323 const Int_t kMaxLeng=20;
324 Int_t idX[kMaxLeng], idZ[kMaxLeng];
329 for(Int_t iDigit=0; iDigit<fMulDigit; iDigit++) {
330 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
331 Int_t absId = digit->GetId();
332 phosgeom->AbsToRelNumbering(absId, relid) ;
336 for (i=0; i<fLengX; i++) if (relid[2]==idX[i]) { dejavu=kTRUE; break; }
338 idX[fLengX]=relid[2];
340 fLengX = TMath::Min(fLengX,kMaxLeng);
344 for (i=0; i<fLengZ; i++) if (relid[3]==idZ[i]) { dejavu=kTRUE; break; }
346 idZ[fLengZ]=relid[3];
348 fLengZ = TMath::Min(fLengZ,kMaxLeng);
355 //____________________________________________________________________________
356 void AliPHOSCpvRecPoint::Print(const Option_t *) const
358 // Print the list of digits belonging to the cluster
361 message = "AliPHOSCpvRecPoint: " ;
362 message += "Digits # " ;
363 AliInfo(Form(message.Data())) ;
367 for(iDigit=0; iDigit<fMulDigit; iDigit++)
368 printf(" %d \n", fDigitsList[iDigit]) ;
370 printf("Energies: \n") ;
371 for(iDigit=0; iDigit<fMulDigit; iDigit++)
372 printf(" %f ", fEnergyList[iDigit]) ;
374 message = " Multiplicity = %d\n" ;
375 message += " Cluster Energy = %f\n" ;
376 message += " Stored at position %d\n" ;
378 printf(message.Data(), fMulDigit, fAmp, GetIndexInList() ) ;