]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TRD/AliTRDmatrix.cxx
updated from TPC-dev1
[u/mrichter/AliRoot.git] / TRD / AliTRDmatrix.cxx
CommitLineData
4c039060 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
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 **************************************************************************/
15
16/*
17$Log$
851d3db9 18Revision 1.3.4.1 2000/02/28 17:57:47 cblume
19GetTrack returns now -1 if no track is found
20
21Revision 1.3 1999/10/04 14:48:07 fca
22Avoid warnings on non-ansi compiler HP-UX CC
23
7d83513f 24Revision 1.2 1999/09/29 09:24:35 fca
25Introduction of the Copyright and cvs Log
26
4c039060 27*/
28
99d5402e 29///////////////////////////////////////////////////////////////////////////////
30// //
31// Contains the pixel information for one TRD chamber //
32// //
33///////////////////////////////////////////////////////////////////////////////
34
35#include "AliTRDmatrix.h"
36
37ClassImp(AliTRDmatrix)
38
39//_____________________________________________________________________________
40AliTRDmatrix::AliTRDmatrix():TObject()
41{
42 //
43 // Create a TRD detector matrix
44 //
45
46 fRow = 0;
47 fCol = 0;
48 fTime = 0;
49 fPixel = 0;
50 fSector = 0;
51 fChamber = 0;
52 fPlane = 0;
53 fPixelArray = NULL;
54
55}
56
57//_____________________________________________________________________________
58AliTRDmatrix::AliTRDmatrix(Int_t nRow, Int_t nCol, Int_t nTime
59 , Int_t iSec, Int_t iCha, Int_t iPla)
60 :TObject()
61{
62 //
63 // Create a TRD detector matrix with a given size
64 //
65
66 fRow = nRow;
67 fCol = nCol;
68 fTime = nTime;
69 fPixel = nRow * nCol * nTime;
70 fSector = iSec;
71 fChamber = iCha;
72 fPlane = iPla;
73 fPixelArray = new TObjArray(fPixel);
74 for (Int_t iPixel = 0; iPixel < fPixel; iPixel++) {
75 AliTRDpixel *pixel = new AliTRDpixel();
76 fPixelArray->Add(pixel);
77 }
78
79}
80
81//_____________________________________________________________________________
82AliTRDmatrix::~AliTRDmatrix()
83{
84
85 if (fPixelArray) {
86 fPixelArray->Delete();
87 delete fPixelArray;
88 }
89
90}
91
92//_____________________________________________________________________________
93void AliTRDmatrix::AddSignal(Int_t iRow, Int_t iCol, Int_t iTime, Float_t signal)
94{
95 //
96 // Add a value to the amplitude of the signal for one specific pixel
97 //
98
99 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
100 if (pixel) {
101 signal += pixel->GetSignal();
102 pixel->SetSignal(signal);
103 }
104
105}
106
107//_____________________________________________________________________________
7d83513f 108void AliTRDmatrix::Draw(Option_t *)
99d5402e 109{
110 //
111 // Draws a 3D view of the detector matrix
112 //
113
114 Char_t ctitle[50];
115 sprintf(ctitle,"Matrix (Sector:%d Chamber:%d Plane:%d)"
116 ,fSector,fChamber,fPlane);
117 TH3F *hMatrix = new TH3F("hMatrix",ctitle,fRow ,-0.5,fRow +0.5
118 ,fCol ,-0.5,fCol +0.5
119 ,fTime,-0.5,fTime+0.5);
120
121 for (Int_t iRow = 0; iRow < fRow; iRow++ ) {
122 for (Int_t iCol = 0; iCol < fCol; iCol++ ) {
123 for (Int_t iTime = 0; iTime < fTime; iTime++) {
124 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
125 if (pixel) hMatrix->Fill3(iRow,iCol,iTime,pixel->GetSignal());
126 }
127 }
128 }
129
130 gStyle->SetOptStat(0);
131 TCanvas *cMatrix = new TCanvas("cMatrix","Detector matrix 3D-view"
132 ,50,50,600,400);
133 cMatrix->ToggleEventStatus();
134 hMatrix->SetXTitle("Pad-row (z)");
135 hMatrix->SetYTitle("Pad-column (rphi)");
136 hMatrix->SetZTitle("Timebucket");
137 hMatrix->Draw("BOX");
138
139}
140
141//_____________________________________________________________________________
142void AliTRDmatrix::DrawRow(Int_t iRow)
143{
144 //
145 // Draws a 2D slice of the detector matrix along one row
146 //
147
148 if ((iRow < 0) || (iRow >= fRow)) {
851d3db9 149 printf("AliTRDmatrix::DrawRow -- ");
99d5402e 150 printf("Index out of bounds (%d/%d)\n",iRow,fRow);
151 return;
152 }
153
154 Char_t ctitle[50];
155 sprintf(ctitle,"Pad-row %d (Sector:%d Chamber:%d Plane:%d)"
156 ,iRow,fSector,fChamber,fPlane);
157 TH2F *hSliceRow = new TH2F("hSliceRow",ctitle,fCol ,-0.5,fCol +0.5
158 ,fTime,-0.5,fTime+0.5);
159
160 for (Int_t iCol = 0; iCol < fCol; iCol++ ) {
161 for (Int_t iTime = 0; iTime < fTime; iTime++) {
162 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
163 if (pixel) hSliceRow->Fill(iCol,iTime,pixel->GetSignal());
164 }
165 }
166
167 gStyle->SetOptStat(0);
168 TCanvas *cSliceRow = new TCanvas("cSliceRow","Detector matrix 2D-slice"
169 ,50,50,600,400);
170 cSliceRow->ToggleEventStatus();
171 hSliceRow->SetXTitle("Pad-column (rphi)");
172 hSliceRow->SetYTitle("Timebucket");
173 hSliceRow->Draw("COLZ");
174
175}
176
177//_____________________________________________________________________________
178void AliTRDmatrix::DrawCol(Int_t iCol)
179{
180 //
181 // Draws a 2D slice of the detector matrix along one column
182 //
183
184 if ((iCol < 0) || (iCol >= fCol)) {
851d3db9 185 printf("AliTRDmatrix::DrawCol -- ");
99d5402e 186 printf("Index out of bounds (%d/%d)\n",iCol,fCol);
187 return;
188 }
189
190 Char_t ctitle[50];
191 sprintf(ctitle,"Pad-column %d (Sector:%d Chamber:%d Plane:%d)"
192 ,iCol,fSector,fChamber,fPlane);
193 TH2F *hSliceCol = new TH2F("hSliceCol",ctitle,fRow ,-0.5,fRow +0.5
194 ,fTime,-0.5,fTime+0.5);
195
196 for (Int_t iRow = 0; iRow < fRow; iRow++ ) {
197 for (Int_t iTime = 0; iTime < fTime; iTime++) {
198 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
199 if (pixel) hSliceCol->Fill(iRow,iTime,pixel->GetSignal());
200 }
201 }
202
203 gStyle->SetOptStat(0);
204 TCanvas *cSliceCol = new TCanvas("cSliceCol","Detector matrix 2D-slice"
205 ,50,50,600,400);
206 cSliceCol->ToggleEventStatus();
207 hSliceCol->SetXTitle("Pad-row (z)");
208 hSliceCol->SetYTitle("Timebucket");
209 hSliceCol->Draw("COLZ");
210
211}
212
213//_____________________________________________________________________________
214void AliTRDmatrix::DrawTime(Int_t iTime)
215{
216 //
217 // Draws a 2D slice of the detector matrix along one time slice
218 //
219
220 if ((iTime < 0) || (iTime >= fTime)) {
851d3db9 221 printf("AliTRDmatrix::DrawTime -- ");
99d5402e 222 printf("Index out of bounds (%d/%d)\n",iTime,fTime);
223 return;
224 }
225
226 Char_t ctitle[50];
227 sprintf(ctitle,"Time-slice %d (Sector:%d Chamber:%d Plane:%d)"
228 ,iTime,fSector,fChamber,fPlane);
229 TH2F *hSliceTime = new TH2F("hSliceTime",ctitle,fRow,-0.5,fRow+0.5
230 ,fCol,-0.5,fCol+0.5);
231
232 for (Int_t iRow = 0; iRow < fRow; iRow++) {
233 for (Int_t iCol = 0; iCol < fCol; iCol++) {
234 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
235 if (pixel) hSliceTime->Fill(iRow,iCol,pixel->GetSignal());
236 }
237 }
238
239 gStyle->SetOptStat(0);
240 TCanvas *cSliceTime = new TCanvas("cSliceTime","Detector matrix 2D-slice"
241 ,50,50,600,400);
242 cSliceTime->ToggleEventStatus();
243 hSliceTime->SetXTitle("Pad-row (z)");
244 hSliceTime->SetYTitle("Pad-column (rphi)");
245 hSliceTime->Draw("COLZ");
246
247}
248
249//_____________________________________________________________________________
250void AliTRDmatrix::SetSignal(Int_t iRow, Int_t iCol, Int_t iTime, Float_t signal)
251{
252 //
253 // Set the amplitude of the signal for one specific pixel
254 //
255
256 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
257 if (pixel) {
258 pixel->SetSignal(signal);
259 }
260
261}
262
263//_____________________________________________________________________________
264Bool_t AliTRDmatrix::AddTrack(Int_t iRow, Int_t iCol, Int_t iTime, Int_t track)
265{
266 //
267 // Add this track number to the stored tracks passing through this pixel.
268 // If there are already three stored the return status is FALSE.
269 //
270
271 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
272 if (!(pixel)) return kTRUE;
273
274 Bool_t trackSet = kFALSE;
275 for (Int_t i = 0; i < kTrackPixel; i++) {
276 if (pixel->GetTrack(i) == track) {
277 trackSet = kTRUE;
278 break;
279 }
851d3db9 280 if (pixel->GetTrack(i) == -1) {
99d5402e 281 pixel->SetTrack(i,track);
282 trackSet = kTRUE;
283 break;
284 }
285 }
286
287 return trackSet;
288
289}
290
291//_____________________________________________________________________________
292void AliTRDmatrix::SetTrack(Int_t iRow, Int_t iCol, Int_t iTime
293 , Int_t iTrack, Int_t track)
294{
295 //
296 // Store the number of a track which is passing through this pixel
297 //
298
299 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
300 if (pixel) {
301 pixel->SetTrack(iTrack,track);
302 }
303
304}
305
306//_____________________________________________________________________________
307Float_t AliTRDmatrix::GetSignal(Int_t iRow, Int_t iCol, Int_t iTime)
308{
309 //
310 // Returns the amplitude of the signal for one specific pixel
311 //
312
313 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
314 if (pixel) {
315 return (pixel->GetSignal());
316 }
317 else {
318 return 0;
319 }
320
321}
322
323//_____________________________________________________________________________
324Int_t AliTRDmatrix::GetTrack(Int_t iRow, Int_t iCol, Int_t iTime, Int_t iTrack)
325{
326 //
327 // Returns the numbers of the tracks passing through one specific pixel
328 //
329
330 if ((iTrack < 0) || (iTrack >= kTrackPixel)) {
851d3db9 331 printf("AliTRDmatrix::GetTrack -- ");
332 printf("Index out of bounds (%d)\n",iTrack);
333 return -1;
99d5402e 334 }
335
336 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
337 if (pixel) {
338 return (pixel->GetTrack(iTrack));
339 }
340 else {
851d3db9 341 return -1;
99d5402e 342 }
343
344}
345
346//_____________________________________________________________________________
347Int_t AliTRDmatrix::GetIndex(Int_t iRow, Int_t iCol, Int_t iTime)
348{
349
350 if ((iRow >= 0) && (iRow < fRow ) &&
351 (iCol >= 0) && (iCol < fCol ) &&
352 (iTime >= 0) && (iTime < fTime)) {
353 return (iTime + iCol * fTime + iRow * fTime * fCol);
354 }
355 else {
356 return -1;
357 }
358
359}
360
361//_____________________________________________________________________________
362AliTRDpixel *AliTRDmatrix::GetPixel(Int_t iRow, Int_t iCol, Int_t iTime)
363{
364
365 Int_t iPixel = GetIndex(iRow,iCol,iTime);
366 if (iPixel < 0) {
367 return NULL;
368 }
369 else {
370 return ((AliTRDpixel *) fPixelArray->At(iPixel));
371 }
372
373}