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 //_________________________________________________________________________
19 // RecPoint implementation for PHOS-CPV
20 // An CpvRecPoint is a cluster of digits
21 //*-- Author: Yuri Kharlov
22 // (after Dmitri Peressounko (RRC KI & SUBATECH))
25 // --- ROOT system ---
28 #include "TClonesArray.h"
30 // --- Standard library ---
32 // --- AliRoot header files ---
34 #include "AliPHOSGeometry.h"
35 #include "AliPHOSDigit.h"
36 #include "AliPHOSCpvRecPoint.h"
37 #include "AliPHOSLoader.h"
39 ClassImp(AliPHOSCpvRecPoint)
41 //____________________________________________________________________________
42 AliPHOSCpvRecPoint::AliPHOSCpvRecPoint() : AliPHOSEmcRecPoint()
50 //____________________________________________________________________________
51 AliPHOSCpvRecPoint::AliPHOSCpvRecPoint(const char * opt) : AliPHOSEmcRecPoint(opt)
59 //____________________________________________________________________________
60 AliPHOSCpvRecPoint::~AliPHOSCpvRecPoint()
65 //____________________________________________________________________________
66 Bool_t AliPHOSCpvRecPoint::AreNeighbours(AliPHOSDigit * digit1, AliPHOSDigit * digit2 ) const
68 // Tells if (true) or not (false) two digits are neighbors)
70 Bool_t aren = kFALSE ;
72 AliPHOSGeometry * phosgeom = AliPHOSLoader::GetPHOSGeometry();
75 phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ;
78 phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ;
80 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
81 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
83 if (( coldiff <= 1 ) && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0))
89 //____________________________________________________________________________
90 Int_t AliPHOSCpvRecPoint::Compare(const TObject * obj) const
92 // Compares two RecPoints according to their position in the PHOS modules
94 Float_t delta = 1 ; //Width of "Sorting row". If you changibg this
95 //value (what is senseless) change as vell delta in
96 //AliPHOSTrackSegmentMakerv* and other RecPoints...
100 AliPHOSCpvRecPoint * clu = (AliPHOSCpvRecPoint *) obj ;
102 Int_t phosmod1 = GetPHOSMod() ;
103 Int_t phosmod2 = clu->GetPHOSMod() ;
106 GetLocalPosition(locpos1) ;
108 clu->GetLocalPosition(locpos2) ;
110 if(phosmod1 == phosmod2 ) {
111 Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
116 else if(locpos1.Z()>locpos2.Z())
123 if(phosmod1 < phosmod2 )
133 //______________________________________________________________________________
134 void AliPHOSCpvRecPoint::ExecuteEvent(Int_t, Int_t, Int_t ) const
136 // // Execute action corresponding to one event
137 // // This member function is called when a AliPHOSRecPoint is clicked with the locator
139 // // If Left button is clicked on AliPHOSRecPoint, the digits are switched on
140 // // and switched off when the mouse button is released.
143 // // static Int_t pxold, pyold;
145 // AliPHOSLoader * gime = AliPHOSLoader::GetInstance() ;
147 // static TGraph * digitgraph = 0 ;
149 // if (!gPad->IsEditable()) return;
151 // TH2F * histo = 0 ;
152 // TCanvas * histocanvas ;
156 // case kButton1Down: {
157 // AliPHOSDigit * digit ;
158 // AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
162 // const Int_t kMulDigit = AliPHOSCpvRecPoint::GetDigitsMultiplicity() ;
163 // Float_t * xi = new Float_t[kMulDigit] ;
164 // Float_t * zi = new Float_t[kMulDigit] ;
166 // // create the histogram for the single cluster
167 // // 1. gets histogram boundaries
168 // Float_t ximax = -999. ;
169 // Float_t zimax = -999. ;
170 // Float_t ximin = 999. ;
171 // Float_t zimin = 999. ;
173 // for(iDigit=0; iDigit<kMulDigit; iDigit++) {
174 // digit = (AliPHOSDigit *) ( gime->Digit(fDigitsList[iDigit]) ) ;
175 // phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
176 // phosgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
177 // if ( xi[iDigit] > ximax )
178 // ximax = xi[iDigit] ;
179 // if ( xi[iDigit] < ximin )
180 // ximin = xi[iDigit] ;
181 // if ( zi[iDigit] > zimax )
182 // zimax = zi[iDigit] ;
183 // if ( zi[iDigit] < zimin )
184 // zimin = zi[iDigit] ;
186 // ximax += phosgeom->GetCrystalSize(0) / 2. ;
187 // zimax += phosgeom->GetCrystalSize(2) / 2. ;
188 // ximin -= phosgeom->GetCrystalSize(0) / 2. ;
189 // zimin -= phosgeom->GetCrystalSize(2) / 2. ;
190 // Int_t xdim = (int)( (ximax - ximin ) / phosgeom->GetCrystalSize(0) + 0.5 ) ;
191 // Int_t zdim = (int)( (zimax - zimin ) / phosgeom->GetCrystalSize(2) + 0.5 ) ;
193 // // 2. gets the histogram title
195 // Text_t title[100] ;
196 // sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
202 // histo = new TH2F("cluster3D", title, xdim, ximin, ximax, zdim, zimin, zimax) ;
205 // for(iDigit=0; iDigit<kMulDigit; iDigit++) {
206 // digit = (AliPHOSDigit *) ( gime->Digit(fDigitsList[iDigit]) ) ;
207 // phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
208 // phosgeom->RelPosInModule(relid, x, z);
209 // histo->Fill(x, z, fEnergyList[iDigit] ) ;
212 // if (!digitgraph) {
213 // digitgraph = new TGraph(kMulDigit,xi,zi);
214 // digitgraph-> SetMarkerStyle(5) ;
215 // digitgraph-> SetMarkerSize(1.) ;
216 // digitgraph-> SetMarkerColor(1) ;
217 // digitgraph-> Paint("P") ;
221 // histocanvas = new TCanvas("cluser", "a single cluster", 600, 500) ;
222 // histocanvas->Draw() ;
223 // histo->Draw("lego1") ;
233 // delete digitgraph ;
241 //____________________________________________________________________________
242 void AliPHOSCpvRecPoint::EvalAll(Float_t logWeight,TClonesArray * digits)
244 // wraps other methods
245 AliPHOSEmcRecPoint::EvalAll(logWeight,digits) ;
246 EvalClusterLengths(digits) ;
248 //____________________________________________________________________________
249 void AliPHOSCpvRecPoint::EvalLocalPosition(Float_t logWeight,TClonesArray * digits)
251 // Calculates the center of gravity in the local PHOS-module coordinates
260 AliPHOSDigit * digit ;
262 AliPHOSGeometry * phosgeom = AliPHOSLoader::GetPHOSGeometry();
266 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
267 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]);
271 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
272 phosgeom->RelPosInModule(relid, xi, zi);
273 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
286 AliWarning(Form("Too low log weight factor to evaluate cluster's center" )) ;
295 //____________________________________________________________________________
296 void AliPHOSCpvRecPoint::EvalClusterLengths(TClonesArray * digits)
298 //Modified 15.03.2001 by Dmitri Peressounko
300 // Calculates the cluster lengths along X and Z axes
301 // These characteristics are needed for CPV to tune
302 // digitization+reconstruction to experimental data
303 // Yuri Kharlov. 24 October 2000
307 AliPHOSDigit * digit ;
309 AliPHOSGeometry * phosgeom = AliPHOSLoader::GetPHOSGeometry();
311 const Int_t kMaxLeng=20;
312 Int_t idX[kMaxLeng], idZ[kMaxLeng];
317 for(Int_t iDigit=0; iDigit<fMulDigit; iDigit++) {
318 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
319 Int_t absId = digit->GetId();
320 phosgeom->AbsToRelNumbering(absId, relid) ;
324 for (i=0; i<fLengX; i++) if (relid[2]==idX[i]) { dejavu=kTRUE; break; }
326 idX[fLengX]=relid[2];
328 fLengX = TMath::Min(fLengX,kMaxLeng);
332 for (i=0; i<fLengZ; i++) if (relid[3]==idZ[i]) { dejavu=kTRUE; break; }
334 idZ[fLengZ]=relid[3];
336 fLengZ = TMath::Min(fLengZ,kMaxLeng);
343 //____________________________________________________________________________
344 void AliPHOSCpvRecPoint::Print()
346 // Print the list of digits belonging to the cluster
349 message = "AliPHOSCpvRecPoint: " ;
350 message += "Digits # " ;
351 AliInfo(Form(message.Data())) ;
355 for(iDigit=0; iDigit<fMulDigit; iDigit++)
356 printf(" %d \n", fDigitsList[iDigit]) ;
358 printf("Energies: \n") ;
359 for(iDigit=0; iDigit<fMulDigit; iDigit++)
360 printf(" %f ", fEnergyList[iDigit]) ;
362 message = " Multiplicity = %d\n" ;
363 message += " Cluster Energy = %f\n" ;
364 message += " Stored at position %d\n" ;
366 printf(message.Data(), fMulDigit, fAmp, GetIndexInList() ) ;