]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TRD/AliTRDmatrix.cxx
Forgot this one last time...
[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$
0a29d0f1 18Revision 1.10 2000/11/20 08:56:51 cblume
19Fix the binning of the histograms
20
580727b0 21Revision 1.9 2000/11/01 14:53:21 cblume
22Merge with TRD-develop
23
793ff80c 24Revision 1.4.2.5 2000/10/17 02:27:34 cblume
25Get rid of global constants
26
27Revision 1.4.2.4 2000/10/06 16:49:46 cblume
28Made Getters const
29
30Revision 1.4.2.3 2000/10/04 16:34:58 cblume
31Replace include files by forward declarations
32
33Revision 1.8 2000/06/09 11:10:07 cblume
34Compiler warnings and coding conventions, next round
35
dd9a6ee3 36Revision 1.7 2000/06/08 18:32:58 cblume
37Make code compliant to coding conventions
38
8230f242 39Revision 1.6 2000/05/08 15:48:30 cblume
40Resolved merge conflict
41
9583ab1c 42Revision 1.4.2.2 2000/05/08 14:50:58 cblume
43Add functions ProjRow(), ProjCol(), and ProjTime()
44
45Revision 1.4.2.1 2000/04/27 12:47:02 cblume
46Replace Fill3() by Fill()
47
1befd3b2 48Revision 1.4 2000/02/28 19:10:26 cblume
49Include the new TRD classes
50
851d3db9 51Revision 1.3.4.1 2000/02/28 17:57:47 cblume
52GetTrack returns now -1 if no track is found
53
54Revision 1.3 1999/10/04 14:48:07 fca
55Avoid warnings on non-ansi compiler HP-UX CC
56
7d83513f 57Revision 1.2 1999/09/29 09:24:35 fca
58Introduction of the Copyright and cvs Log
59
4c039060 60*/
61
99d5402e 62///////////////////////////////////////////////////////////////////////////////
63// //
64// Contains the pixel information for one TRD chamber //
65// //
0a29d0f1 66// //
99d5402e 67///////////////////////////////////////////////////////////////////////////////
68
793ff80c 69#include <TObjArray.h>
70#include <TH2.h>
71#include <TH3.h>
72#include <TStyle.h>
73#include <TCanvas.h>
74
99d5402e 75#include "AliTRDmatrix.h"
793ff80c 76#include "AliTRDpixel.h"
99d5402e 77
78ClassImp(AliTRDmatrix)
79
80//_____________________________________________________________________________
81AliTRDmatrix::AliTRDmatrix():TObject()
82{
83 //
84 // Create a TRD detector matrix
85 //
86
87 fRow = 0;
88 fCol = 0;
89 fTime = 0;
90 fPixel = 0;
91 fSector = 0;
92 fChamber = 0;
93 fPlane = 0;
94 fPixelArray = NULL;
95
96}
97
98//_____________________________________________________________________________
99AliTRDmatrix::AliTRDmatrix(Int_t nRow, Int_t nCol, Int_t nTime
100 , Int_t iSec, Int_t iCha, Int_t iPla)
101 :TObject()
102{
103 //
104 // Create a TRD detector matrix with a given size
105 //
106
107 fRow = nRow;
108 fCol = nCol;
109 fTime = nTime;
110 fPixel = nRow * nCol * nTime;
111 fSector = iSec;
112 fChamber = iCha;
113 fPlane = iPla;
114 fPixelArray = new TObjArray(fPixel);
115 for (Int_t iPixel = 0; iPixel < fPixel; iPixel++) {
116 AliTRDpixel *pixel = new AliTRDpixel();
117 fPixelArray->Add(pixel);
118 }
119
120}
121
8230f242 122//_____________________________________________________________________________
dd9a6ee3 123AliTRDmatrix::AliTRDmatrix(const AliTRDmatrix &m)
8230f242 124{
125 //
126 // AliTRDmatrix copy constructor
127 //
128
dd9a6ee3 129 ((AliTRDmatrix &) m).Copy(*this);
8230f242 130
131}
132
99d5402e 133//_____________________________________________________________________________
134AliTRDmatrix::~AliTRDmatrix()
135{
8230f242 136 //
137 // AliTRDmatrix destructor
138 //
99d5402e 139
140 if (fPixelArray) {
141 fPixelArray->Delete();
142 delete fPixelArray;
143 }
144
145}
146
dd9a6ee3 147//_____________________________________________________________________________
148AliTRDmatrix &AliTRDmatrix::operator=(const AliTRDmatrix &m)
149{
150 //
151 // Assignment operator
152 //
153
154 if (this != &m) ((AliTRDmatrix &) m).Copy(*this);
155 return *this;
156
157}
158
99d5402e 159//_____________________________________________________________________________
160void AliTRDmatrix::AddSignal(Int_t iRow, Int_t iCol, Int_t iTime, Float_t signal)
161{
162 //
163 // Add a value to the amplitude of the signal for one specific pixel
164 //
165
166 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
167 if (pixel) {
168 signal += pixel->GetSignal();
169 pixel->SetSignal(signal);
170 }
171
172}
173
8230f242 174//_____________________________________________________________________________
dd9a6ee3 175void AliTRDmatrix::Copy(TObject &m)
8230f242 176{
177 //
178 // Copy function
179 //
180
dd9a6ee3 181 ((AliTRDmatrix &) m).fRow = fRow;
182 ((AliTRDmatrix &) m).fCol = fCol;
183 ((AliTRDmatrix &) m).fTime = fTime;
184 ((AliTRDmatrix &) m).fPixel = fPixel;
185 ((AliTRDmatrix &) m).fSector = fSector;
186 ((AliTRDmatrix &) m).fChamber = fChamber;
187 ((AliTRDmatrix &) m).fPlane = fPlane;
8230f242 188
dd9a6ee3 189 ((AliTRDmatrix &) m).fPixelArray = new TObjArray(*fPixelArray);
8230f242 190
191}
192
99d5402e 193//_____________________________________________________________________________
7d83513f 194void AliTRDmatrix::Draw(Option_t *)
99d5402e 195{
196 //
197 // Draws a 3D view of the detector matrix
198 //
199
200 Char_t ctitle[50];
201 sprintf(ctitle,"Matrix (Sector:%d Chamber:%d Plane:%d)"
202 ,fSector,fChamber,fPlane);
580727b0 203 TH3F *hMatrix = new TH3F("hMatrix",ctitle,fRow ,-0.5,fRow -0.5
204 ,fCol ,-0.5,fCol -0.5
205 ,fTime,-0.5,fTime-0.5);
99d5402e 206
207 for (Int_t iRow = 0; iRow < fRow; iRow++ ) {
208 for (Int_t iCol = 0; iCol < fCol; iCol++ ) {
209 for (Int_t iTime = 0; iTime < fTime; iTime++) {
210 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
1befd3b2 211 if (pixel) hMatrix->Fill(iRow,iCol,iTime,pixel->GetSignal());
99d5402e 212 }
213 }
214 }
215
216 gStyle->SetOptStat(0);
217 TCanvas *cMatrix = new TCanvas("cMatrix","Detector matrix 3D-view"
218 ,50,50,600,400);
219 cMatrix->ToggleEventStatus();
220 hMatrix->SetXTitle("Pad-row (z)");
221 hMatrix->SetYTitle("Pad-column (rphi)");
222 hMatrix->SetZTitle("Timebucket");
223 hMatrix->Draw("BOX");
224
225}
226
227//_____________________________________________________________________________
228void AliTRDmatrix::DrawRow(Int_t iRow)
229{
230 //
231 // Draws a 2D slice of the detector matrix along one row
232 //
233
234 if ((iRow < 0) || (iRow >= fRow)) {
851d3db9 235 printf("AliTRDmatrix::DrawRow -- ");
99d5402e 236 printf("Index out of bounds (%d/%d)\n",iRow,fRow);
237 return;
238 }
239
240 Char_t ctitle[50];
241 sprintf(ctitle,"Pad-row %d (Sector:%d Chamber:%d Plane:%d)"
242 ,iRow,fSector,fChamber,fPlane);
580727b0 243 TH2F *hSliceRow = new TH2F("hSliceRow",ctitle,fCol ,-0.5,fCol -0.5
244 ,fTime,-0.5,fTime-0.5);
99d5402e 245
246 for (Int_t iCol = 0; iCol < fCol; iCol++ ) {
247 for (Int_t iTime = 0; iTime < fTime; iTime++) {
248 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
249 if (pixel) hSliceRow->Fill(iCol,iTime,pixel->GetSignal());
250 }
251 }
252
253 gStyle->SetOptStat(0);
254 TCanvas *cSliceRow = new TCanvas("cSliceRow","Detector matrix 2D-slice"
255 ,50,50,600,400);
256 cSliceRow->ToggleEventStatus();
257 hSliceRow->SetXTitle("Pad-column (rphi)");
258 hSliceRow->SetYTitle("Timebucket");
259 hSliceRow->Draw("COLZ");
260
261}
262
263//_____________________________________________________________________________
264void AliTRDmatrix::DrawCol(Int_t iCol)
265{
266 //
267 // Draws a 2D slice of the detector matrix along one column
268 //
269
270 if ((iCol < 0) || (iCol >= fCol)) {
851d3db9 271 printf("AliTRDmatrix::DrawCol -- ");
99d5402e 272 printf("Index out of bounds (%d/%d)\n",iCol,fCol);
273 return;
274 }
275
276 Char_t ctitle[50];
277 sprintf(ctitle,"Pad-column %d (Sector:%d Chamber:%d Plane:%d)"
278 ,iCol,fSector,fChamber,fPlane);
580727b0 279 TH2F *hSliceCol = new TH2F("hSliceCol",ctitle,fRow ,-0.5,fRow -0.5
280 ,fTime,-0.5,fTime-0.5);
99d5402e 281
282 for (Int_t iRow = 0; iRow < fRow; iRow++ ) {
283 for (Int_t iTime = 0; iTime < fTime; iTime++) {
284 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
285 if (pixel) hSliceCol->Fill(iRow,iTime,pixel->GetSignal());
286 }
287 }
288
289 gStyle->SetOptStat(0);
290 TCanvas *cSliceCol = new TCanvas("cSliceCol","Detector matrix 2D-slice"
291 ,50,50,600,400);
292 cSliceCol->ToggleEventStatus();
293 hSliceCol->SetXTitle("Pad-row (z)");
294 hSliceCol->SetYTitle("Timebucket");
295 hSliceCol->Draw("COLZ");
296
297}
298
299//_____________________________________________________________________________
300void AliTRDmatrix::DrawTime(Int_t iTime)
301{
302 //
303 // Draws a 2D slice of the detector matrix along one time slice
304 //
305
306 if ((iTime < 0) || (iTime >= fTime)) {
851d3db9 307 printf("AliTRDmatrix::DrawTime -- ");
99d5402e 308 printf("Index out of bounds (%d/%d)\n",iTime,fTime);
309 return;
310 }
311
312 Char_t ctitle[50];
313 sprintf(ctitle,"Time-slice %d (Sector:%d Chamber:%d Plane:%d)"
314 ,iTime,fSector,fChamber,fPlane);
580727b0 315 TH2F *hSliceTime = new TH2F("hSliceTime",ctitle,fRow,-0.5,fRow-0.5
316 ,fCol,-0.5,fCol-0.5);
99d5402e 317
318 for (Int_t iRow = 0; iRow < fRow; iRow++) {
319 for (Int_t iCol = 0; iCol < fCol; iCol++) {
320 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
321 if (pixel) hSliceTime->Fill(iRow,iCol,pixel->GetSignal());
322 }
323 }
324
325 gStyle->SetOptStat(0);
326 TCanvas *cSliceTime = new TCanvas("cSliceTime","Detector matrix 2D-slice"
327 ,50,50,600,400);
328 cSliceTime->ToggleEventStatus();
329 hSliceTime->SetXTitle("Pad-row (z)");
330 hSliceTime->SetYTitle("Pad-column (rphi)");
331 hSliceTime->Draw("COLZ");
332
333}
334
9583ab1c 335//_____________________________________________________________________________
336void AliTRDmatrix::ProjRow()
337{
338 //
339 // Projects the detector matrix along the row-axis
340 //
341
342 Char_t ctitle[60];
343 sprintf(ctitle,"Row-projection (Sector:%d Chamber:%d Plane:%d)"
344 ,fSector,fChamber,fPlane);
580727b0 345 TH2F *hProjRow = new TH2F("hProjRow",ctitle,fCol ,-0.5,fCol -0.5
346 ,fTime,-0.5,fTime-0.5);
9583ab1c 347
348 for (Int_t iRow = 0; iRow < fRow; iRow++ ) {
349 for (Int_t iCol = 0; iCol < fCol; iCol++ ) {
350 for (Int_t iTime = 0; iTime < fTime; iTime++) {
351 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
352 if (pixel) hProjRow->Fill(iCol,iTime,pixel->GetSignal());
353 }
354 }
355 }
356
357 gStyle->SetOptStat(0);
358 TCanvas *cProjRow = new TCanvas("cProjRow","Detector matrix 2D-projection"
359 ,50,50,600,400);
360 cProjRow->ToggleEventStatus();
361 hProjRow->SetXTitle("Pad-column (rphi)");
362 hProjRow->SetYTitle("Timebucket");
363 hProjRow->Draw("COLZ");
364
365}
366
367//_____________________________________________________________________________
368void AliTRDmatrix::ProjCol()
369{
370 //
371 // Projects the detector matrix along the column-axis
372 //
373
374 Char_t ctitle[60];
375 sprintf(ctitle,"Column-projection (Sector:%d Chamber:%d Plane:%d)"
376 ,fSector,fChamber,fPlane);
580727b0 377 TH2F *hProjCol = new TH2F("hProjCol",ctitle,fRow ,-0.5,fRow -0.5
378 ,fTime,-0.5,fTime-0.5);
9583ab1c 379
380 for (Int_t iRow = 0; iRow < fRow; iRow++ ) {
381 for (Int_t iCol = 0; iCol < fCol; iCol++ ) {
382 for (Int_t iTime = 0; iTime < fTime; iTime++) {
383 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
384 if (pixel) hProjCol->Fill(iRow,iTime,pixel->GetSignal());
385 }
386 }
387 }
388
389 gStyle->SetOptStat(0);
390 TCanvas *cProjCol = new TCanvas("cProjCol","Detector matrix 2D-projection"
391 ,50,50,600,400);
392 cProjCol->ToggleEventStatus();
393 hProjCol->SetXTitle("Pad-row (z)");
394 hProjCol->SetYTitle("Timebucket");
395 hProjCol->Draw("COLZ");
396
397}
398
399//_____________________________________________________________________________
400void AliTRDmatrix::ProjTime()
401{
402 //
403 // Projects the detector matrix along the time-axis
404 //
405
406 Char_t ctitle[60];
407 sprintf(ctitle,"Time-projection (Sector:%d Chamber:%d Plane:%d)"
408 ,fSector,fChamber,fPlane);
580727b0 409 TH2F *hProjTime = new TH2F("hProjTime",ctitle,fRow,-0.5,fRow-0.5
410 ,fCol,-0.5,fCol-0.5);
9583ab1c 411
412 for (Int_t iRow = 0; iRow < fRow; iRow++) {
413 for (Int_t iCol = 0; iCol < fCol; iCol++) {
414 for (Int_t iTime = 0; iTime < fTime; iTime++) {
415 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
416 if (pixel) hProjTime->Fill(iRow,iCol,pixel->GetSignal());
417 }
418 }
419 }
420
421 gStyle->SetOptStat(0);
422 TCanvas *cProjTime = new TCanvas("cProjTime","Detector matrix 2D-projection"
423 ,50,50,600,400);
424 cProjTime->ToggleEventStatus();
425 hProjTime->SetXTitle("Pad-row (z)");
426 hProjTime->SetYTitle("Pad-column (rphi)");
427 hProjTime->Draw("COLZ");
428
429}
430
99d5402e 431//_____________________________________________________________________________
432void AliTRDmatrix::SetSignal(Int_t iRow, Int_t iCol, Int_t iTime, Float_t signal)
433{
434 //
435 // Set the amplitude of the signal for one specific pixel
436 //
437
438 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
439 if (pixel) {
440 pixel->SetSignal(signal);
441 }
442
443}
444
445//_____________________________________________________________________________
446Bool_t AliTRDmatrix::AddTrack(Int_t iRow, Int_t iCol, Int_t iTime, Int_t track)
447{
448 //
449 // Add this track number to the stored tracks passing through this pixel.
450 // If there are already three stored the return status is FALSE.
451 //
452
453 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
454 if (!(pixel)) return kTRUE;
455
456 Bool_t trackSet = kFALSE;
793ff80c 457 for (Int_t i = 0; i < AliTRDpixel::NTrackPixel(); i++) {
99d5402e 458 if (pixel->GetTrack(i) == track) {
459 trackSet = kTRUE;
460 break;
461 }
851d3db9 462 if (pixel->GetTrack(i) == -1) {
99d5402e 463 pixel->SetTrack(i,track);
464 trackSet = kTRUE;
465 break;
466 }
467 }
468
469 return trackSet;
470
471}
472
473//_____________________________________________________________________________
474void AliTRDmatrix::SetTrack(Int_t iRow, Int_t iCol, Int_t iTime
475 , Int_t iTrack, Int_t track)
476{
477 //
478 // Store the number of a track which is passing through this pixel
479 //
480
481 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
482 if (pixel) {
483 pixel->SetTrack(iTrack,track);
484 }
485
486}
487
488//_____________________________________________________________________________
793ff80c 489Float_t AliTRDmatrix::GetSignal(Int_t iRow, Int_t iCol, Int_t iTime) const
99d5402e 490{
491 //
492 // Returns the amplitude of the signal for one specific pixel
493 //
494
495 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
496 if (pixel) {
497 return (pixel->GetSignal());
498 }
499 else {
500 return 0;
501 }
502
503}
504
505//_____________________________________________________________________________
793ff80c 506Int_t AliTRDmatrix::GetTrack(Int_t iRow, Int_t iCol, Int_t iTime
507 , Int_t iTrack) const
99d5402e 508{
509 //
510 // Returns the numbers of the tracks passing through one specific pixel
511 //
512
793ff80c 513 if ((iTrack < 0) || (iTrack >= AliTRDpixel::NTrackPixel())) {
851d3db9 514 printf("AliTRDmatrix::GetTrack -- ");
515 printf("Index out of bounds (%d)\n",iTrack);
516 return -1;
99d5402e 517 }
518
519 AliTRDpixel *pixel = GetPixel(iRow,iCol,iTime);
520 if (pixel) {
521 return (pixel->GetTrack(iTrack));
522 }
523 else {
851d3db9 524 return -1;
99d5402e 525 }
526
527}
528
529//_____________________________________________________________________________
793ff80c 530Int_t AliTRDmatrix::GetIndex(Int_t iRow, Int_t iCol, Int_t iTime) const
99d5402e 531{
0a29d0f1 532 //
533 // Get the index of a given pixel
534 //
99d5402e 535
536 if ((iRow >= 0) && (iRow < fRow ) &&
537 (iCol >= 0) && (iCol < fCol ) &&
538 (iTime >= 0) && (iTime < fTime)) {
539 return (iTime + iCol * fTime + iRow * fTime * fCol);
540 }
541 else {
542 return -1;
543 }
544
545}
546
547//_____________________________________________________________________________
793ff80c 548AliTRDpixel *AliTRDmatrix::GetPixel(Int_t iRow, Int_t iCol, Int_t iTime) const
99d5402e 549{
0a29d0f1 550 //
551 // Get one pixel
552 //
99d5402e 553
554 Int_t iPixel = GetIndex(iRow,iCol,iTime);
555 if (iPixel < 0) {
556 return NULL;
557 }
558 else {
559 return ((AliTRDpixel *) fPixelArray->At(iPixel));
560 }
561
562}