1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
16 /** @file AliFMDRawReader.cxx
17 @author Christian Holm Christensen <cholm@nbi.dk>
18 @date Mon Mar 27 12:45:23 2006
19 @brief Class to read raw data
21 //____________________________________________________________________
23 // Class to read ADC values from a AliRawReader object.
25 // This class uses the AliFMDRawStreamer class to read the ALTRO
33 // +-----------------+ <<references>> +--------------+
34 // | AliFMDRawReader |<>----------------| AliRawReader |
35 // +-----------------+ +--------------+
39 // +-----------------+ <<uses>> |
40 // | AliFMDRawStream |------------------------+
41 // +-----------------+
48 #include <AliLog.h> // ALILOG_H
49 #include "AliFMDParameters.h" // ALIFMDPARAMETERS_H
50 #include "AliFMDDigit.h" // ALIFMDDIGIT_H
51 #include "AliFMDRawStream.h" // ALIFMDRAWSTREAM_H
52 #include "AliRawReader.h" // ALIRAWREADER_H
53 #include "AliFMDRawReader.h" // ALIFMDRAWREADER_H
54 // #include "AliFMDAltroIO.h" // ALIFMDALTROIO_H
55 #include <TArrayI.h> // ROOT_TArrayI
56 #include <TTree.h> // ROOT_TTree
57 #include <TClonesArray.h> // ROOT_TClonesArray
61 #define PRETTY_HEX(N,X) \
62 " 0x" << std::setfill('0') << std::setw(N) << std::hex << X \
63 << std::setfill(' ') << std::dec
65 //____________________________________________________________________
66 ClassImp(AliFMDRawReader)
68 ; // This is here to keep Emacs for indenting the next line
71 //____________________________________________________________________
72 AliFMDRawReader::AliFMDRawReader(AliRawReader* reader, TTree* tree)
73 : TTask("FMDRawReader", "Reader of Raw ADC values from the FMD"),
81 //____________________________________________________________________
83 AliFMDRawReader::Exec(Option_t*)
85 TClonesArray* array = new TClonesArray("AliFMDDigit");
90 fTree->Branch("FMD", &array);
92 Int_t nWrite = fTree->Fill();
93 AliDebug(1, Form("Got a grand total of %d digits, wrote %d bytes to tree",
94 array->GetEntries(), nWrite));
99 //____________________________________________________________________
101 AliFMDRawReader::ReadAdcs(TClonesArray* array)
103 // Read raw data into the digits array, using AliFMDAltroReader.
105 AliError("No TClonesArray passed");
108 if (!fReader->ReadHeader()) {
109 AliError("Couldn't read header");
113 AliFMDParameters* pars = AliFMDParameters::Instance();
114 AliFMDRawStream input(fReader);
116 fReader->Select(AliFMDParameters::kBaseDDL>>8);
118 UShort_t stripMin = 0;
119 UShort_t stripMax = 127;
120 UShort_t preSamp = 0;
126 // Data array is approx twice the size needed.
128 while (input.ReadChannel(ddl, hwaddr, last, data)) {
129 AliDebug(5, Form("Read channel 0x%x of size %d", hwaddr, last));
130 UShort_t det, sec, str;
132 if (!pars->Hardware2Detector(ddl, hwaddr, det, ring, sec, str)) {
133 AliError(Form("Failed to get detector id from DDL %d "
134 "and hardware address 0x%x", ddl, hwaddr));
137 rate = pars->GetSampleRate(det, ring, sec, str);
138 stripMin = pars->GetMinStrip(det, ring, sec, str);
139 stripMax = pars->GetMaxStrip(det, ring, sec, str);
140 AliDebug(5, Form("DDL 0x%04x, address 0x%03x maps to FMD%d%c[%2d,%3d]",
141 ddl, hwaddr, det, ring, sec, str));
143 // Loop over the `timebins', and make the digits
144 for (size_t i = 0; i < last; i++) {
145 if (i < preSamp) continue;
146 Int_t n = array->GetEntries();
147 UShort_t curStr = str + stripMin + i / rate;
148 if ((curStr-str) > stripMax) {
149 AliError(Form("Current strip is %d but DB says max is %d",
152 AliDebug(5, Form("making digit for FMD%d%c[%2d,%3d] from sample %4d",
153 det, ring, sec, curStr, i));
154 new ((*array)[n]) AliFMDDigit(det, ring, sec, curStr, data[i],
155 (rate >= 2 ? data[i+1] : 0),
156 (rate >= 3 ? data[i+2] : 0));
164 //____________________________________________________________________
166 AliFMDRawReader::ReadAdcs(TClonesArray* array)
168 // Read raw data into the digits array, using AliFMDAltroReader.
170 AliError("No TClonesArray passed");
173 if (!fReader->ReadHeader()) {
174 AliError("Couldn't read header");
178 AliFMDParameters* pars = AliFMDParameters::Instance();
181 fReader->Select(AliFMDParameters::kBaseDDL>>8);
183 UShort_t stripMin = 0;
184 UShort_t stripMax = 127;
185 UShort_t preSamp = 0;
189 if (!fReader->ReadNextData(cdata)) break;
190 size_t nchar = fReader->GetDataSize();
191 UShort_t ddl = AliFMDParameters::kBaseDDL + fReader->GetDDLID();
193 AliDebug(1, Form("Reading %d bytes (%d 10bit words) from %d",
194 nchar, nchar * 8 / 10, ddl));
195 // Make a stream to read from
196 std::string str((char*)(cdata), nchar);
197 std::istringstream s(str);
198 // Prep the reader class.
199 AliFMDAltroReader r(s);
200 // Data array is approx twice the size needed.
201 UShort_t data[2048], hwaddr, last;
202 while (r.ReadChannel(hwaddr, last, data) > 0) {
203 AliDebug(5, Form("Read channel 0x%x of size %d", hwaddr, last));
204 UShort_t det, sec, str;
206 if (!pars->Hardware2Detector(ddl, hwaddr, det, ring, sec, str)) {
207 AliError(Form("Failed to detector id from DDL %d "
208 "and hardware address 0x%x", ddl, hwaddr));
211 rate = pars->GetSampleRate(det, ring, sec, str);
212 stripMin = pars->GetMinStrip(det, ring, sec, str);
213 stripMax = pars->GetMaxStrip(det, ring, sec, str);
214 AliDebug(5, Form("DDL 0x%04x, address 0x%03x maps to FMD%d%c[%2d,%3d]",
215 ddl, hwaddr, det, ring, sec, str));
217 // Loop over the `timebins', and make the digits
218 for (size_t i = 0; i < last; i++) {
219 if (i < preSamp) continue;
220 Int_t n = array->GetEntries();
221 UShort_t curStr = str + stripMin + i / rate;
222 if ((curStr-str) > stripMax) {
223 AliError(Form("Current strip is %d but DB says max is %d",
226 AliDebug(5, Form("making digit for FMD%d%c[%2d,%3d] from sample %4d",
227 det, ring, sec, curStr, i));
228 new ((*array)[n]) AliFMDDigit(det, ring, sec, curStr, data[i],
229 (rate >= 2 ? data[i+1] : 0),
230 (rate >= 3 ? data[i+2] : 0));
234 if (r.IsBof()) break;
242 // This is the old method, for comparison. It's really ugly, and far
244 //____________________________________________________________________
246 AliFMDRawReader::Exec(Option_t*)
248 // Read raw data into the digits array
249 if (!fReader->ReadHeader()) {
250 Error("ReadAdcs", "Couldn't read header");
255 TClonesArray* array = new TClonesArray("AliFMDDigit");
256 fTree->Branch("FMD", &array);
259 AliFMDParameters* pars = AliFMDParameters::Instance();
260 fSampleRate = pars->GetSampleRate(AliFMDParameters::kBaseDDL);
262 // Use AliAltroRawStream to read the ALTRO format. No need to
263 // reinvent the wheel :-)
264 AliFMDRawStream input(fReader, fSampleRate);
266 fReader->Select(AliFMDParameters::kBaseDDL);
270 UShort_t detector = 1; // Must be one here
271 UShort_t oldDetector = 0;
278 // Loop over data in file
283 Int_t ddl = fReader->GetDDLID();
284 AliDebug(10, Form("Current DDL is %d", ddl));
285 if (ddl != oldDDL || input.IsNewStrip() || !next) {
286 // Make a new digit, if we have some data (oldDetector == 0,
287 // means that we haven't really read anything yet - that is,
288 // it's the first time we get here).
289 if (oldDetector > 0) {
291 AliDebug(10, Form("Add a new strip: FMD%d%c[%2d,%3d] "
292 "(current: FMD%d%c[%2d,%3d])",
293 oldDetector, input.PrevRing(),
294 input.PrevSector() , input.PrevStrip(),
295 detector , input.Ring(), input.Sector(),
297 new ((*array)[n]) AliFMDDigit(oldDetector,
301 counts[0], counts[1], counts[2]);
305 static_cast<AliFMDDigit*>(fFMD->Digits()->
306 UncheckedAt(fFMD->GetNdigits()-1));
311 AliDebug(10, Form("Read %d channels for FMD%d",
312 count + 1, detector));
317 // If we got a new DDL, it means we have a new detector.
320 AliDebug(10, Form("Read %d channels for FMD%d", count + 1, detector));
321 // Reset counts, and update the DDL cache
324 // Check that we're processing a FMD detector
325 Int_t detId = fReader->GetDetectorID();
326 if (detId != (AliFMDParameters::kBaseDDL >> 8)) {
327 AliError(Form("Detector ID %d != %d",
328 detId, (AliFMDParameters::kBaseDDL >> 8)));
331 // Figure out what detector we're deling with
332 oldDetector = detector;
334 case 0: detector = 1; break;
335 case 1: detector = 2; break;
336 case 2: detector = 3; break;
338 AliError(Form("Unknown DDL 0x%x for FMD", ddl));
341 AliDebug(10, Form("Reading ADCs for 0x%x - That is FMD%d",
342 fReader->GetEquipmentId(), detector));
347 counts[input.Sample()] = input.Count();
349 AliDebug(10, Form("ADC of FMD%d%c[%2d,%3d] += %d",
350 detector, input.Ring(), input.Sector(),
351 input.Strip(), input.Count()));
352 oldDetector = detector;
360 //____________________________________________________________________