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 Revision 1.4.2.2 2000/05/08 14:50:58 cblume
19 Add functions ProjRow(), ProjCol(), and ProjTime()
21 Revision 1.4.2.1 2000/04/27 12:47:02 cblume
22 Replace Fill3() by Fill()
24 Revision 1.4 2000/02/28 19:10:26 cblume
25 Include the new TRD classes
27 Revision 1.3.4.1 2000/02/28 17:57:47 cblume
28 GetTrack returns now -1 if no track is found
30 Revision 1.3 1999/10/04 14:48:07 fca
31 Avoid warnings on non-ansi compiler HP-UX CC
33 Revision 1.2 1999/09/29 09:24:35 fca
34 Introduction of the Copyright and cvs Log
38 ///////////////////////////////////////////////////////////////////////////////
40 // Contains the pixel information for one TRD chamber //
42 ///////////////////////////////////////////////////////////////////////////////
44 #include "AliTRDmatrix.h"
46 ClassImp(AliTRDmatrix)
48 //_____________________________________________________________________________
49 AliTRDmatrix::AliTRDmatrix():TObject()
52 // Create a TRD detector matrix
66 //_____________________________________________________________________________
67 AliTRDmatrix::AliTRDmatrix(Int_t nRow, Int_t nCol, Int_t nTime
68 , Int_t iSec, Int_t iCha, Int_t iPla)
72 // Create a TRD detector matrix with a given size
78 fPixel = nRow * nCol * nTime;
82 fPixelArray = new TObjArray(fPixel);
83 for (Int_t iPixel = 0; iPixel < fPixel; iPixel++) {
84 AliTRDpixel *pixel = new AliTRDpixel();
85 fPixelArray->Add(pixel);
90 //_____________________________________________________________________________
91 AliTRDmatrix::~AliTRDmatrix()
95 fPixelArray->Delete();
101 //_____________________________________________________________________________
102 void AliTRDmatrix::AddSignal(Int_t iRow, Int_t iCol, Int_t iTime, Float_t signal)
105 // Add a value to the amplitude of the signal for one specific pixel
108 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
110 signal += pixel->GetSignal();
111 pixel->SetSignal(signal);
116 //_____________________________________________________________________________
117 void AliTRDmatrix::Draw(Option_t *)
120 // Draws a 3D view of the detector matrix
124 sprintf(ctitle,"Matrix (Sector:%d Chamber:%d Plane:%d)"
125 ,fSector,fChamber,fPlane);
126 TH3F *hMatrix = new TH3F("hMatrix",ctitle,fRow ,-0.5,fRow +0.5
127 ,fCol ,-0.5,fCol +0.5
128 ,fTime,-0.5,fTime+0.5);
130 for (Int_t iRow = 0; iRow < fRow; iRow++ ) {
131 for (Int_t iCol = 0; iCol < fCol; iCol++ ) {
132 for (Int_t iTime = 0; iTime < fTime; iTime++) {
133 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
134 if (pixel) hMatrix->Fill(iRow,iCol,iTime,pixel->GetSignal());
139 gStyle->SetOptStat(0);
140 TCanvas *cMatrix = new TCanvas("cMatrix","Detector matrix 3D-view"
142 cMatrix->ToggleEventStatus();
143 hMatrix->SetXTitle("Pad-row (z)");
144 hMatrix->SetYTitle("Pad-column (rphi)");
145 hMatrix->SetZTitle("Timebucket");
146 hMatrix->Draw("BOX");
150 //_____________________________________________________________________________
151 void AliTRDmatrix::DrawRow(Int_t iRow)
154 // Draws a 2D slice of the detector matrix along one row
157 if ((iRow < 0) || (iRow >= fRow)) {
158 printf("AliTRDmatrix::DrawRow -- ");
159 printf("Index out of bounds (%d/%d)\n",iRow,fRow);
164 sprintf(ctitle,"Pad-row %d (Sector:%d Chamber:%d Plane:%d)"
165 ,iRow,fSector,fChamber,fPlane);
166 TH2F *hSliceRow = new TH2F("hSliceRow",ctitle,fCol ,-0.5,fCol +0.5
167 ,fTime,-0.5,fTime+0.5);
169 for (Int_t iCol = 0; iCol < fCol; iCol++ ) {
170 for (Int_t iTime = 0; iTime < fTime; iTime++) {
171 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
172 if (pixel) hSliceRow->Fill(iCol,iTime,pixel->GetSignal());
176 gStyle->SetOptStat(0);
177 TCanvas *cSliceRow = new TCanvas("cSliceRow","Detector matrix 2D-slice"
179 cSliceRow->ToggleEventStatus();
180 hSliceRow->SetXTitle("Pad-column (rphi)");
181 hSliceRow->SetYTitle("Timebucket");
182 hSliceRow->Draw("COLZ");
186 //_____________________________________________________________________________
187 void AliTRDmatrix::DrawCol(Int_t iCol)
190 // Draws a 2D slice of the detector matrix along one column
193 if ((iCol < 0) || (iCol >= fCol)) {
194 printf("AliTRDmatrix::DrawCol -- ");
195 printf("Index out of bounds (%d/%d)\n",iCol,fCol);
200 sprintf(ctitle,"Pad-column %d (Sector:%d Chamber:%d Plane:%d)"
201 ,iCol,fSector,fChamber,fPlane);
202 TH2F *hSliceCol = new TH2F("hSliceCol",ctitle,fRow ,-0.5,fRow +0.5
203 ,fTime,-0.5,fTime+0.5);
205 for (Int_t iRow = 0; iRow < fRow; iRow++ ) {
206 for (Int_t iTime = 0; iTime < fTime; iTime++) {
207 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
208 if (pixel) hSliceCol->Fill(iRow,iTime,pixel->GetSignal());
212 gStyle->SetOptStat(0);
213 TCanvas *cSliceCol = new TCanvas("cSliceCol","Detector matrix 2D-slice"
215 cSliceCol->ToggleEventStatus();
216 hSliceCol->SetXTitle("Pad-row (z)");
217 hSliceCol->SetYTitle("Timebucket");
218 hSliceCol->Draw("COLZ");
222 //_____________________________________________________________________________
223 void AliTRDmatrix::DrawTime(Int_t iTime)
226 // Draws a 2D slice of the detector matrix along one time slice
229 if ((iTime < 0) || (iTime >= fTime)) {
230 printf("AliTRDmatrix::DrawTime -- ");
231 printf("Index out of bounds (%d/%d)\n",iTime,fTime);
236 sprintf(ctitle,"Time-slice %d (Sector:%d Chamber:%d Plane:%d)"
237 ,iTime,fSector,fChamber,fPlane);
238 TH2F *hSliceTime = new TH2F("hSliceTime",ctitle,fRow,-0.5,fRow+0.5
239 ,fCol,-0.5,fCol+0.5);
241 for (Int_t iRow = 0; iRow < fRow; iRow++) {
242 for (Int_t iCol = 0; iCol < fCol; iCol++) {
243 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
244 if (pixel) hSliceTime->Fill(iRow,iCol,pixel->GetSignal());
248 gStyle->SetOptStat(0);
249 TCanvas *cSliceTime = new TCanvas("cSliceTime","Detector matrix 2D-slice"
251 cSliceTime->ToggleEventStatus();
252 hSliceTime->SetXTitle("Pad-row (z)");
253 hSliceTime->SetYTitle("Pad-column (rphi)");
254 hSliceTime->Draw("COLZ");
258 //_____________________________________________________________________________
259 void AliTRDmatrix::ProjRow()
262 // Projects the detector matrix along the row-axis
266 sprintf(ctitle,"Row-projection (Sector:%d Chamber:%d Plane:%d)"
267 ,fSector,fChamber,fPlane);
268 TH2F *hProjRow = new TH2F("hProjRow",ctitle,fCol ,-0.5,fCol +0.5
269 ,fTime,-0.5,fTime+0.5);
271 for (Int_t iRow = 0; iRow < fRow; iRow++ ) {
272 for (Int_t iCol = 0; iCol < fCol; iCol++ ) {
273 for (Int_t iTime = 0; iTime < fTime; iTime++) {
274 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
275 if (pixel) hProjRow->Fill(iCol,iTime,pixel->GetSignal());
280 gStyle->SetOptStat(0);
281 TCanvas *cProjRow = new TCanvas("cProjRow","Detector matrix 2D-projection"
283 cProjRow->ToggleEventStatus();
284 hProjRow->SetXTitle("Pad-column (rphi)");
285 hProjRow->SetYTitle("Timebucket");
286 hProjRow->Draw("COLZ");
290 //_____________________________________________________________________________
291 void AliTRDmatrix::ProjCol()
294 // Projects the detector matrix along the column-axis
298 sprintf(ctitle,"Column-projection (Sector:%d Chamber:%d Plane:%d)"
299 ,fSector,fChamber,fPlane);
300 TH2F *hProjCol = new TH2F("hProjCol",ctitle,fRow ,-0.5,fRow +0.5
301 ,fTime,-0.5,fTime+0.5);
303 for (Int_t iRow = 0; iRow < fRow; iRow++ ) {
304 for (Int_t iCol = 0; iCol < fCol; iCol++ ) {
305 for (Int_t iTime = 0; iTime < fTime; iTime++) {
306 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
307 if (pixel) hProjCol->Fill(iRow,iTime,pixel->GetSignal());
312 gStyle->SetOptStat(0);
313 TCanvas *cProjCol = new TCanvas("cProjCol","Detector matrix 2D-projection"
315 cProjCol->ToggleEventStatus();
316 hProjCol->SetXTitle("Pad-row (z)");
317 hProjCol->SetYTitle("Timebucket");
318 hProjCol->Draw("COLZ");
322 //_____________________________________________________________________________
323 void AliTRDmatrix::ProjTime()
326 // Projects the detector matrix along the time-axis
330 sprintf(ctitle,"Time-projection (Sector:%d Chamber:%d Plane:%d)"
331 ,fSector,fChamber,fPlane);
332 TH2F *hProjTime = new TH2F("hProjTime",ctitle,fRow,-0.5,fRow+0.5
333 ,fCol,-0.5,fCol+0.5);
335 for (Int_t iRow = 0; iRow < fRow; iRow++) {
336 for (Int_t iCol = 0; iCol < fCol; iCol++) {
337 for (Int_t iTime = 0; iTime < fTime; iTime++) {
338 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
339 if (pixel) hProjTime->Fill(iRow,iCol,pixel->GetSignal());
344 gStyle->SetOptStat(0);
345 TCanvas *cProjTime = new TCanvas("cProjTime","Detector matrix 2D-projection"
347 cProjTime->ToggleEventStatus();
348 hProjTime->SetXTitle("Pad-row (z)");
349 hProjTime->SetYTitle("Pad-column (rphi)");
350 hProjTime->Draw("COLZ");
354 //_____________________________________________________________________________
355 void AliTRDmatrix::SetSignal(Int_t iRow, Int_t iCol, Int_t iTime, Float_t signal)
358 // Set the amplitude of the signal for one specific pixel
361 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
363 pixel->SetSignal(signal);
368 //_____________________________________________________________________________
369 Bool_t AliTRDmatrix::AddTrack(Int_t iRow, Int_t iCol, Int_t iTime, Int_t track)
372 // Add this track number to the stored tracks passing through this pixel.
373 // If there are already three stored the return status is FALSE.
376 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
377 if (!(pixel)) return kTRUE;
379 Bool_t trackSet = kFALSE;
380 for (Int_t i = 0; i < kTrackPixel; i++) {
381 if (pixel->GetTrack(i) == track) {
385 if (pixel->GetTrack(i) == -1) {
386 pixel->SetTrack(i,track);
396 //_____________________________________________________________________________
397 void AliTRDmatrix::SetTrack(Int_t iRow, Int_t iCol, Int_t iTime
398 , Int_t iTrack, Int_t track)
401 // Store the number of a track which is passing through this pixel
404 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
406 pixel->SetTrack(iTrack,track);
411 //_____________________________________________________________________________
412 Float_t AliTRDmatrix::GetSignal(Int_t iRow, Int_t iCol, Int_t iTime)
415 // Returns the amplitude of the signal for one specific pixel
418 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
420 return (pixel->GetSignal());
428 //_____________________________________________________________________________
429 Int_t AliTRDmatrix::GetTrack(Int_t iRow, Int_t iCol, Int_t iTime, Int_t iTrack)
432 // Returns the numbers of the tracks passing through one specific pixel
435 if ((iTrack < 0) || (iTrack >= kTrackPixel)) {
436 printf("AliTRDmatrix::GetTrack -- ");
437 printf("Index out of bounds (%d)\n",iTrack);
441 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
443 return (pixel->GetTrack(iTrack));
451 //_____________________________________________________________________________
452 Int_t AliTRDmatrix::GetIndex(Int_t iRow, Int_t iCol, Int_t iTime)
455 if ((iRow >= 0) && (iRow < fRow ) &&
456 (iCol >= 0) && (iCol < fCol ) &&
457 (iTime >= 0) && (iTime < fTime)) {
458 return (iTime + iCol * fTime + iRow * fTime * fCol);
466 //_____________________________________________________________________________
467 AliTRDpixel *AliTRDmatrix::GetPixel(Int_t iRow, Int_t iCol, Int_t iTime)
470 Int_t iPixel = GetIndex(iRow,iCol,iTime);
475 return ((AliTRDpixel *) fPixelArray->At(iPixel));