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 ---
30 #include "TClonesArray.h"
32 // --- Standard library ---
36 // --- AliRoot header files ---
38 #include "AliPHOSCpvRecPoint.h"
39 #include "AliPHOSPpsdRecPoint.h"
41 ClassImp(AliPHOSCpvRecPoint)
43 //____________________________________________________________________________
44 AliPHOSCpvRecPoint::AliPHOSCpvRecPoint() : AliPHOSEmcRecPoint()
52 //____________________________________________________________________________
53 AliPHOSCpvRecPoint::~AliPHOSCpvRecPoint()
59 //____________________________________________________________________________
60 Bool_t AliPHOSCpvRecPoint::AreNeighbours(AliPHOSDigit * digit1, AliPHOSDigit * digit2 ) const
62 // Tells if (true) or not (false) two digits are neighbors)
64 Bool_t aren = kFALSE ;
66 AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
68 phosgeom->AbsToRelNumbering(digit1->GetId(), relid1) ;
71 phosgeom->AbsToRelNumbering(digit2->GetId(), relid2) ;
73 Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;
74 Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;
76 if (( coldiff <= 1 ) && ( rowdiff <= 1 ) && (coldiff + rowdiff > 0))
82 //____________________________________________________________________________
83 Int_t AliPHOSCpvRecPoint::Compare(const TObject * obj) const
85 // Compares two RecPoints according to their position in the PHOS modules
87 Float_t delta = 1 ; //Width of "Sorting row". If you changibg this
88 //value (what is senseless) change as vell delta in
89 //AliPHOSTrackSegmentMakerv* and other RecPoints...
93 if( (strcmp(obj->ClassName() , "AliPHOSPpsdRecPoint" )) == 0) // PPSD Rec Point
95 AliPHOSPpsdRecPoint * clu = (AliPHOSPpsdRecPoint *)obj ;
96 if(this->GetPHOSMod() < clu->GetPHOSMod() )
104 AliPHOSCpvRecPoint * clu = (AliPHOSCpvRecPoint *) obj ;
106 Int_t phosmod1 = GetPHOSMod() ;
107 Int_t phosmod2 = clu->GetPHOSMod() ;
110 GetLocalPosition(locpos1) ;
112 clu->GetLocalPosition(locpos2) ;
114 if(phosmod1 == phosmod2 ) {
115 Int_t rowdif = (Int_t)TMath::Ceil(locpos1.X()/delta)-(Int_t)TMath::Ceil(locpos2.X()/delta) ;
120 else if(locpos1.Z()>locpos2.Z())
127 if(phosmod1 < phosmod2 )
137 //______________________________________________________________________________
138 void AliPHOSCpvRecPoint::ExecuteEvent(Int_t event, Int_t px, Int_t py) const
140 // // Execute action corresponding to one event
141 // // This member function is called when a AliPHOSRecPoint is clicked with the locator
143 // // If Left button is clicked on AliPHOSRecPoint, the digits are switched on
144 // // and switched off when the mouse button is released.
147 // // static Int_t pxold, pyold;
149 // AliPHOSIndexToObject * please = AliPHOSIndexToObject::GetInstance() ;
151 // static TGraph * digitgraph = 0 ;
153 // if (!gPad->IsEditable()) return;
155 // TH2F * histo = 0 ;
156 // TCanvas * histocanvas ;
160 // case kButton1Down: {
161 // AliPHOSDigit * digit ;
162 // AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
166 // const Int_t kMulDigit = AliPHOSCpvRecPoint::GetDigitsMultiplicity() ;
167 // Float_t * xi = new Float_t[kMulDigit] ;
168 // Float_t * zi = new Float_t[kMulDigit] ;
170 // // create the histogram for the single cluster
171 // // 1. gets histogram boundaries
172 // Float_t ximax = -999. ;
173 // Float_t zimax = -999. ;
174 // Float_t ximin = 999. ;
175 // Float_t zimin = 999. ;
177 // for(iDigit=0; iDigit<kMulDigit; iDigit++) {
178 // digit = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigit]) ) ;
179 // phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
180 // phosgeom->RelPosInModule(relid, xi[iDigit], zi[iDigit]);
181 // if ( xi[iDigit] > ximax )
182 // ximax = xi[iDigit] ;
183 // if ( xi[iDigit] < ximin )
184 // ximin = xi[iDigit] ;
185 // if ( zi[iDigit] > zimax )
186 // zimax = zi[iDigit] ;
187 // if ( zi[iDigit] < zimin )
188 // zimin = zi[iDigit] ;
190 // ximax += phosgeom->GetCrystalSize(0) / 2. ;
191 // zimax += phosgeom->GetCrystalSize(2) / 2. ;
192 // ximin -= phosgeom->GetCrystalSize(0) / 2. ;
193 // zimin -= phosgeom->GetCrystalSize(2) / 2. ;
194 // Int_t xdim = (int)( (ximax - ximin ) / phosgeom->GetCrystalSize(0) + 0.5 ) ;
195 // Int_t zdim = (int)( (zimax - zimin ) / phosgeom->GetCrystalSize(2) + 0.5 ) ;
197 // // 2. gets the histogram title
199 // Text_t title[100] ;
200 // sprintf(title,"Energy=%1.2f GeV ; Digits ; %d ", GetEnergy(), GetDigitsMultiplicity()) ;
206 // histo = new TH2F("cluster3D", title, xdim, ximin, ximax, zdim, zimin, zimax) ;
209 // for(iDigit=0; iDigit<kMulDigit; iDigit++) {
210 // digit = (AliPHOSDigit *) ( please->GimeDigit(fDigitsList[iDigit]) ) ;
211 // phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
212 // phosgeom->RelPosInModule(relid, x, z);
213 // histo->Fill(x, z, fEnergyList[iDigit] ) ;
216 // if (!digitgraph) {
217 // digitgraph = new TGraph(kMulDigit,xi,zi);
218 // digitgraph-> SetMarkerStyle(5) ;
219 // digitgraph-> SetMarkerSize(1.) ;
220 // digitgraph-> SetMarkerColor(1) ;
221 // digitgraph-> Paint("P") ;
225 // histocanvas = new TCanvas("cluser", "a single cluster", 600, 500) ;
226 // histocanvas->Draw() ;
227 // histo->Draw("lego1") ;
237 // delete digitgraph ;
245 //____________________________________________________________________________
246 void AliPHOSCpvRecPoint::EvalAll(Float_t logWeight,TClonesArray * digits)
248 // wraps other methods
249 AliPHOSEmcRecPoint::EvalAll(logWeight,digits) ;
250 EvalClusterLengths(digits) ;
252 //____________________________________________________________________________
253 void AliPHOSCpvRecPoint::EvalLocalPosition(Float_t logWeight,TClonesArray * digits)
255 // Calculates the center of gravity in the local PHOS-module coordinates
264 AliPHOSDigit * digit ;
266 AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
270 for(iDigit=0; iDigit<fMulDigit; iDigit++) {
271 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]);
275 phosgeom->AbsToRelNumbering(digit->GetId(), relid) ;
276 phosgeom->RelPosInModule(relid, xi, zi);
277 Float_t w = TMath::Max( 0., logWeight + TMath::Log( fEnergyList[iDigit] / fAmp ) ) ;
289 if (fMulDigit != 0) cout << "AliPHOSCpvRecPoint: too low log weight factor "
290 << "to evaluate cluster's center\n";
298 //____________________________________________________________________________
299 void AliPHOSCpvRecPoint::EvalClusterLengths(TClonesArray * digits)
301 //Modified 15.03.2001 by Dmitri Peressounko
303 // Calculates the cluster lengths along X and Z axes
304 // These characteristics are needed for CPV to tune
305 // digitization+reconstruction to experimental data
306 // Yuri Kharlov. 24 October 2000
310 AliPHOSDigit * digit ;
312 AliPHOSGeometry * phosgeom = (AliPHOSGeometry *) fGeom ;
314 const Int_t kMaxLeng=20;
315 Int_t idX[kMaxLeng], idZ[kMaxLeng];
320 for(Int_t iDigit=0; iDigit<fMulDigit; iDigit++) {
321 digit = (AliPHOSDigit *) digits->At(fDigitsList[iDigit]) ;
322 Int_t absId = digit->GetId();
323 phosgeom->AbsToRelNumbering(absId, relid) ;
327 for (i=0; i<fLengX; i++) if (relid[2]==idX[i]) { dejavu=kTRUE; break; }
329 idX[fLengX]=relid[2];
331 fLengX = TMath::Min(fLengX,kMaxLeng);
335 for (i=0; i<fLengZ; i++) if (relid[3]==idZ[i]) { dejavu=kTRUE; break; }
337 idZ[fLengZ]=relid[3];
339 fLengZ = TMath::Min(fLengZ,kMaxLeng);
346 //____________________________________________________________________________
347 void AliPHOSCpvRecPoint::Print(Option_t * option)
349 // Print the list of digits belonging to the cluster
351 cout << "AliPHOSCpvRecPoint: " << endl ;
355 cout << "Digits # " ;
356 for(iDigit=0; iDigit<fMulDigit; iDigit++)
357 cout << fDigitsList[iDigit] << " " ;
360 cout << "Energies: " ;
361 for(iDigit=0; iDigit<fMulDigit; iDigit++)
362 cout << fEnergyList[iDigit] << " " ;
365 cout << " Multiplicity = " << fMulDigit << endl ;
366 cout << " Cluster Energy = " << fAmp << endl ;
367 cout << " Stored at position " << GetIndexInList() << endl ;