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.26 2007/06/18 07:02:44 kharlov
22 * Change the signature of EvalLocalPosition() to obey the method virtuality from the parent class
24 * Revision 1.25 2007/03/06 06:47:28 kharlov
25 * DP:Possibility to use actual vertex position added
27 * Revision 1.24 2006/08/28 10:01:56 kharlov
28 * Effective C++ warnings fixed (Timur Pocheptsov)
30 * Revision 1.23 2005/12/20 14:28:47 hristov
31 * Additional protection
33 * Revision 1.22 2005/05/28 14:19:04 schutz
34 * Compilation warnings fixed by T.P.
38 //_________________________________________________________________________
39 // RecPoint implementation for PHOS-CPV
40 // An CpvRecPoint is a cluster of digits
41 //*-- Author: Yuri Kharlov
42 // (after Dmitri Peressounko (RRC KI & SUBATECH))
45 // --- ROOT system ---
48 #include "TClonesArray.h"
50 // --- Standard library ---
52 // --- AliRoot header files ---
54 #include "AliPHOSGeometry.h"
55 #include "AliPHOSDigit.h"
56 #include "AliPHOSCpvRecPoint.h"
57 #include "AliPHOSLoader.h"
59 ClassImp(AliPHOSCpvRecPoint)
61 //____________________________________________________________________________
62 AliPHOSCpvRecPoint::AliPHOSCpvRecPoint() :
70 //____________________________________________________________________________
71 AliPHOSCpvRecPoint::AliPHOSCpvRecPoint(const char * opt) :
72 AliPHOSEmcRecPoint(opt),
79 //____________________________________________________________________________
80 AliPHOSCpvRecPoint::~AliPHOSCpvRecPoint()
85 //____________________________________________________________________________
86 Bool_t AliPHOSCpvRecPoint::AreNeighbours(AliPHOSDigit * digit1, AliPHOSDigit * digit2 ) const
88 // Tells if (true) or not (false) two digits are neighbors)
90 Bool_t aren = kFALSE ;
92 AliPHOSGeometry * phosgeom = AliPHOSLoader::GetPHOSGeometry();
95 phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ;
98 phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ;
100 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
101 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
103 if (( coldiff <= 1 ) && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0))
109 //____________________________________________________________________________
110 Int_t AliPHOSCpvRecPoint::Compare(const TObject * obj) const
112 // Compares two RecPoints according to their position in the PHOS modules
114 Float_t delta = 1 ; //Width of "Sorting row". If you changibg this
115 //value (what is senseless) change as vell delta in
116 //AliPHOSTrackSegmentMakerv* and other RecPoints...
120 AliPHOSCpvRecPoint * clu = (AliPHOSCpvRecPoint *) obj ;
122 Int_t phosmod1 = GetPHOSMod() ;
123 Int_t phosmod2 = clu->GetPHOSMod() ;
126 GetLocalPosition(locpos1) ;
128 clu->GetLocalPosition(locpos2) ;
130 if(phosmod1 == phosmod2 ) {
131 Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
136 else if(locpos1.Z()>locpos2.Z())
143 if(phosmod1 < phosmod2 )
153 //______________________________________________________________________________
154 void AliPHOSCpvRecPoint::ExecuteEvent(Int_t, Int_t, Int_t ) /*const*/
156 // // Execute action corresponding to one event
157 // // This member function is called when a AliPHOSRecPoint is clicked with the locator
159 // // If Left button is clicked on AliPHOSRecPoint, the digits are switched on
160 // // and switched off when the mouse button is released.
163 // // static Int_t pxold, pyold;
165 // AliPHOSLoader * gime = AliPHOSLoader::GetInstance() ;
167 // static TGraph * digitgraph = 0 ;
169 // if (!gPad->IsEditable()) return;
171 // TH2F * histo = 0 ;
172 // TCanvas * histocanvas ;
176 // case kButton1Down: {
177 // AliPHOSDigit * digit ;
178 // AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
182 // const Int_t kMulDigit = AliPHOSCpvRecPoint::GetDigitsMultiplicity() ;
183 // Float_t * xi = new Float_t[kMulDigit] ;
184 // Float_t * zi = new Float_t[kMulDigit] ;
186 // // create the histogram for the single cluster
187 // // 1. gets histogram boundaries
188 // Float_t ximax = -999. ;
189 // Float_t zimax = -999. ;
190 // Float_t ximin = 999. ;
191 // Float_t zimin = 999. ;
193 // for(iDigit=0; iDigit<kMulDigit; iDigit++) {
194 // digit = (AliPHOSDigit *) ( gime->Digit(fDigitsList[iDigit]) ) ;
195 // phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
196 // phosgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
197 // if ( xi[iDigit] > ximax )
198 // ximax = xi[iDigit] ;
199 // if ( xi[iDigit] < ximin )
200 // ximin = xi[iDigit] ;
201 // if ( zi[iDigit] > zimax )
202 // zimax = zi[iDigit] ;
203 // if ( zi[iDigit] < zimin )
204 // zimin = zi[iDigit] ;
206 // ximax += phosgeom->GetCrystalSize(0) / 2. ;
207 // zimax += phosgeom->GetCrystalSize(2) / 2. ;
208 // ximin -= phosgeom->GetCrystalSize(0) / 2. ;
209 // zimin -= phosgeom->GetCrystalSize(2) / 2. ;
210 // Int_t xdim = (int)( (ximax - ximin ) / phosgeom->GetCrystalSize(0) + 0.5 ) ;
211 // Int_t zdim = (int)( (zimax - zimin ) / phosgeom->GetCrystalSize(2) + 0.5 ) ;
213 // // 2. gets the histogram title
215 // Text_t title[100] ;
216 // sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
222 // histo = new TH2F("cluster3D", title, xdim, ximin, ximax, zdim, zimin, zimax) ;
225 // for(iDigit=0; iDigit<kMulDigit; iDigit++) {
226 // digit = (AliPHOSDigit *) ( gime->Digit(fDigitsList[iDigit]) ) ;
227 // phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
228 // phosgeom->RelPosInModule(relid, x, z);
229 // histo->Fill(x, z, fEnergyList[iDigit] ) ;
232 // if (!digitgraph) {
233 // digitgraph = new TGraph(kMulDigit,xi,zi);
234 // digitgraph-> SetMarkerStyle(5) ;
235 // digitgraph-> SetMarkerSize(1.) ;
236 // digitgraph-> SetMarkerColor(1) ;
237 // digitgraph-> Paint("P") ;
241 // histocanvas = new TCanvas("cluser", "a single cluster", 600, 500) ;
242 // histocanvas->Draw() ;
243 // histo->Draw("lego1") ;
253 // delete digitgraph ;
261 //____________________________________________________________________________
262 void AliPHOSCpvRecPoint::EvalAll(Float_t logWeight, TClonesArray * digits)
264 // Evaluate local coordinate assuming the vertex in (000) and no inclination
265 TVector3 vtx(0,0,0), vInc(0,1,0);
266 AliPHOSEmcRecPoint::EvalAll(logWeight,digits) ;
267 EvalLocalPosition(logWeight, vtx, digits,vInc) ;
268 EvalClusterLengths(digits) ;
270 //____________________________________________________________________________
271 void AliPHOSCpvRecPoint::EvalAll(Float_t logWeight, TVector3 &vtx, TClonesArray * digits)
273 // wraps other methods
274 AliPHOSEmcRecPoint::EvalAll(logWeight,vtx,digits) ;
276 //____________________________________________________________________________
277 void AliPHOSCpvRecPoint::EvalLocalPosition(Float_t logWeight, TVector3 & /*vtx */, TClonesArray * digits, TVector3 &/* vInc */)
279 // Calculates the center of gravity in the local PHOS-module coordinates
288 AliPHOSDigit * digit ;
290 AliPHOSGeometry * phosgeom = AliPHOSLoader::GetPHOSGeometry();
294 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
295 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]);
299 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
300 phosgeom->RelPosInModule(relid, xi, zi);
301 if (fAmp>0 && fEnergyList[iDigit]>0) {
302 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
308 AliError(Form("Wrong energy %f and/or amplitude %f\n", fEnergyList[iDigit], fAmp));
318 AliWarning(Form("Too low log weight factor to evaluate cluster's center" )) ;
327 //____________________________________________________________________________
328 void AliPHOSCpvRecPoint::EvalClusterLengths(TClonesArray * digits)
330 //Modified 15.03.2001 by Dmitri Peressounko
332 // Calculates the cluster lengths along X and Z axes
333 // These characteristics are needed for CPV to tune
334 // digitization+reconstruction to experimental data
335 // Yuri Kharlov. 24 October 2000
339 AliPHOSDigit * digit ;
341 AliPHOSGeometry * phosgeom = AliPHOSLoader::GetPHOSGeometry();
343 const Int_t kMaxLeng=20;
344 Int_t idX[kMaxLeng], idZ[kMaxLeng];
349 for(Int_t iDigit=0; iDigit<fMulDigit; iDigit++) {
350 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
351 Int_t absId = digit->GetId();
352 phosgeom->AbsToRelNumbering(absId, relid) ;
356 for (i=0; i<fLengX; i++) if (relid[2]==idX[i]) { dejavu=kTRUE; break; }
358 idX[fLengX]=relid[2];
360 fLengX = TMath::Min(fLengX,kMaxLeng);
364 for (i=0; i<fLengZ; i++) if (relid[3]==idZ[i]) { dejavu=kTRUE; break; }
366 idZ[fLengZ]=relid[3];
368 fLengZ = TMath::Min(fLengZ,kMaxLeng);
373 //____________________________________________________________________________
374 void AliPHOSCpvRecPoint::Print(const Option_t *) const
376 // Print the list of digits belonging to the cluster
379 message = "AliPHOSCpvRecPoint: " ;
380 message += "Digits # " ;
381 AliInfo(Form(message.Data())) ;
385 for(iDigit=0; iDigit<fMulDigit; iDigit++)
386 printf(" %d \n", fDigitsList[iDigit]) ;
388 printf("Energies: \n") ;
389 for(iDigit=0; iDigit<fMulDigit; iDigit++)
390 printf(" %f ", fEnergyList[iDigit]) ;
392 message = " Multiplicity = %d\n" ;
393 message += " Cluster Energy = %f\n" ;
394 message += " Stored at position %d\n" ;
396 printf(message.Data(), fMulDigit, fAmp, GetIndexInList() ) ;