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