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