]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDanaDigits.C
Pyhtia comparison extended
[u/mrichter/AliRoot.git] / TRD / AliTRDanaDigits.C
1 void AliTRDanaDigits() 
2 {
3
4 /////////////////////////////////////////////////////////////////////////
5 //
6 // Example macro for the analysis of the TRD digits and the use
7 // of the AliTRDmatrix class.
8 //
9 /////////////////////////////////////////////////////////////////////////
10
11   // Dynamically link some shared libs
12   if (gClassTable->GetID("AliRun") < 0) {
13     gROOT->LoadMacro("loadlibs.C");
14     loadlibs();
15   }
16
17   // Input file name
18   Char_t *alifile = "galice.root"; 
19
20   // Event number
21   Int_t   nEvent  = 0;
22
23   // Define the objects
24   AliTRDv1        *trd;
25   AliTRDgeometry  *geo;
26   AliTRDdigit     *digit;
27   AliTRDparameter *par;
28
29   Int_t           track;
30
31   AliRunLoader* rl = AliRunLoader::Open(alifile);
32   AliLoader* loader = rl->GetLoader("TRDLoader");
33   rl->LoadDigits();
34   
35   rl->LoadgAlice();
36   gAlice = rl->GetAliRun();
37   
38   // Get the pointer to the detector object
39   trd = (AliTRDv1*) gAlice->GetDetector("TRD");
40
41   // Get the pointer to the geometry object
42   if (trd) {
43     geo = trd->GetGeometry();
44   }
45   else {
46     cout << "Cannot find the geometry" << endl;
47     break;
48   }
49
50   // The parameter object
51   par = new AliTRDparameter("TRDparameter","TRD parameter class");
52
53   // Create the digits manager
54   AliTRDdigitsManager *digitsManager = new AliTRDdigitsManager();
55   digitsManager->SetDebug(1);
56
57   // Read the digits from the file
58   
59   digitsManager->ReadDigits(loader->TreeD());
60
61   // Get the detector number
62   Int_t iDet = 514;
63   cout << " iDet = " << iDet << endl;
64
65   // Define the detector matrix for one chamber
66   const Int_t iSec = geo->GetSector(iDet);
67   const Int_t iCha = geo->GetChamber(iDet);
68   const Int_t iPla = geo->GetPlane(iDet);
69   Int_t  rowMax = par->GetRowMax(iPla,iCha,iSec);
70   Int_t  colMax = par->GetColMax(iPla);
71   Int_t timeMax = par->GetTimeMax();
72   cout << "Geometry: rowMax = "  <<  rowMax
73                 << " colMax = "  <<  colMax
74                 << " timeMax = " << timeMax << endl;
75   AliTRDmatrix *matrix = new AliTRDmatrix(rowMax,colMax,timeMax,iSec,iCha,iPla);
76
77   // Loop through the detector pixel
78   for (Int_t time = 0; time < timeMax; time++) {
79     for (Int_t  col = 0;  col <  colMax;  col++) {
80       for (Int_t  row = 0;  row <  rowMax;  row++) {
81
82         digit = digitsManager->GetDigit(row,col,time,iDet);
83         track = digitsManager->GetTrack(0,row,col,time,iDet);
84         
85         matrix->SetSignal(row,col,time,digit->GetAmp());
86
87         delete digit;
88
89       }
90     }
91   }
92
93   // Display the detector matrix
94   matrix->Draw();
95   matrix->ProjRow();
96   matrix->ProjCol();
97   matrix->ProjTime();
98   delete rl;
99 }