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.24 2006/08/28 10:01:56 kharlov
22 * Effective C++ warnings fixed (Timur Pocheptsov)
24 * Revision 1.23 2005/12/20 14:28:47 hristov
25 * Additional protection
27 * Revision 1.22 2005/05/28 14:19:04 schutz
28 * Compilation warnings fixed by T.P.
32 //_________________________________________________________________________
33 // RecPoint implementation for PHOS-CPV
34 // An CpvRecPoint is a cluster of digits
35 //*-- Author: Yuri Kharlov
36 // (after Dmitri Peressounko (RRC KI & SUBATECH))
39 // --- ROOT system ---
42 #include "TClonesArray.h"
44 // --- Standard library ---
46 // --- AliRoot header files ---
48 #include "AliPHOSGeometry.h"
49 #include "AliPHOSDigit.h"
50 #include "AliPHOSCpvRecPoint.h"
51 #include "AliPHOSLoader.h"
53 ClassImp(AliPHOSCpvRecPoint)
55 //____________________________________________________________________________
56 AliPHOSCpvRecPoint::AliPHOSCpvRecPoint() :
64 //____________________________________________________________________________
65 AliPHOSCpvRecPoint::AliPHOSCpvRecPoint(const char * opt) :
66 AliPHOSEmcRecPoint(opt),
73 //____________________________________________________________________________
74 AliPHOSCpvRecPoint::~AliPHOSCpvRecPoint()
79 //____________________________________________________________________________
80 Bool_t AliPHOSCpvRecPoint::AreNeighbours(AliPHOSDigit * digit1, AliPHOSDigit * digit2 ) const
82 // Tells if (true) or not (false) two digits are neighbors)
84 Bool_t aren = kFALSE ;
86 AliPHOSGeometry * phosgeom = AliPHOSLoader::GetPHOSGeometry();
89 phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ;
92 phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ;
94 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
95 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
97 if (( coldiff <= 1 ) && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0))
103 //____________________________________________________________________________
104 Int_t AliPHOSCpvRecPoint::Compare(const TObject * obj) const
106 // Compares two RecPoints according to their position in the PHOS modules
108 Float_t delta = 1 ; //Width of "Sorting row". If you changibg this
109 //value (what is senseless) change as vell delta in
110 //AliPHOSTrackSegmentMakerv* and other RecPoints...
114 AliPHOSCpvRecPoint * clu = (AliPHOSCpvRecPoint *) obj ;
116 Int_t phosmod1 = GetPHOSMod() ;
117 Int_t phosmod2 = clu->GetPHOSMod() ;
120 GetLocalPosition(locpos1) ;
122 clu->GetLocalPosition(locpos2) ;
124 if(phosmod1 == phosmod2 ) {
125 Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
130 else if(locpos1.Z()>locpos2.Z())
137 if(phosmod1 < phosmod2 )
147 //______________________________________________________________________________
148 void AliPHOSCpvRecPoint::ExecuteEvent(Int_t, Int_t, Int_t ) /*const*/
150 // // Execute action corresponding to one event
151 // // This member function is called when a AliPHOSRecPoint is clicked with the locator
153 // // If Left button is clicked on AliPHOSRecPoint, the digits are switched on
154 // // and switched off when the mouse button is released.
157 // // static Int_t pxold, pyold;
159 // AliPHOSLoader * gime = AliPHOSLoader::GetInstance() ;
161 // static TGraph * digitgraph = 0 ;
163 // if (!gPad->IsEditable()) return;
165 // TH2F * histo = 0 ;
166 // TCanvas * histocanvas ;
170 // case kButton1Down: {
171 // AliPHOSDigit * digit ;
172 // AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
176 // const Int_t kMulDigit = AliPHOSCpvRecPoint::GetDigitsMultiplicity() ;
177 // Float_t * xi = new Float_t[kMulDigit] ;
178 // Float_t * zi = new Float_t[kMulDigit] ;
180 // // create the histogram for the single cluster
181 // // 1. gets histogram boundaries
182 // Float_t ximax = -999. ;
183 // Float_t zimax = -999. ;
184 // Float_t ximin = 999. ;
185 // Float_t zimin = 999. ;
187 // for(iDigit=0; iDigit<kMulDigit; iDigit++) {
188 // digit = (AliPHOSDigit *) ( gime->Digit(fDigitsList[iDigit]) ) ;
189 // phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
190 // phosgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
191 // if ( xi[iDigit] > ximax )
192 // ximax = xi[iDigit] ;
193 // if ( xi[iDigit] < ximin )
194 // ximin = xi[iDigit] ;
195 // if ( zi[iDigit] > zimax )
196 // zimax = zi[iDigit] ;
197 // if ( zi[iDigit] < zimin )
198 // zimin = zi[iDigit] ;
200 // ximax += phosgeom->GetCrystalSize(0) / 2. ;
201 // zimax += phosgeom->GetCrystalSize(2) / 2. ;
202 // ximin -= phosgeom->GetCrystalSize(0) / 2. ;
203 // zimin -= phosgeom->GetCrystalSize(2) / 2. ;
204 // Int_t xdim = (int)( (ximax - ximin ) / phosgeom->GetCrystalSize(0) + 0.5 ) ;
205 // Int_t zdim = (int)( (zimax - zimin ) / phosgeom->GetCrystalSize(2) + 0.5 ) ;
207 // // 2. gets the histogram title
209 // Text_t title[100] ;
210 // sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
216 // histo = new TH2F("cluster3D", title, xdim, ximin, ximax, zdim, zimin, zimax) ;
219 // for(iDigit=0; iDigit<kMulDigit; iDigit++) {
220 // digit = (AliPHOSDigit *) ( gime->Digit(fDigitsList[iDigit]) ) ;
221 // phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
222 // phosgeom->RelPosInModule(relid, x, z);
223 // histo->Fill(x, z, fEnergyList[iDigit] ) ;
226 // if (!digitgraph) {
227 // digitgraph = new TGraph(kMulDigit,xi,zi);
228 // digitgraph-> SetMarkerStyle(5) ;
229 // digitgraph-> SetMarkerSize(1.) ;
230 // digitgraph-> SetMarkerColor(1) ;
231 // digitgraph-> Paint("P") ;
235 // histocanvas = new TCanvas("cluser", "a single cluster", 600, 500) ;
236 // histocanvas->Draw() ;
237 // histo->Draw("lego1") ;
247 // delete digitgraph ;
255 //____________________________________________________________________________
256 void AliPHOSCpvRecPoint::EvalAll(Float_t logWeight, TClonesArray * digits)
258 AliPHOSEmcRecPoint::EvalAll(logWeight,digits) ;
259 EvalClusterLengths(digits) ;
261 //____________________________________________________________________________
262 void AliPHOSCpvRecPoint::EvalAll(Float_t logWeight, TVector3 &vtx, TClonesArray * digits)
264 // wraps other methods
265 AliPHOSEmcRecPoint::EvalAll(logWeight,vtx,digits) ;
267 //____________________________________________________________________________
268 void AliPHOSCpvRecPoint::EvalLocalPosition(Float_t logWeight, TVector3 & /*vtx */, TClonesArray * digits)
270 // Calculates the center of gravity in the local PHOS-module coordinates
279 AliPHOSDigit * digit ;
281 AliPHOSGeometry * phosgeom = AliPHOSLoader::GetPHOSGeometry();
285 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
286 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]);
290 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
291 phosgeom->RelPosInModule(relid, xi, zi);
292 if (fAmp>0 && fEnergyList[iDigit]>0) {
293 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
299 AliError(Form("Wrong energy %f and/or amplitude %f\n", fEnergyList[iDigit], fAmp));
309 AliWarning(Form("Too low log weight factor to evaluate cluster's center" )) ;
318 //____________________________________________________________________________
319 void AliPHOSCpvRecPoint::EvalClusterLengths(TClonesArray * digits)
321 //Modified 15.03.2001 by Dmitri Peressounko
323 // Calculates the cluster lengths along X and Z axes
324 // These characteristics are needed for CPV to tune
325 // digitization+reconstruction to experimental data
326 // Yuri Kharlov. 24 October 2000
330 AliPHOSDigit * digit ;
332 AliPHOSGeometry * phosgeom = AliPHOSLoader::GetPHOSGeometry();
334 const Int_t kMaxLeng=20;
335 Int_t idX[kMaxLeng], idZ[kMaxLeng];
340 for(Int_t iDigit=0; iDigit<fMulDigit; iDigit++) {
341 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
342 Int_t absId = digit->GetId();
343 phosgeom->AbsToRelNumbering(absId, relid) ;
347 for (i=0; i<fLengX; i++) if (relid[2]==idX[i]) { dejavu=kTRUE; break; }
349 idX[fLengX]=relid[2];
351 fLengX = TMath::Min(fLengX,kMaxLeng);
355 for (i=0; i<fLengZ; i++) if (relid[3]==idZ[i]) { dejavu=kTRUE; break; }
357 idZ[fLengZ]=relid[3];
359 fLengZ = TMath::Min(fLengZ,kMaxLeng);
364 //____________________________________________________________________________
365 void AliPHOSCpvRecPoint::Print(const Option_t *) const
367 // Print the list of digits belonging to the cluster
370 message = "AliPHOSCpvRecPoint: " ;
371 message += "Digits # " ;
372 AliInfo(Form(message.Data())) ;
376 for(iDigit=0; iDigit<fMulDigit; iDigit++)
377 printf(" %d \n", fDigitsList[iDigit]) ;
379 printf("Energies: \n") ;
380 for(iDigit=0; iDigit<fMulDigit; iDigit++)
381 printf(" %f ", fEnergyList[iDigit]) ;
383 message = " Multiplicity = %d\n" ;
384 message += " Cluster Energy = %f\n" ;
385 message += " Stored at position %d\n" ;
387 printf(message.Data(), fMulDigit, fAmp, GetIndexInList() ) ;