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