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
28 //_________________________________________________________________________
29 // RecPoint implementation for PHOS-CPV
30 // An CpvRecPoint is a cluster of digits
31 //-- Author: Yuri Kharlov
32 // (after Dmitri Peressounko (RRC KI & SUBATECH))
35 // --- ROOT system ---
38 #include "TClonesArray.h"
40 // --- Standard library ---
42 // --- AliRoot header files ---
44 #include "AliPHOSGeometry.h"
45 #include "AliPHOSDigit.h"
46 #include "AliPHOSCpvRecPoint.h"
47 #include "AliPHOSLoader.h"
49 ClassImp(AliPHOSCpvRecPoint)
51 //____________________________________________________________________________
52 AliPHOSCpvRecPoint::AliPHOSCpvRecPoint() :
60 //____________________________________________________________________________
61 AliPHOSCpvRecPoint::AliPHOSCpvRecPoint(const char * opt) :
62 AliPHOSEmcRecPoint(opt),
69 //____________________________________________________________________________
70 AliPHOSCpvRecPoint::~AliPHOSCpvRecPoint()
75 //____________________________________________________________________________
76 Bool_t AliPHOSCpvRecPoint::AreNeighbours(AliPHOSDigit * digit1, AliPHOSDigit * digit2 ) const
78 // Tells if (true) or not (false) two digits are neighbors)
80 Bool_t aren = kFALSE ;
82 AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance() ;
85 phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ;
88 phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ;
90 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
91 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
93 if (( coldiff <= 1 ) && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0))
99 //____________________________________________________________________________
100 Int_t AliPHOSCpvRecPoint::Compare(const TObject * obj) const
102 // Compares two RecPoints according to their position in the PHOS modules
104 Float_t delta = 1 ; //Width of "Sorting row". If you changibg this
105 //value (what is senseless) change as vell delta in
106 //AliPHOSTrackSegmentMakerv* and other RecPoints...
110 AliPHOSCpvRecPoint * clu = (AliPHOSCpvRecPoint *) obj ;
112 Int_t phosmod1 = GetPHOSMod() ;
113 Int_t phosmod2 = clu->GetPHOSMod() ;
116 GetLocalPosition(locpos1) ;
118 clu->GetLocalPosition(locpos2) ;
120 if(phosmod1 == phosmod2 ) {
121 Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
126 else if(locpos1.Z()>locpos2.Z())
133 if(phosmod1 < phosmod2 )
143 //______________________________________________________________________________
144 void AliPHOSCpvRecPoint::ExecuteEvent(Int_t, Int_t, Int_t ) /*const*/
146 // // Execute action corresponding to one event
147 // // This member function is called when a AliPHOSRecPoint is clicked with the locator
149 // // If Left button is clicked on AliPHOSRecPoint, the digits are switched on
150 // // and switched off when the mouse button is released.
153 // // static Int_t pxold, pyold;
155 // AliPHOSLoader * gime = AliPHOSLoader::GetInstance() ;
157 // static TGraph * digitgraph = 0 ;
159 // if (!gPad->IsEditable()) return;
161 // TH2F * histo = 0 ;
162 // TCanvas * histocanvas ;
166 // case kButton1Down: {
167 // AliPHOSDigit * digit ;
168 // AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
172 // const Int_t kMulDigit = AliPHOSCpvRecPoint::GetDigitsMultiplicity() ;
173 // Float_t * xi = new Float_t[kMulDigit] ;
174 // Float_t * zi = new Float_t[kMulDigit] ;
176 // // create the histogram for the single cluster
177 // // 1. gets histogram boundaries
178 // Float_t ximax = -999. ;
179 // Float_t zimax = -999. ;
180 // Float_t ximin = 999. ;
181 // Float_t zimin = 999. ;
183 // for(iDigit=0; iDigit<kMulDigit; iDigit++) {
184 // digit = (AliPHOSDigit *) ( gime->Digit(fDigitsList[iDigit]) ) ;
185 // phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
186 // phosgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
187 // if ( xi[iDigit] > ximax )
188 // ximax = xi[iDigit] ;
189 // if ( xi[iDigit] < ximin )
190 // ximin = xi[iDigit] ;
191 // if ( zi[iDigit] > zimax )
192 // zimax = zi[iDigit] ;
193 // if ( zi[iDigit] < zimin )
194 // zimin = zi[iDigit] ;
196 // ximax += phosgeom->GetCrystalSize(0) / 2. ;
197 // zimax += phosgeom->GetCrystalSize(2) / 2. ;
198 // ximin -= phosgeom->GetCrystalSize(0) / 2. ;
199 // zimin -= phosgeom->GetCrystalSize(2) / 2. ;
200 // Int_t xdim = (int)( (ximax - ximin ) / phosgeom->GetCrystalSize(0) + 0.5 ) ;
201 // Int_t zdim = (int)( (zimax - zimin ) / phosgeom->GetCrystalSize(2) + 0.5 ) ;
203 // // 2. gets the histogram title
205 // Text_t title[100] ;
206 // sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
212 // histo = new TH2F("cluster3D", title, xdim, ximin, ximax, zdim, zimin, zimax) ;
215 // for(iDigit=0; iDigit<kMulDigit; iDigit++) {
216 // digit = (AliPHOSDigit *) ( gime->Digit(fDigitsList[iDigit]) ) ;
217 // phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
218 // phosgeom->RelPosInModule(relid, x, z);
219 // histo->Fill(x, z, fEnergyList[iDigit] ) ;
222 // if (!digitgraph) {
223 // digitgraph = new TGraph(kMulDigit,xi,zi);
224 // digitgraph-> SetMarkerStyle(5) ;
225 // digitgraph-> SetMarkerSize(1.) ;
226 // digitgraph-> SetMarkerColor(1) ;
227 // digitgraph-> Paint("P") ;
231 // histocanvas = new TCanvas("cluser", "a single cluster", 600, 500) ;
232 // histocanvas->Draw() ;
233 // histo->Draw("lego1") ;
243 // delete digitgraph ;
251 //____________________________________________________________________________
252 void AliPHOSCpvRecPoint::EvalAll(TClonesArray * digits)
254 // Evaluate local coordinate assuming the vertex in (000) and no inclination
255 AliPHOSEmcRecPoint::EvalAll(digits) ;
257 //____________________________________________________________________________
258 void AliPHOSCpvRecPoint::EvalAll(Float_t logWeight, TVector3 &vtx, TClonesArray * digits)
260 // wraps other methods
261 TVector3 vInc(0,1,0);
262 AliPHOSEmcRecPoint::EvalAll(logWeight,vtx,digits) ;
263 EvalLocalPosition(logWeight, vtx, digits,vInc) ;
264 EvalClusterLengths(digits) ;
266 //____________________________________________________________________________
267 void AliPHOSCpvRecPoint::EvalLocalPosition(Float_t logWeight, TVector3 & /*vtx */, TClonesArray * digits, TVector3 &/* vInc */)
269 // Calculates the center of gravity in the local PHOS-module coordinates
278 AliPHOSDigit * digit ;
280 AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance();
284 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
285 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]);
289 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
290 phosgeom->RelPosInModule(relid, xi, zi);
291 if (fAmp>0 && fEnergyList[iDigit]>0) {
292 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
298 AliError(Form("Wrong energy %f and/or amplitude %f\n", fEnergyList[iDigit], fAmp));
308 AliWarning(Form("Too low log weight factor to evaluate cluster's center" )) ;
316 //____________________________________________________________________________
317 void AliPHOSCpvRecPoint::EvalClusterLengths(TClonesArray * digits)
319 //Modified 15.03.2001 by Dmitri Peressounko
321 // Calculates the cluster lengths along X and Z axes
322 // These characteristics are needed for CPV to tune
323 // digitization+reconstruction to experimental data
324 // Yuri Kharlov. 24 October 2000
328 AliPHOSDigit * digit ;
330 AliPHOSGeometry * phosgeom = AliPHOSGeometry::GetInstance();
332 const Int_t kMaxLeng=20;
333 Int_t idX[kMaxLeng], idZ[kMaxLeng];
338 for(Int_t iDigit=0; iDigit<fMulDigit; iDigit++) {
339 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
340 Int_t absId = digit->GetId();
341 phosgeom->AbsToRelNumbering(absId, relid) ;
345 for (i=0; i<fLengX; i++) if (relid[2]==idX[i]) { dejavu=kTRUE; break; }
347 idX[fLengX]=relid[2];
349 fLengX = TMath::Min(fLengX,kMaxLeng);
353 for (i=0; i<fLengZ; i++) if (relid[3]==idZ[i]) { dejavu=kTRUE; break; }
355 idZ[fLengZ]=relid[3];
357 fLengZ = TMath::Min(fLengZ,kMaxLeng);
362 //____________________________________________________________________________
363 void AliPHOSCpvRecPoint::Print(const Option_t *) const
365 // Print the list of digits belonging to the cluster
368 message = "AliPHOSCpvRecPoint: " ;
369 message += "Digits # " ;
370 AliInfo(message.Data()) ;
374 for(iDigit=0; iDigit<fMulDigit; iDigit++)
375 printf(" %d \n", fDigitsList[iDigit]) ;
377 printf("Energies: \n") ;
378 for(iDigit=0; iDigit<fMulDigit; iDigit++)
379 printf(" %f ", fEnergyList[iDigit]) ;
381 message = " Multiplicity = %d\n" ;
382 message += " Cluster Energy = %f\n" ;
383 message += " Stored at position %d\n" ;
385 printf(message.Data(), fMulDigit, fAmp, GetIndexInList() ) ;