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