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
22 //____________________________________________________________________
24 // Class to read ADC values from a AliRawReader object.
26 // This class uses the AliFMDRawStreamer class to read the ALTRO
34 // +-----------------+ <<references>> +--------------+
35 // | AliFMDRawReader |<>----------------| AliRawReader |
36 // +-----------------+ +--------------+
40 // +-----------------+ <<uses>> |
41 // | AliFMDRawStream |------------------------+
42 // +-----------------+
49 // #include <AliLog.h> // ALILOG_H
50 #include "AliFMDDebug.h" // Better debug macros
51 #include "AliFMDParameters.h" // ALIFMDPARAMETERS_H
52 #include "AliFMDDigit.h" // ALIFMDDIGIT_H
53 #include "AliFMDRawStream.h" // ALIFMDRAWSTREAM_H
54 // #include "AliRawReader.h" // ALIRAWREADER_H
55 #include "AliFMDRawReader.h" // ALIFMDRAWREADER_H
56 // #include "AliFMDAltroIO.h" // ALIFMDALTROIO_H
57 // #include <TArrayI.h> // ROOT_TArrayI
58 #include <TTree.h> // ROOT_TTree
59 #include <TClonesArray.h> // ROOT_TClonesArray
60 // #include <iostream>
63 //____________________________________________________________________
64 ClassImp(AliFMDRawReader)
66 ; // This is here to keep Emacs for indenting the next line
69 //____________________________________________________________________
70 AliFMDRawReader::AliFMDRawReader(AliRawReader* reader, TTree* tree)
71 : TTask("FMDRawReader", "Reader of Raw ADC values from the FMD"),
79 //____________________________________________________________________
81 AliFMDRawReader::Exec(Option_t*)
84 TClonesArray* array = new TClonesArray("AliFMDDigit");
89 fTree->Branch("FMD", &array);
91 Int_t nWrite = fTree->Fill();
92 AliFMDDebug(1, ("Got a grand total of %d digits, wrote %d bytes to tree",
93 array->GetEntries(), nWrite));
98 //____________________________________________________________________
100 AliFMDRawReader::ReadAdcs(TClonesArray* array)
102 // Read raw data into the digits array, using AliFMDAltroReader.
104 AliError("No TClonesArray passed");
107 // if (!fReader->ReadHeader()) {
108 // AliError("Couldn't read header");
112 AliFMDParameters* pars = AliFMDParameters::Instance();
113 AliFMDRawStream input(fReader);
115 UShort_t stripMin = 0;
116 UShort_t stripMax = 127;
117 UShort_t preSamp = 0;
123 // Data array is approx twice the size needed.
126 Bool_t isGood = kTRUE;
128 isGood = input.ReadChannel(ddl, hwaddr, last, data);
130 AliFMDDebug(5, ("Read channel 0x%x of size %d", hwaddr, last));
131 UShort_t det, sec, str;
133 if (!pars->Hardware2Detector(ddl, hwaddr, det, ring, sec, str)) {
134 AliError(Form("Failed to get detector id from DDL %d "
135 "and hardware address 0x%x", ddl, hwaddr));
138 rate = pars->GetSampleRate(det, ring, sec, str);
139 stripMin = pars->GetMinStrip(det, ring, sec, str);
140 stripMax = pars->GetMaxStrip(det, ring, sec, str);
141 AliFMDDebug(5, ("DDL 0x%04x, address 0x%03x maps to FMD%d%c[%2d,%3d]",
142 ddl, hwaddr, det, ring, sec, str));
144 // Loop over the `timebins', and make the digits
145 for (size_t i = 0; i < last; i++) {
146 if (i < preSamp) continue;
147 Int_t n = array->GetEntriesFast();
148 UShort_t curStr = str + stripMin + i / rate;
149 if ((curStr-str) > stripMax) {
150 AliError(Form("Current strip is %d but DB says max is %d",
153 AliFMDDebug(5, ("making digit for FMD%d%c[%2d,%3d] from sample %4d",
154 det, ring, sec, curStr, i));
155 new ((*array)[n]) AliFMDDigit(det, ring, sec, curStr, data[i],
156 (rate >= 2 ? data[i+1] : 0),
157 (rate >= 3 ? data[i+2] : 0),
158 (rate >= 4 ? data[i+3] : 0));
167 //____________________________________________________________________
169 AliFMDRawReader::ReadAdcs(TClonesArray* array)
171 // Read raw data into the digits array, using AliFMDAltroReader.
173 AliError("No TClonesArray passed");
176 // if (!fReader->ReadHeader()) {
177 // AliError("Couldn't read header");
181 AliFMDParameters* pars = AliFMDParameters::Instance();
184 fReader->Select("FMD");
186 UShort_t stripMin = 0;
187 UShort_t stripMax = 127;
188 UShort_t preSamp = 0;
192 if (!fReader->ReadNextData(cdata)) break;
193 size_t nchar = fReader->GetDataSize();
194 UShort_t ddl = fReader->GetDDLID();
196 AliFMDDebug(1, ("Reading %d bytes (%d 10bit words) from %d",
197 nchar, nchar * 8 / 10, ddl));
198 // Make a stream to read from
199 std::string str((char*)(cdata), nchar);
200 std::istringstream s(str);
201 // Prep the reader class.
202 AliFMDAltroReader r(s);
203 // Data array is approx twice the size needed.
204 UShort_t data[2048], hwaddr, last;
205 while (r.ReadChannel(hwaddr, last, data) > 0) {
206 AliFMDDebug(5, ("Read channel 0x%x of size %d", hwaddr, last));
207 UShort_t det, sec, str;
209 if (!pars->Hardware2Detector(ddl, hwaddr, det, ring, sec, str)) {
210 AliError(Form("Failed to detector id from DDL %d "
211 "and hardware address 0x%x", ddl, hwaddr));
214 rate = pars->GetSampleRate(det, ring, sec, str);
215 stripMin = pars->GetMinStrip(det, ring, sec, str);
216 stripMax = pars->GetMaxStrip(det, ring, sec, str);
217 AliFMDDebug(5, ("DDL 0x%04x, address 0x%03x maps to FMD%d%c[%2d,%3d]",
218 ddl, hwaddr, det, ring, sec, str));
220 // Loop over the `timebins', and make the digits
221 for (size_t i = 0; i < last; i++) {
222 if (i < preSamp) continue;
223 Int_t n = array->GetEntries();
224 UShort_t curStr = str + stripMin + i / rate;
225 if ((curStr-str) > stripMax) {
226 AliError(Form("Current strip is %d but DB says max is %d",
229 AliFMDDebug(5, ("making digit for FMD%d%c[%2d,%3d] from sample %4d",
230 det, ring, sec, curStr, i));
231 new ((*array)[n]) AliFMDDigit(det, ring, sec, curStr, data[i],
232 (rate >= 2 ? data[i+1] : 0),
233 (rate >= 3 ? data[i+2] : 0));
237 if (r.IsBof()) break;
245 // This is the old method, for comparison. It's really ugly, and far
247 //____________________________________________________________________
249 AliFMDRawReader::Exec(Option_t*)
251 // Read raw data into the digits array
252 // if (!fReader->ReadHeader()) {
253 // Error("ReadAdcs", "Couldn't read header");
258 TClonesArray* array = new TClonesArray("AliFMDDigit");
259 fTree->Branch("FMD", &array);
262 AliFMDParameters* pars = AliFMDParameters::Instance();
263 fSampleRate = pars->GetSampleRate(0);
265 // Use AliAltroRawStream to read the ALTRO format. No need to
266 // reinvent the wheel :-)
267 AliFMDRawStream input(fReader, fSampleRate);
269 fReader->Select("FMD");
273 UShort_t detector = 1; // Must be one here
274 UShort_t oldDetector = 0;
281 // Loop over data in file
286 Int_t ddl = fReader->GetDDLID();
287 AliFMDDebug(10, ("Current DDL is %d", ddl));
288 if (ddl != oldDDL || input.IsNewStrip() || !next) {
289 // Make a new digit, if we have some data (oldDetector == 0,
290 // means that we haven't really read anything yet - that is,
291 // it's the first time we get here).
292 if (oldDetector > 0) {
294 AliFMDDebug(10, ("Add a new strip: FMD%d%c[%2d,%3d] "
295 "(current: FMD%d%c[%2d,%3d])",
296 oldDetector, input.PrevRing(),
297 input.PrevSector() , input.PrevStrip(),
298 detector , input.Ring(), input.Sector(),
300 new ((*array)[n]) AliFMDDigit(oldDetector,
304 counts[0], counts[1], counts[2]);
308 static_cast<AliFMDDigit*>(fFMD->Digits()->
309 UncheckedAt(fFMD->GetNdigits()-1));
314 AliFMDDebug(10, ("Read %d channels for FMD%d",
315 count + 1, detector));
320 // If we got a new DDL, it means we have a new detector.
323 AliFMDDebug(10, ("Read %d channels for FMD%d", count + 1, detector));
324 // Reset counts, and update the DDL cache
327 // Check that we're processing a FMD detector
328 Int_t detId = fReader->GetDetectorID();
329 if (detId != (AliDAQ::DetectorID("FMD"))) {
330 AliError(Form("Detector ID %d != %d",
331 detId, (AliDAQ::DetectorID("FMD"))));
334 // Figure out what detector we're deling with
335 oldDetector = detector;
337 case 0: detector = 1; break;
338 case 1: detector = 2; break;
339 case 2: detector = 3; break;
341 AliError(Form("Unknown DDL 0x%x for FMD", ddl));
344 AliFMDDebug(10, ("Reading ADCs for 0x%x - That is FMD%d",
345 fReader->GetEquipmentId(), detector));
350 counts[input.Sample()] = input.Count();
352 AliFMDDebug(10, ("ADC of FMD%d%c[%2d,%3d] += %d",
353 detector, input.Ring(), input.Sector(),
354 input.Strip(), input.Count()));
355 oldDetector = detector;
363 //____________________________________________________________________