]> git.uio.no Git - u/mrichter/AliRoot.git/blame - T0/AliT0RawReader.cxx
AliHMPIDDigitN no longer needed
[u/mrichter/AliRoot.git] / T0 / AliT0RawReader.cxx
CommitLineData
dc7ca31d 1#include "AliT0RawReader.h"
2#include "AliT0RawData.h"
3#include "AliT0digit.h"
4#include "AliBitPacking.h"
5#include "TBits.h"
6
7#include <Riostream.h>
8#include "TMath.h"
9#include "TH1F.h"
10#include "TArrayI.h"
11#include "AliLog.h"
12
13ClassImp(AliT0RawReader)
14
15 AliT0RawReader::AliT0RawReader (AliRawReader *rawReader, TTree* tree)
16 : TTask("T0RawReader","read raw T0 data"),
17 fDigits(NULL),
18 fTree(tree),
19 fRawReader(rawReader),
20 fData(NULL),
21 fPosition(0)
22{
23 //
24// create an object to read T0raw digits
25 AliDebug(1,"Start ");
26 if (fDigits == 0x0) fDigits = new AliT0digit();
27 fTree->Branch("T0","AliT0digit",&fDigits,405,1);
28
29 fRawReader->Reset();
30 fRawReader->Select("T0");
31
32
33}
34 AliT0RawReader::~AliT0RawReader ()
35{
36 //
37}
38
39Bool_t AliT0RawReader::Next()
40{
41// read the next raw digit
42// returns kFALSE if there is no digit left
43//"LookUpTable":
44// Amplitude LED TRM=0; chain=0; TDC 0 -5 channel 0,2,4,6
45// Time CFD TRM=0; chain=0; TDC 6 - 11 channel 0,2,4,6
46// mean time TRM=0; chain=0; TDC 12 channel 0
47// T0A TRM=0; chain=0; TDC 12 channel 2
48// T0C TRM=0; chain=0; TDC 12 channel 4
49// vertex TRM=0; chain=0; TDC 12 channel 6
50// mult QTC0 TRM=0; chain=0; TDC 13 channel 0
51// mult QTC1 TRM=0; chain=0; TDC 13 channel 2
52
53// Charge QTC0 TRM=1; chain=0; TDC 0 -5 channel 0,2,4,6
54// Charge QTC1 TRM=1; chain=0; TDC 6 - 11 channel 0,2,4,6
55// T0A trigger TRM=1; chain=0; TDC 12 channel 0
56// T0C trigger TRM=1; chain=0; TDC 12 channel 2
57// vertex trigger TRM=1; chain=0; TDC 12 channel 4
58// trigger central TRM=1; chain=0; TDC 13 channel 0
59// tigger semicenral TRM=1; chain=0; TDC 13 channel 2
60//
61// allData array collect data from all channels in one :
62// allData[0] - allData[23] 24 CFD channels
63// allData[24] - allData[47] 24 LED channels
64// allData[48] mean (T0) signal
65// allData[49] time difference (vertex)
66
67 UInt_t word;
68 Int_t time=0, itdc=0, ichannel=0;
69 Int_t numberOfWordsInTRM=0, iTRM=0;
70 Int_t tdcTime, koef, meanTime, timeDiff ;
71 Int_t allData[107];
72
73 TArrayI *timeTDC1 = new TArrayI(24);
74 TArrayI * chargeTDC1 = new TArrayI(24);
75 TArrayI *timeTDC2 = new TArrayI(24);
76 TArrayI *chargeTDC2 = new TArrayI(24);
77
78 for ( Int_t k=0; k<107; k++) allData[k]=0;
79 do {
80 if (!fRawReader->ReadNextData(fData)) return kFALSE;
81 } while (fRawReader->GetDataSize() == 0);
82
83 // fPosition = GetPosition();
84 fPosition = 0;
85
86 //DRM header
87 for (Int_t i=0; i<4; i++) {
88 word = GetNextWord();
89 }
90 //TRMheader
91 word = GetNextWord();
92 numberOfWordsInTRM=AliBitPacking::UnpackWord(word,4,16);
93 iTRM=AliBitPacking::UnpackWord(word,0,3);
94
95 //chain header
96 word = GetNextWord();
97
98 for (Int_t i=0; i<numberOfWordsInTRM; i++) {
99 word = GetNextWord();
100 tdcTime = AliBitPacking::UnpackWord(word,31,31);
101
102 if ( tdcTime == 1)
103 {
104 itdc=AliBitPacking::UnpackWord(word,24,27);
105 ichannel=AliBitPacking::UnpackWord(word,21,23);
106 time=AliBitPacking::UnpackWord(word,0,20);
107 koef = itdc*4 + ichannel/2;
108 allData[koef]=time;
109 }
110 }
111 word = GetNextWord(); //chain trailer
112 word = GetNextWord(); //TRM trailer
113
114 //TRMheader
115 word = GetNextWord();
116 numberOfWordsInTRM=AliBitPacking::UnpackWord(word,4,16);
117 iTRM=AliBitPacking::UnpackWord(word,0,3);
118 //chain header
119 word = GetNextWord();
120
121 for (Int_t iword=0; iword<numberOfWordsInTRM; iword++) {
122 word = GetNextWord();
123 tdcTime = AliBitPacking::UnpackWord(word,31,31);
124
125 if ( tdcTime == 1)
126 {
127 itdc=AliBitPacking::UnpackWord(word,24,27);
128 ichannel=AliBitPacking::UnpackWord(word,21,23);
129 time=AliBitPacking::UnpackWord(word,0,20);
130 koef = itdc*4 + ichannel/2;
131 allData[koef+54]=time;
132 }
133 }
134
135 for (Int_t in=0; in<24; in++)
136 {
137 timeTDC1->AddAt(allData[in],in);
138 timeTDC2->AddAt(allData[in+24],in);
139 chargeTDC1->AddAt(allData[in+54],in);
140 chargeTDC2->AddAt(allData[in+78],in);
141 }
142
143 meanTime = allData[48]; // T0 !!!!!!
144 timeDiff = allData[49];
145
146 word = GetNextWord();
147 word = GetNextWord();
148
149 fDigits->SetTime(*timeTDC2);
150 fDigits->SetADC(*chargeTDC1);
151
152 fDigits->SetTimeAmp(*timeTDC1);
153 fDigits->SetADCAmp(*chargeTDC2);
154
155 fDigits->SetMeanTime(meanTime);
156 fDigits->SetDiffTime(timeDiff);
157 fTree->Fill();
158
159 delete timeTDC1 ;
160 delete chargeTDC1;
161 delete timeTDC2 ;
162 delete chargeTDC2;
163
164 return kTRUE;
165}
166//_____________________________________________________________________________
167/*
168void AliT0RawReader::UnpackTime(Int_t outTime, Int_t outCh)
169{
170 UInt_t word=0;
171 UInt_t unpackword=0;
172
173 word = GetNextWord();
174 unpackword=AliBitPacking::UnpackWord(word,0,12);
175 outTime=unpackword;
176 unpackword= AliBitPacking::UnpackWord(word,21,27);
177 outCh=unpackword;
178 }
179 */
180//_____________________________________________________________________________
181Int_t AliT0RawReader::GetPosition()
182{
183 // Sets the position in the
184 // input stream
185 if (((fRawReader->GetDataSize() * 8) % 32) != 0)
186 AliFatal(Form("Incorrect raw data size ! %d words are found !",fRawReader->GetDataSize()));
187 return (fRawReader->GetDataSize() * 8) / 32;
188}
189//_____________________________________________________________________________
190UInt_t AliT0RawReader::GetNextWord()
191{
192 // Read the next 32 bit word in backward direction
193 // The input stream access is given by fData and fPosition
194
195
196 // fPosition--;
197 Int_t iBit = fPosition * 32;
198 Int_t iByte = iBit / 8;
199
200 UInt_t word = 0;
201 word = fData[iByte+3]<<24;
202 word |= fData[iByte+2]<<16;
203 word |= fData[iByte+1]<<8;
204 word |= fData[iByte];
205 fPosition++;
206
207 return word;
208
209}
210