]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ACORDE/macros/readACORDERawData.C
The volumes that were still coded in the Geant3 framework have been redefined within...
[u/mrichter/AliRoot.git] / ACORDE / macros / readACORDERawData.C
1 //////////////////////////////////////////////////////////////////
2 //                                                              //
3 //      readACORDERawData.C macro                               //
4 //                                                              //
5 //      Reads the information of ACORDE from raw data file      //
6 //      from the 4-32-words (SL0 & MCN mode)                    //
7 //      and draws four histograms (2 by each mode)              //
8 //                                                              //
9 //      Author: Mario Rodriguez Cahuantzi                       //
10 //      E-MaiL: mario.rocah@gmail.com, mrodrigu@mail.cern.ch    //
11 //                                                              //
12 //      Created: July 1st. 2010 @ FCFM -BUAP, Puebla, MX        //
13 //      Last update: created from old AcoReco.C macro           //
14 //                                                              //
15 //////////////////////////////////////////////////////////////////
16
17 void readACORDERawData(char* fileName)
18 {
19
20         TStopwatch timer;
21         timer.Start();
22         // Pointer to rawReader class
23
24         AliRawReader* rawReader = new AliRawReaderRoot(fileName); 
25         
26         // Pointer to rawData of Acorde
27
28         AliACORDERawStream* rawStream  = new AliACORDERawStream(rawReader);    
29         
30         // Create some histograms
31
32         TH1F *h1 = new TH1F("h1","ACORDE - Single Muon Hits (SL0)",60,-0.5,59.5);
33         TH1F *h2 = new TH1F("h2","ACORDE - Multiplicity of Acorde Modules (SL0)",61,-1,60);
34         TH1F *h3 = new TH1F("h3","ACORDE - Single Muon Hits (MCN)",60,-0.5,59.5);
35         TH1F *h4 = new TH1F("h4","ACORDE - Multiplicity of Acorde Modules (MCN)",61,-1,60);
36         
37         // Declare some counters
38         
39         size_t contSingle=0;
40         size_t contMulti=0;
41         UInt_t acorde_word[4]; // array to store the 4-words
42         bool word_sl0[60],word_mcn[60]; // boolean array if some hit in module
43         UInt_t shiftword; // shift word
44
45         for(Int_t m=0;m<60;m++) {word_sl0[m]=0;word_mcn[m]=0;}
46   
47  
48         Int_t nEvents = rawStream->GetNEvents(fileName);
49   
50         printf("File: %s, Number of events: %d \n",fileName,nEvents);
51
52
53         // Loop over all the events
54
55         for (Int_t i=1; i<=nEvents; i++) 
56         {
57
58                 if (!rawReader->NextEvent()) break;
59                 rawStream->Reset();
60                 if (!rawStream->Next()) continue;
61
62                 acorde_word[0] = rawStream->GetWord(0);
63                 acorde_word[1] = rawStream->GetWord(1);
64                 acorde_word[2] = rawStream->GetWord(2);
65                 acorde_word[3] = rawStream->GetWord(3);
66
67                 shiftword = acorde_word[0];
68                 for(Int_t iaco=0;iaco<30;iaco++)
69                 {
70                         word_sl0[iaco] = shiftword & 1;
71                         shiftword>>=1;
72                 }
73         
74                 shiftword = acorde_word[1];
75                 for(Int_t iaco=30;iaco<60;iaco++)
76                 {
77                         word_sl0[iaco] = shiftword & 1;
78                         shiftword>>=1;
79                 }
80
81                 shiftword=acorde_word[2];
82                 for(Int_t iaco=0;iaco<30;iaco++)
83                 {
84                         word_mcn[iaco] = shiftword & 1;
85                         shiftword>>=1;
86                 }
87
88                 shiftword=acorde_word[3];
89                 for(Int_t iaco=30;iaco<60;iaco++)
90                 {
91                         word_mcn[iaco] = shiftword & 1;
92                         shiftword>>=1;
93                 }
94
95                 contSingle=0;   
96                 for(Int_t iaco=0;iaco<60;iaco++) 
97                 {
98                         if(word_sl0[iaco]==1) 
99                         {
100                                 h1->Fill(iaco);
101                                 contSingle++;
102                         }
103                 
104                 }h2->Fill(contSingle);
105                 contMulti=0;
106                 for(Int_t iaco=0;iaco<60;iaco++) 
107                 {
108                         if(word_mcn[iaco]==1) 
109                         {
110                                 h3->Fill(iaco);
111                                 contMulti++;
112                         }
113                 
114                 }h4->Fill(contMulti);
115
116         }
117
118         TCanvas *acorde = new TCanvas("ACORDE","ACORDE-Histograms from Raw-Data",1);
119         acorde->Divide(2,2);
120         acorde->cd(1);
121         h1->GetXaxis()->SetTitle("No. of module");
122         h1->GetYaxis()->SetTitle("No. of Hits");
123         h1->SetFillColor(kRed);
124         h1->Draw();
125
126         acorde->cd(2);
127         gPad->SetLogy();
128         h2->GetXaxis()->SetTitle("No. of fired modules");
129         h2->GetYaxis()->SetTitle("No. of events");
130         h2->SetFillColor(kBlue);
131         h2->Draw();
132
133         acorde->cd(3);
134         h3->GetXaxis()->SetTitle("No. of module");
135         h3->GetYaxis()->SetTitle("No. of events");
136         h3->SetFillColor(kRed);
137         h3->Draw();
138
139         acorde->cd(4);
140         gPad->SetLogy();
141         h4->GetXaxis()->SetTitle("No. of fired modules");
142         h4->GetYaxis()->SetTitle("No. of events");
143         h4->SetFillColor(kBlue);
144         h4->Draw();
145         
146
147         delete rawReader;
148         delete rawStream;
149
150         timer.Stop();
151         timer.Print();
152 }