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.25 2007/03/06 06:47:28 kharlov
22 * DP:Possibility to use actual vertex position added
24 * Revision 1.24 2006/08/28 10:01:56 kharlov
25 * Effective C++ warnings fixed (Timur Pocheptsov)
27 * Revision 1.23 2005/12/20 14:28:47 hristov
28 * Additional protection
30 * Revision 1.22 2005/05/28 14:19:04 schutz
31 * Compilation warnings fixed by T.P.
35 //_________________________________________________________________________
36 // RecPoint implementation for PHOS-CPV
37 // An CpvRecPoint is a cluster of digits
38 //*-- Author: Yuri Kharlov
39 // (after Dmitri Peressounko (RRC KI & SUBATECH))
42 // --- ROOT system ---
45 #include "TClonesArray.h"
47 // --- Standard library ---
49 // --- AliRoot header files ---
51 #include "AliPHOSGeometry.h"
52 #include "AliPHOSDigit.h"
53 #include "AliPHOSCpvRecPoint.h"
54 #include "AliPHOSLoader.h"
56 ClassImp(AliPHOSCpvRecPoint)
58 //____________________________________________________________________________
59 AliPHOSCpvRecPoint::AliPHOSCpvRecPoint() :
67 //____________________________________________________________________________
68 AliPHOSCpvRecPoint::AliPHOSCpvRecPoint(const char * opt) :
69 AliPHOSEmcRecPoint(opt),
76 //____________________________________________________________________________
77 AliPHOSCpvRecPoint::~AliPHOSCpvRecPoint()
82 //____________________________________________________________________________
83 Bool_t AliPHOSCpvRecPoint::AreNeighbours(AliPHOSDigit * digit1, AliPHOSDigit * digit2 ) const
85 // Tells if (true) or not (false) two digits are neighbors)
87 Bool_t aren = kFALSE ;
89 AliPHOSGeometry * phosgeom = AliPHOSLoader::GetPHOSGeometry();
92 phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ;
95 phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ;
97 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
98 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
100 if (( coldiff <= 1 ) && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0))
106 //____________________________________________________________________________
107 Int_t AliPHOSCpvRecPoint::Compare(const TObject * obj) const
109 // Compares two RecPoints according to their position in the PHOS modules
111 Float_t delta = 1 ; //Width of "Sorting row". If you changibg this
112 //value (what is senseless) change as vell delta in
113 //AliPHOSTrackSegmentMakerv* and other RecPoints...
117 AliPHOSCpvRecPoint * clu = (AliPHOSCpvRecPoint *) obj ;
119 Int_t phosmod1 = GetPHOSMod() ;
120 Int_t phosmod2 = clu->GetPHOSMod() ;
123 GetLocalPosition(locpos1) ;
125 clu->GetLocalPosition(locpos2) ;
127 if(phosmod1 == phosmod2 ) {
128 Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
133 else if(locpos1.Z()>locpos2.Z())
140 if(phosmod1 < phosmod2 )
150 //______________________________________________________________________________
151 void AliPHOSCpvRecPoint::ExecuteEvent(Int_t, Int_t, Int_t ) /*const*/
153 // // Execute action corresponding to one event
154 // // This member function is called when a AliPHOSRecPoint is clicked with the locator
156 // // If Left button is clicked on AliPHOSRecPoint, the digits are switched on
157 // // and switched off when the mouse button is released.
160 // // static Int_t pxold, pyold;
162 // AliPHOSLoader * gime = AliPHOSLoader::GetInstance() ;
164 // static TGraph * digitgraph = 0 ;
166 // if (!gPad->IsEditable()) return;
168 // TH2F * histo = 0 ;
169 // TCanvas * histocanvas ;
173 // case kButton1Down: {
174 // AliPHOSDigit * digit ;
175 // AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
179 // const Int_t kMulDigit = AliPHOSCpvRecPoint::GetDigitsMultiplicity() ;
180 // Float_t * xi = new Float_t[kMulDigit] ;
181 // Float_t * zi = new Float_t[kMulDigit] ;
183 // // create the histogram for the single cluster
184 // // 1. gets histogram boundaries
185 // Float_t ximax = -999. ;
186 // Float_t zimax = -999. ;
187 // Float_t ximin = 999. ;
188 // Float_t zimin = 999. ;
190 // for(iDigit=0; iDigit<kMulDigit; iDigit++) {
191 // digit = (AliPHOSDigit *) ( gime->Digit(fDigitsList[iDigit]) ) ;
192 // phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
193 // phosgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
194 // if ( xi[iDigit] > ximax )
195 // ximax = xi[iDigit] ;
196 // if ( xi[iDigit] < ximin )
197 // ximin = xi[iDigit] ;
198 // if ( zi[iDigit] > zimax )
199 // zimax = zi[iDigit] ;
200 // if ( zi[iDigit] < zimin )
201 // zimin = zi[iDigit] ;
203 // ximax += phosgeom->GetCrystalSize(0) / 2. ;
204 // zimax += phosgeom->GetCrystalSize(2) / 2. ;
205 // ximin -= phosgeom->GetCrystalSize(0) / 2. ;
206 // zimin -= phosgeom->GetCrystalSize(2) / 2. ;
207 // Int_t xdim = (int)( (ximax - ximin ) / phosgeom->GetCrystalSize(0) + 0.5 ) ;
208 // Int_t zdim = (int)( (zimax - zimin ) / phosgeom->GetCrystalSize(2) + 0.5 ) ;
210 // // 2. gets the histogram title
212 // Text_t title[100] ;
213 // sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
219 // histo = new TH2F("cluster3D", title, xdim, ximin, ximax, zdim, zimin, zimax) ;
222 // for(iDigit=0; iDigit<kMulDigit; iDigit++) {
223 // digit = (AliPHOSDigit *) ( gime->Digit(fDigitsList[iDigit]) ) ;
224 // phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
225 // phosgeom->RelPosInModule(relid, x, z);
226 // histo->Fill(x, z, fEnergyList[iDigit] ) ;
229 // if (!digitgraph) {
230 // digitgraph = new TGraph(kMulDigit,xi,zi);
231 // digitgraph-> SetMarkerStyle(5) ;
232 // digitgraph-> SetMarkerSize(1.) ;
233 // digitgraph-> SetMarkerColor(1) ;
234 // digitgraph-> Paint("P") ;
238 // histocanvas = new TCanvas("cluser", "a single cluster", 600, 500) ;
239 // histocanvas->Draw() ;
240 // histo->Draw("lego1") ;
250 // delete digitgraph ;
258 //____________________________________________________________________________
259 void AliPHOSCpvRecPoint::EvalAll(Float_t logWeight, TClonesArray * digits)
261 AliPHOSEmcRecPoint::EvalAll(logWeight,digits) ;
262 EvalClusterLengths(digits) ;
264 //____________________________________________________________________________
265 void AliPHOSCpvRecPoint::EvalAll(Float_t logWeight, TVector3 &vtx, TClonesArray * digits)
267 // wraps other methods
268 AliPHOSEmcRecPoint::EvalAll(logWeight,vtx,digits) ;
270 //____________________________________________________________________________
271 void AliPHOSCpvRecPoint::EvalLocalPosition(Float_t logWeight, TVector3 & /*vtx */, TClonesArray * digits, TVector3 &/* vInc */)
273 // Calculates the center of gravity in the local PHOS-module coordinates
282 AliPHOSDigit * digit ;
284 AliPHOSGeometry * phosgeom = AliPHOSLoader::GetPHOSGeometry();
288 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
289 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]);
293 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
294 phosgeom->RelPosInModule(relid, xi, zi);
295 if (fAmp>0 && fEnergyList[iDigit]>0) {
296 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
302 AliError(Form("Wrong energy %f and/or amplitude %f\n", fEnergyList[iDigit], fAmp));
312 AliWarning(Form("Too low log weight factor to evaluate cluster's center" )) ;
321 //____________________________________________________________________________
322 void AliPHOSCpvRecPoint::EvalClusterLengths(TClonesArray * digits)
324 //Modified 15.03.2001 by Dmitri Peressounko
326 // Calculates the cluster lengths along X and Z axes
327 // These characteristics are needed for CPV to tune
328 // digitization+reconstruction to experimental data
329 // Yuri Kharlov. 24 October 2000
333 AliPHOSDigit * digit ;
335 AliPHOSGeometry * phosgeom = AliPHOSLoader::GetPHOSGeometry();
337 const Int_t kMaxLeng=20;
338 Int_t idX[kMaxLeng], idZ[kMaxLeng];
343 for(Int_t iDigit=0; iDigit<fMulDigit; iDigit++) {
344 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
345 Int_t absId = digit->GetId();
346 phosgeom->AbsToRelNumbering(absId, relid) ;
350 for (i=0; i<fLengX; i++) if (relid[2]==idX[i]) { dejavu=kTRUE; break; }
352 idX[fLengX]=relid[2];
354 fLengX = TMath::Min(fLengX,kMaxLeng);
358 for (i=0; i<fLengZ; i++) if (relid[3]==idZ[i]) { dejavu=kTRUE; break; }
360 idZ[fLengZ]=relid[3];
362 fLengZ = TMath::Min(fLengZ,kMaxLeng);
367 //____________________________________________________________________________
368 void AliPHOSCpvRecPoint::Print(const Option_t *) const
370 // Print the list of digits belonging to the cluster
373 message = "AliPHOSCpvRecPoint: " ;
374 message += "Digits # " ;
375 AliInfo(Form(message.Data())) ;
379 for(iDigit=0; iDigit<fMulDigit; iDigit++)
380 printf(" %d \n", fDigitsList[iDigit]) ;
382 printf("Energies: \n") ;
383 for(iDigit=0; iDigit<fMulDigit; iDigit++)
384 printf(" %f ", fEnergyList[iDigit]) ;
386 message = " Multiplicity = %d\n" ;
387 message += " Cluster Energy = %f\n" ;
388 message += " Stored at position %d\n" ;
390 printf(message.Data(), fMulDigit, fAmp, GetIndexInList() ) ;