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