]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDRawReader.cxx
Minor fixes
[u/mrichter/AliRoot.git] / FMD / AliFMDRawReader.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15 /* $Id$ */
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 
20 */
21 //____________________________________________________________________
22 //
23 // Class to read ADC values from a AliRawReader object. 
24 //
25 // This class uses the AliFMDRawStreamer class to read the ALTRO
26 // formatted data. 
27 // 
28 //          +-------+
29 //          | TTask |
30 //          +-------+
31 //              ^
32 //              |
33 //      +-----------------+  <<references>>  +--------------+
34 //      | AliFMDRawReader |<>----------------| AliRawReader |
35 //      +-----------------+                  +--------------+
36 //              |                                  ^
37 //              | <<uses>>                         |
38 //              V                                  |
39 //      +-----------------+      <<uses>>          |
40 //      | AliFMDRawStream |------------------------+
41 //      +-----------------+
42 //              |
43 //              V
44 //      +----------------+
45 //      | AliAltroStream |
46 //      +----------------+
47 //
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
58 #include <iostream>
59 #include <iomanip>
60 #include <sstream>
61 #define PRETTY_HEX(N,X) \
62   "  0x" << std::setfill('0') << std::setw(N) << std::hex << X \
63          << std::setfill(' ') << std::dec
64
65 //____________________________________________________________________
66 ClassImp(AliFMDRawReader)
67 #if 0
68   ; // This is here to keep Emacs for indenting the next line
69 #endif
70
71 //____________________________________________________________________
72 AliFMDRawReader::AliFMDRawReader(AliRawReader* reader, TTree* tree) 
73   : TTask("FMDRawReader", "Reader of Raw ADC values from the FMD"),
74     fTree(tree),
75     fReader(reader), 
76     fSampleRate(1)
77 {
78   // Default CTOR
79 }
80
81 //____________________________________________________________________
82 void
83 AliFMDRawReader::Exec(Option_t*) 
84 {
85   TClonesArray* array = new TClonesArray("AliFMDDigit");
86   if (!fTree) {
87     AliError("No tree");
88     return;
89   }
90   fTree->Branch("FMD", &array);
91   ReadAdcs(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));
95 }
96
97
98 #if 1
99 //____________________________________________________________________
100 Bool_t
101 AliFMDRawReader::ReadAdcs(TClonesArray* array) 
102 {
103   // Read raw data into the digits array, using AliFMDAltroReader. 
104   if (!array) {
105     AliError("No TClonesArray passed");
106     return kFALSE;
107   }
108   if (!fReader->ReadHeader()) {
109     AliError("Couldn't read header");
110     return kFALSE;
111   }
112   // Get sample rate 
113   AliFMDParameters* pars = AliFMDParameters::Instance();
114   AliFMDRawStream input(fReader);
115
116   UShort_t stripMin = 0;
117   UShort_t stripMax = 127;
118   UShort_t preSamp  = 0;
119   
120   UInt_t ddl    = 0;
121   UInt_t rate   = 0;
122   UInt_t last   = 0;
123   UInt_t hwaddr = 0;
124   // Data array is approx twice the size needed. 
125   UShort_t data[2048];
126   while (input.ReadChannel(ddl, hwaddr, last, data)) {
127     AliDebug(5, Form("Read channel 0x%x of size %d", hwaddr, last));
128     UShort_t det, sec, str;
129     Char_t   ring;
130     if (!pars->Hardware2Detector(ddl, hwaddr, det, ring, sec, str)) {
131       AliError(Form("Failed to get detector id from DDL %d "
132                     "and hardware address 0x%x", ddl, hwaddr));
133       continue;
134     }
135     rate     = pars->GetSampleRate(det, ring, sec, str);
136     stripMin = pars->GetMinStrip(det, ring, sec, str);
137     stripMax = pars->GetMaxStrip(det, ring, sec, str);
138     AliDebug(5, Form("DDL 0x%04x, address 0x%03x maps to FMD%d%c[%2d,%3d]", 
139                        ddl, hwaddr, det, ring, sec, str));
140     
141     // Loop over the `timebins', and make the digits
142     for (size_t i = 0; i < last; i++) {
143       if (i < preSamp) continue;
144       Int_t    n      = array->GetEntries();
145       UShort_t curStr = str + stripMin + i / rate;
146       if ((curStr-str) > stripMax) {
147         AliError(Form("Current strip is %d but DB says max is %d", 
148                       curStr, stripMax));
149       }
150       AliDebug(5, Form("making digit for FMD%d%c[%2d,%3d] from sample %4d", 
151                        det, ring, sec, curStr, i));
152       new ((*array)[n]) AliFMDDigit(det, ring, sec, curStr, data[i], 
153                                     (rate >= 2 ? data[i+1] : 0),
154                                     (rate >= 3 ? data[i+2] : 0));
155       if (rate >= 2) i++;
156       if (rate >= 3) i++;
157     }
158   }
159   return kTRUE;
160 }
161 #else
162 //____________________________________________________________________
163 Bool_t
164 AliFMDRawReader::ReadAdcs(TClonesArray* array) 
165 {
166   // Read raw data into the digits array, using AliFMDAltroReader. 
167   if (!array) {
168     AliError("No TClonesArray passed");
169     return kFALSE;
170   }
171   if (!fReader->ReadHeader()) {
172     AliError("Couldn't read header");
173     return kFALSE;
174   }
175   // Get sample rate 
176   AliFMDParameters* pars = AliFMDParameters::Instance();
177
178   // Select FMD DDL's 
179   fReader->Select(AliFMDParameters::kBaseDDL>>8);
180
181   UShort_t stripMin = 0;
182   UShort_t stripMax = 127;
183   UShort_t preSamp  = 0;
184   
185   do {
186     UChar_t* cdata;
187     if (!fReader->ReadNextData(cdata)) break;
188     size_t   nchar = fReader->GetDataSize();
189     UShort_t ddl   = AliFMDParameters::kBaseDDL + fReader->GetDDLID();
190     UShort_t rate  = 0;
191     AliDebug(1, Form("Reading %d bytes (%d 10bit words) from %d", 
192                      nchar, nchar * 8 / 10, ddl));
193     // Make a stream to read from 
194     std::string str((char*)(cdata), nchar);
195     std::istringstream s(str);
196     // Prep the reader class.
197     AliFMDAltroReader r(s);
198     // Data array is approx twice the size needed. 
199     UShort_t data[2048], hwaddr, last;
200     while (r.ReadChannel(hwaddr, last, data) > 0) {
201       AliDebug(5, Form("Read channel 0x%x of size %d", hwaddr, last));
202       UShort_t det, sec, str;
203       Char_t   ring;
204       if (!pars->Hardware2Detector(ddl, hwaddr, det, ring, sec, str)) {
205         AliError(Form("Failed to detector id from DDL %d "
206                       "and hardware address 0x%x", ddl, hwaddr));
207         continue;
208       }
209       rate     = pars->GetSampleRate(det, ring, sec, str);
210       stripMin = pars->GetMinStrip(det, ring, sec, str);
211       stripMax = pars->GetMaxStrip(det, ring, sec, str);
212       AliDebug(5, Form("DDL 0x%04x, address 0x%03x maps to FMD%d%c[%2d,%3d]", 
213                        ddl, hwaddr, det, ring, sec, str));
214
215       // Loop over the `timebins', and make the digits
216       for (size_t i = 0; i < last; i++) {
217         if (i < preSamp) continue;
218         Int_t    n      = array->GetEntries();
219         UShort_t curStr = str + stripMin + i / rate;
220         if ((curStr-str) > stripMax) {
221           AliError(Form("Current strip is %d but DB says max is %d", 
222                         curStr, stripMax));
223         }
224         AliDebug(5, Form("making digit for FMD%d%c[%2d,%3d] from sample %4d", 
225                          det, ring, sec, curStr, i));
226         new ((*array)[n]) AliFMDDigit(det, ring, sec, curStr, data[i], 
227                                       (rate >= 2 ? data[i+1] : 0),
228                                       (rate >= 3 ? data[i+2] : 0));
229         if (rate >= 2) i++;
230         if (rate >= 3) i++;
231         }
232         if (r.IsBof()) break;
233     }
234   } while (true);
235   return kTRUE;
236 }
237
238   
239
240 // This is the old method, for comparison.   It's really ugly, and far
241 // too convoluted. 
242 //____________________________________________________________________
243 void
244 AliFMDRawReader::Exec(Option_t*) 
245 {
246   // Read raw data into the digits array
247   if (!fReader->ReadHeader()) {
248     Error("ReadAdcs", "Couldn't read header");
249     return;
250   }
251
252   Int_t n = 0;
253   TClonesArray* array = new TClonesArray("AliFMDDigit");
254   fTree->Branch("FMD", &array);
255
256   // Get sample rate 
257   AliFMDParameters* pars = AliFMDParameters::Instance();
258   fSampleRate = pars->GetSampleRate(AliFMDParameters::kBaseDDL);
259
260   // Use AliAltroRawStream to read the ALTRO format.  No need to
261   // reinvent the wheel :-) 
262   AliFMDRawStream input(fReader, fSampleRate);
263   // Select FMD DDL's 
264   fReader->Select(AliFMDParameters::kBaseDDL);
265   
266   Int_t    oldDDL      = -1;
267   Int_t    count       = 0;
268   UShort_t detector    = 1; // Must be one here
269   UShort_t oldDetector = 0;
270   Bool_t   next        = kTRUE;
271
272   // local Cache 
273   TArrayI counts(10);
274   counts.Reset(-1);
275   
276   // Loop over data in file 
277   while (next) {
278     next = input.Next();
279
280     count++; 
281     Int_t ddl = fReader->GetDDLID();
282     AliDebug(10, Form("Current DDL is %d", ddl));
283     if (ddl != oldDDL || input.IsNewStrip() || !next) {
284       // Make a new digit, if we have some data (oldDetector == 0,
285       // means that we haven't really read anything yet - that is,
286       // it's the first time we get here). 
287       if (oldDetector > 0) {
288         // Got a new strip. 
289         AliDebug(10, Form("Add a new strip: FMD%d%c[%2d,%3d] "
290                           "(current: FMD%d%c[%2d,%3d])", 
291                           oldDetector, input.PrevRing(), 
292                           input.PrevSector() , input.PrevStrip(),
293                           detector , input.Ring(), input.Sector(), 
294                           input.Strip()));
295         new ((*array)[n]) AliFMDDigit(oldDetector, 
296                                       input.PrevRing(), 
297                                       input.PrevSector(), 
298                                       input.PrevStrip(), 
299                                       counts[0], counts[1], counts[2]);
300         n++;
301 #if 0
302         AliFMDDigit* digit = 
303           static_cast<AliFMDDigit*>(fFMD->Digits()->
304                                     UncheckedAt(fFMD->GetNdigits()-1));
305 #endif 
306       }
307         
308       if (!next) { 
309         AliDebug(10, Form("Read %d channels for FMD%d", 
310                           count + 1, detector));
311         break;
312       }
313     
314     
315       // If we got a new DDL, it means we have a new detector. 
316       if (ddl != oldDDL) {
317         if (detector != 0) 
318           AliDebug(10, Form("Read %d channels for FMD%d", count + 1, detector));
319         // Reset counts, and update the DDL cache 
320         count       = 0;
321         oldDDL      = ddl;
322         // Check that we're processing a FMD detector 
323         Int_t detId = fReader->GetDetectorID();
324         if (detId != (AliFMDParameters::kBaseDDL >> 8)) {
325           AliError(Form("Detector ID %d != %d",
326                         detId, (AliFMDParameters::kBaseDDL >> 8)));
327           break;
328         }
329         // Figure out what detector we're deling with 
330         oldDetector = detector;
331         switch (ddl) {
332         case 0: detector = 1; break;
333         case 1: detector = 2; break;
334         case 2: detector = 3; break;
335         default:
336           AliError(Form("Unknown DDL 0x%x for FMD", ddl));
337           return;
338         }
339         AliDebug(10, Form("Reading ADCs for 0x%x  - That is FMD%d",
340                           fReader->GetEquipmentId(), detector));
341       }
342       counts.Reset(-1);
343     }
344     
345     counts[input.Sample()] = input.Count();
346     
347     AliDebug(10, Form("ADC of FMD%d%c[%2d,%3d] += %d",
348                       detector, input.Ring(), input.Sector(), 
349                       input.Strip(), input.Count()));
350     oldDetector = detector;
351   }
352   fTree->Fill();
353   return;
354
355 }
356 #endif
357
358 //____________________________________________________________________
359 // 
360 // EOF
361 //