]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/digitsAna.C
Avoid non ansi warnings on HP compilers
[u/mrichter/AliRoot.git] / TRD / digitsAna.C
1 {
2
3 /////////////////////////////////////////////////////////////////////////
4 //
5 // Example macro for the analysis of the TRD digits and the use
6 // of the AliTRDmatrix class.
7 //
8 /////////////////////////////////////////////////////////////////////////
9
10   // Dynamically link some shared libs
11   if (gClassTable->GetID("AliRun") < 0) {
12     gROOT->LoadMacro("loadlibs.C");
13     loadlibs();
14   }
15
16   // Input file name
17   Char_t *alifile = "galice_v2.root"; 
18
19   // Event number
20   Int_t   nEvent  = 0;
21
22   // Define the objects
23   AliTRDv2     *TRD;
24   TClonesArray *TRDDigits;
25   AliTRDdigit  *OneTRDDigit;
26
27   // Connect the AliRoot file containing Geometry, Kine, Hits, and Digits
28   TFile *gafl = (TFile*) gROOT->GetListOfFiles()->FindObject(alifile);
29   if (!gafl) {
30     cout << "Open the ALIROOT-file " << alifile << endl;
31     gafl = new TFile(alifile);
32   }
33   else {
34     cout << alifile << " is already open" << endl;
35   }
36
37   // Get AliRun object from file or create it if not on file
38   if (!gAlice) {
39     gAlice = (AliRun*) gafl->Get("gAlice");
40     if (gAlice)  
41       cout << "AliRun object found on file" << endl;
42     else
43       gAlice = new AliRun("gAlice","Alice test program");
44   }
45
46   // Import the Trees for the event nEvent in the file
47   Int_t nparticles = gAlice->GetEvent(nEvent);
48   if (nparticles <= 0) break;
49   
50   // Get the pointer to the hit-tree
51   TTree *DigitsTree = gAlice->TreeD();
52
53   // Get the pointer to the detector classes
54   TRD = (AliTRDv2*) gAlice->GetDetector("TRD");
55   // Get the pointer to the hit container
56   if (TRD) TRDDigits = TRD->Digits();
57
58   // Define the detector matrix for one chamber (Sector 6, Chamber 3, Plane 1)
59   const Int_t iSec = 6;
60   const Int_t iCha = 3;
61   const Int_t iPla = 1;
62   Int_t  rowMax = TRD->GetRowMax(iPla,iCha);
63   Int_t  colMax = TRD->GetColMax(iPla);
64   Int_t timeMax = TRD->GetTimeMax();
65   cout << " rowMax = "  << rowMax
66        << " colMax = "  << colMax
67        << " timeMax = " << timeMax << endl;
68   AliTRDmatrix *TRDMatrix = new AliTRDmatrix(rowMax,colMax,timeMax,iSec,iCha,iPla);
69
70   Int_t nEntries = DigitsTree->GetEntries();
71   cout << "Number of entries in digits tree = " << nEntries << endl; 
72
73   // Loop through all entries in the tree
74   Int_t nbytes;
75   for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
76
77     cout << "iEntry = " << iEntry << endl;
78
79     // Import the tree
80     gAlice->ResetDigits();
81     nbytes += DigitsTree->GetEvent(iEntry);
82
83     // Get the number of digits in the detector 
84     Int_t nTRDDigits = TRDDigits->GetEntriesFast();
85     cout << " nTRDDigits = " << nTRDDigits << endl;    
86
87     // Loop through all TRD digits
88     for (Int_t iTRDDigits = 0; iTRDDigits < nTRDDigits; iTRDDigits++) {
89
90       // Get the information for this digit
91       OneTRDDigit = (AliTRDdigit*) TRDDigits->UncheckedAt(iTRDDigits);
92       Int_t signal    = OneTRDDigit->fSignal;
93       Int_t   sector  = OneTRDDigit->fSector;
94       Int_t   chamber = OneTRDDigit->fChamber;
95       Int_t   plane   = OneTRDDigit->fPlane;
96       Int_t   row     = OneTRDDigit->fRow;
97       Int_t   col     = OneTRDDigit->fCol;
98       Int_t   time    = OneTRDDigit->fTime;
99
100       // Fill the detector matrix
101       if (signal > 1) {
102         TRDMatrix->SetSignal(row,col,time,signal);
103       }
104
105     }
106
107   }
108
109   // Display the detector matrix
110   TRDMatrix->Draw();
111   TRDMatrix->DrawRow(18);
112   TRDMatrix->DrawCol(58);
113   TRDMatrix->DrawTime(20);
114
115 }