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