]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDRawReader.cxx
Introducing the fixes needed in order to read the old RCU firmware format (v<3) with...
[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 "AliFMDSDigit.h"       // ALIFMDSDIGIT_H
54 #include "AliFMDRawStream.h"    // ALIFMDRAWSTREAM_H 
55 #include "AliRawReader.h"       // ALIRAWREADER_H 
56 #include "AliFMDRawReader.h"    // ALIFMDRAWREADER_H 
57 #include "AliFMDDebug.h"
58 #include "AliFMDCalibSampleRate.h"
59 #include "AliFMDCalibStripRange.h"
60 #include "AliFMDAltroMapping.h"
61 // #include "AliFMDAltroIO.h"   // ALIFMDALTROIO_H 
62 #include <TArrayS.h>            // ROOT_TArrayS
63 #include <TTree.h>              // ROOT_TTree
64 #include <TClonesArray.h>       // ROOT_TClonesArray
65 #include <TString.h>
66 #include <iostream>
67 #include <climits>
68 // #include <iomanip>
69
70 //____________________________________________________________________
71 ClassImp(AliFMDRawReader)
72 #if 0
73   ; // This is here to keep Emacs for indenting the next line
74 #endif
75
76 //____________________________________________________________________
77 AliFMDRawReader::AliFMDRawReader(AliRawReader* reader, TTree* tree) 
78   : TTask("FMDRawReader", "Reader of Raw ADC values from the FMD"),
79     fTree(tree),
80     fReader(reader), 
81     // fSampleRate(1),
82     fData(0),
83     fNbytes(0), 
84     fSeen()
85 {
86   // Default CTOR
87   for (Int_t i = 0; i < 3; i++) { 
88     fSampleRate[i]   = 0;
89     fZeroSuppress[i] = kFALSE;
90     fNoiseFactor[i]  = 1;
91   }
92 }
93
94 //____________________________________________________________________
95 void
96 AliFMDRawReader::Exec(Option_t*) 
97 {
98   // Read the data 
99   TClonesArray* array = new TClonesArray("AliFMDDigit");
100   if (!fTree) {
101     AliError("No tree");
102     return;
103   }
104   fTree->Branch("FMD", &array);
105   
106   
107   ReadAdcs(array);
108   Int_t nWrite = fTree->Fill();
109   AliFMDDebug(1, ("Got a grand total of %d digits, wrote %d bytes to tree", 
110                    array->GetEntriesFast(), nWrite));
111 }
112
113 //____________________________________________________________________
114 Bool_t
115 AliFMDRawReader::NextSample(UShort_t& det, Char_t&   rng, UShort_t& sec, 
116                             UShort_t& str, UShort_t& sam, UShort_t& rat, 
117                             Short_t&  adc, Bool_t&   zs,  UShort_t& fac)
118 {
119   // Scan current event for next signal.   It returns kFALSE when
120   // there's no more data in the event. 
121   static AliAltroRawStream   stream(fReader); //    = 0;
122   static AliFMDParameters*   pars     = 0;
123   static AliFMDAltroMapping* map      = 0;
124   static Int_t               ddl      = -1;
125   // static UInt_t              rate     = 0;
126   static UShort_t            tdet     = 0;
127   static Char_t              trng     = '\0';
128   static UShort_t            tsec     = 0;
129   static Short_t             tstr     = 0;   
130   static Short_t             bstr     = -1;
131   static Short_t             tsam     = -1;   
132   static UInt_t              trate    = 0;
133   static Int_t               hwaddr   = -1;
134   static UShort_t            stripMin = 0;
135   static UShort_t            stripMax = 0; // 127;
136   static UShort_t            preSamp  = 0; // 14+5;
137   if (stream.GetDDLNumber() < 0) { 
138     fReader->Select("FMD");
139
140     // Reset "seen" array
141     // const UShort_t kUShortMax = (1 << 16) - 1;
142     // fSeen.Reset(kUShortMax);
143
144     pars   = AliFMDParameters::Instance();
145     map    = pars->GetAltroMap();
146     // stream = new AliAltroRawStream(fReader);
147
148     AliFMDDebug(5, ("Setting %d word headers", 
149                     pars->HasCompleteHeader() ? 8 : 7));
150     stream.SetShortDataHeader(!pars->HasCompleteHeader());
151
152     // Reset variables
153     ddl    = -1;  
154     trate  = 0;   
155     tdet   = 0;   
156     trng   = '\0';
157     tsec   = 0;   
158     tstr   = 0;  
159     tsam   = -1;
160     hwaddr = -1;
161   }
162   do { 
163     Bool_t next = stream.Next();
164     if (!next) { 
165       // if (stream) delete stream;
166       // stream = 0;
167       return kFALSE;
168     }
169     Int_t thisDDL = stream.GetDDLNumber();
170     AliFMDDebug(10, ("RCU @ DDL %d", thisDDL));
171     if (thisDDL != ddl) { 
172       ddl   = thisDDL;
173       fZeroSuppress[ddl] = zs    = stream.GetZeroSupp();
174       fNoiseFactor[ddl]  = fac   = stream.GetNPostsamples();
175       fSampleRate[ddl]   = trate = 0; // stream.GetNPresamples();
176       tdet  = map->DDL2Detector(ddl);
177       AliFMDDebug(10, ("RCU @ DDL %d zero suppression: %s",ddl, zs?"yes":"no"));
178       AliFMDDebug(10, ("RCU @ DDL %d noise factor: %d", ddl,fac));
179       AliFMDDebug(10, ("RCU @ DDL %d sample rate: %d", ddl, trate));
180     }
181     Int_t thisAddr = stream.GetHWAddress();
182     AliFMDDebug(10, ("RCU @ DDL %d, Address 0x%03x", ddl, thisAddr));    
183     if (thisAddr != hwaddr) { 
184       hwaddr = thisAddr;
185       UShort_t  board, chip, channel;
186       map->ChannelAddress(hwaddr, board, chip, channel);
187       map->Channel2StripBase(board, chip, channel, trng, tsec, bstr);
188       AliFMDDebug(10, ("0x%04x/0x%03x maps to FMD%d%c[%2d]-%3d", 
189                       ddl, hwaddr, tdet, trng, tsec, bstr));
190       
191       stripMin = pars->GetMinStrip(tdet, trng, tsec, bstr);
192       stripMax = pars->GetMaxStrip(tdet, trng, tsec, bstr);
193       preSamp  = pars->GetPreSamples(tdet, trng, tsec, bstr);
194       if (trate == 0) 
195         fSampleRate[ddl] = trate = pars->GetSampleRate(tdet, trng, tsec, bstr);
196       AliFMDDebug(10, ("RCU @ DDL %d, Address 0x%03x sample rate: %d", 
197                       ddl, hwaddr, trate));
198       
199       Int_t nChAddrMismatch = stream.GetNChAddrMismatch();
200       Int_t nChLenMismatch  = stream.GetNChLengthMismatch();
201       if (nChAddrMismatch != 0) 
202         AliWarning(Form("Got %d channels with address mis-matches for 0x%03x",
203                         nChAddrMismatch, hwaddr));
204       if (nChLenMismatch != 0) 
205         AliWarning(Form("Got %d channels with length mis-matches for 0x%03x",
206                         nChLenMismatch, hwaddr));
207     }
208     // Get the signal
209     adc = stream.GetSignal();
210
211     // Sanity check - if the total bunch length is less than 1, then
212     // read until we get the next bunch. 
213     Int_t b  = stream.GetTimeLength();
214     if (b < 1) { 
215       AliWarning(Form("Bunch length %0d is less than 0 for "
216                       "DDL %4d address 0x%03x", b, ddl, hwaddr));
217       continue;
218     }
219     // Sanity check - if the current time is less than 0, then read
220     // until we get a new bunch. 
221     Int_t t  = stream.GetTime();
222     if (t < 0) {
223       AliWarning(Form("Time %0d is less than 0 for DDL %4d address 0x%03x", 
224                       t, ddl, hwaddr));
225       continue;
226     }
227
228     Short_t  strOff = 0;
229     UShort_t samp   = 0;
230     map->Timebin2Strip(tsec, t, preSamp, trate, strOff, samp);
231     tstr = bstr + strOff;
232     tsam = samp;
233     AliFMDDebug(20, ("0x%04x/0x%03x/%04d maps to FMD%d%c[%2d,%3d]-%d", 
234                      ddl, hwaddr, t, tdet, trng, tsec, tstr, samp));
235     
236     // Local strip number to channel
237     Short_t l = (tstr > bstr ? tstr - bstr : bstr - tstr);
238     AliFMDDebug(10, ("Checking if strip %d in range [%d,%d]", 
239                      l, stripMin, stripMax));
240     if (l < stripMin || l > stripMax) { 
241       AliFMDDebug(5, ("VA channel %3d (t: %4d) of DDL %4d address 0x%03x "
242                       "is out of range (%3d->%3d)", 
243                       l, t, ddl, hwaddr, stripMin, stripMax));
244       continue;
245     }
246
247
248     det = tdet;
249     rng = trng;
250     sec = tsec;
251     str = tstr;
252     sam = tsam;
253     rat = trate;
254     // adc = stream.GetSignal();
255     
256     break;
257   } while (true);
258   return kTRUE;
259 }
260
261 //____________________________________________________________________
262 Bool_t
263 AliFMDRawReader::NextSignal(UShort_t& det, Char_t&   rng, 
264                             UShort_t& sec, UShort_t& str, 
265                             Short_t&  adc, Bool_t&   zs, 
266                             UShort_t& fac)
267 {
268   
269   do { 
270     UShort_t samp, rate;
271     if (!NextSample(det, rng, sec, str, samp, rate, adc, zs, fac)) 
272       return kFALSE;
273
274     Bool_t take = kFALSE;
275     switch (rate) { 
276     case 1:                      take = kTRUE; break;
277     case 2:  if (samp == 1)      take = kTRUE; break;
278     case 3:  if (samp == 1)      take = kTRUE; break; 
279     case 4:  if (samp == 2)      take = kTRUE; break;
280     default: if (samp == rate-2) take = kTRUE; break;
281     }
282     if (!take) continue;
283     break;
284   } while (true);
285   return kTRUE;
286 }
287
288 //____________________________________________________________________
289 Bool_t
290 AliFMDRawReader::SelectSample(UShort_t samp, UShort_t rate) 
291 {
292   Bool_t take = kFALSE;
293   switch (rate) { 
294   case 1:                      take = kTRUE; break;
295   case 2:  if (samp == 1)      take = kTRUE; break;
296   case 3:  if (samp == 1)      take = kTRUE; break; 
297   case 4:  if (samp == 2)      take = kTRUE; break;
298   default: if (samp == rate-2) take = kTRUE; break;
299   }
300   
301   return take;
302 }
303   
304 #if 1
305 //____________________________________________________________________
306 Bool_t
307 AliFMDRawReader::ReadAdcs(TClonesArray* array) 
308 {
309   
310   fReader->Select("FMD");
311   // Read raw data into the digits array, using AliFMDAltroReader. 
312   if (!array) {
313     AliError("No TClonesArray passed");
314     return kFALSE;
315   }
316   const UShort_t kUShortMax = (1 << 16) - 1;
317   fSeen.Reset(kUShortMax);
318   
319   //  if (!fReader->ReadHeader()) {
320   //    AliError("Couldn't read header");
321   //    return kFALSE;
322   //  }
323   // Get sample rate 
324   AliFMDParameters*     pars = AliFMDParameters::Instance();
325   AliFMDAltroMapping*   map  = pars->GetAltroMap();
326   AliFMDRawStream       input(fReader);
327   AliFMDDebug(5, ("Setting 7 word headers"));
328   input.SetShortDataHeader(!pars->HasCompleteHeader());
329
330   UShort_t stripMin = 0;
331   UShort_t stripMax = 0; // 127;
332   UShort_t preSamp  = 0; // 14+5;
333
334   Int_t  oldddl = -1;
335   UInt_t ddl    = 0;
336   // UInt_t rate   = 0;
337   UInt_t last   = 0;
338   UInt_t hwaddr = 0;
339   // Data array is approx twice the size needed. 
340   UShort_t data[2048];
341   for (size_t i = 0; i < 2048; i++) data[i] = 0; // kUShortMax;
342
343   Bool_t isGood = kTRUE;
344   while (isGood) {
345     isGood = input.ReadChannel(ddl, hwaddr, last, data);
346     // if (!isGood) break;
347     if (ddl >= UINT_MAX /* UInt_t(-1) */) { 
348       AliFMDDebug(5, ("At end of event with %d digits", 
349                       array->GetEntriesFast()));
350       break;
351     }
352     if (UInt_t(oldddl) != ddl) { 
353       fZeroSuppress[ddl] = input.GetZeroSupp();
354       AliFMDDebug(20, ("RCU @ DDL %d zero suppression: %s", 
355                        ddl, (fZeroSuppress[ddl] ? "yes" : "no")));
356
357       // WARNING: We store the noise factor in the 2nd baseline
358       // filters excluded post samples, since we'll never use that
359       // mode. 
360       fNoiseFactor[ddl]  = input.GetNPostsamples();
361       AliFMDDebug(20, ("RCU @ DDL %d noise factor: %d", ddl,fNoiseFactor[ddl]));
362
363       // WARNING: We store the noise factor in the 2nd baseline
364       // filters excluded post samples, since we'll never use that
365       // mode. 
366       fSampleRate[ddl]     = input.GetNPretriggerSamples();
367       AliFMDDebug(20, ("RCU @ DDL %d Sample rate: %d", ddl,fNoiseFactor[ddl]));
368
369       Int_t nChAddrMismatch = input.GetNChAddrMismatch();
370       Int_t nChLenMismatch  = input.GetNChLengthMismatch();
371       if (nChAddrMismatch != 0) 
372         AliWarning(Form("Got %d channels with address mis-matches for 0x%03x",
373                         nChAddrMismatch, hwaddr));
374       if (nChLenMismatch != 0) 
375         AliWarning(Form("Got %d channels with length mis-matches for 0x%03x",
376                         nChLenMismatch, hwaddr));
377       oldddl = ddl;
378     }
379     // AliFMDDebug(5, ("Read channel 0x%x of size %d", hwaddr, last));
380
381     UShort_t det, sec, samp, board, chip, channel;
382     Short_t strbase;
383     Char_t   ring;
384     
385     
386     if (map->DDL2Detector(ddl) < 0) break;
387     det = map->DDL2Detector(ddl);
388     map->ChannelAddress(hwaddr, board, chip, channel);
389     if (!map->Channel2StripBase(board, chip, channel, ring, sec, strbase)) {
390       AliError(Form("Failed to get detector id from DDL %d, "
391                     "hardware address 0x%03x", ddl, hwaddr));
392       continue;
393     }
394     AliFMDDebug(5, ("Board: 0x%02x, Altro: 0x%x, Channel: 0x%x, Length: %4d", 
395                     board, chip, channel, last));
396
397     stripMin = pars->GetMinStrip(det, ring, sec, strbase);
398     stripMax = pars->GetMaxStrip(det, ring, sec, strbase);
399     preSamp  = pars->GetPreSamples(det, ring, sec, strbase);
400     if (fSampleRate[ddl] == 0) 
401       fSampleRate[ddl] = pars->GetSampleRate(det, ring, sec, strbase);
402     
403     // Loop over the `timebins', and make the digits
404     for (size_t i = 0; i < last; i++) {
405       // if (i < preSamp) continue;
406       AliFMDDebug(15, ("0x%04x/0x%03x/%04d %4d", ddl, hwaddr, i, data[i]));
407
408       Short_t  stroff = 0;
409       map->Timebin2Strip(sec, i, preSamp, fSampleRate[ddl], stroff, samp);
410       Short_t  str    = strbase + stroff;
411       
412       AliFMDDebug(10, ("0x%04x/0x%03x/%04d maps to FMD%d%c[%2d,%3d]-%d", 
413                       ddl, hwaddr, i, det, ring, sec, str, samp));
414       if (str < 0) { 
415         AliFMDDebug(8, ("Got presamples at timebin %d", i));
416         data[i] = 0; // Reset cache 
417         continue;
418       }
419       
420       Short_t lstrip = (i - preSamp) / fSampleRate[ddl] + stripMin;
421       
422       AliFMDDebug(15, ("Checking if strip %d (%d) in range [%d,%d]", 
423                       lstrip, str, stripMin, stripMax));
424       if (lstrip < stripMin || lstrip > stripMax) {
425         AliFMDDebug(5, ("FMD%d%c[%02d,%03d]-%d out of range (%3d->%3d)", 
426                         det, ring, sec, samp, str, stripMin, stripMax));
427         data[i] = 0; // Reset cache 
428         continue;
429       }
430       // Possibly do pedestal subtraction of signal 
431       Int_t counts = data[i];
432         
433       
434       // Check the cache of indicies
435       Int_t idx = fSeen(det, ring, sec, str);
436       if (idx == kUShortMax) { 
437         // We haven't seen this strip yet. 
438         fSeen(det, ring, sec, str) = idx = array->GetEntriesFast();
439         AliFMDDebug(7,("making digit for FMD%d%c[%2d,%3d]-%d from timebin %4d", 
440                        det, ring, sec, str, samp, i));
441         new ((*array)[idx]) AliFMDDigit(det, ring, sec, str);
442       }
443       AliFMDBaseDigit* digit = static_cast<AliFMDBaseDigit*>(array->At(idx));
444       AliFMDDebug(10,
445                   ("Setting FMD%d%c[%2d,%3d]-%d from timebin %4d=%4d (%4d)", 
446                    det, ring, sec, str, samp, i, counts, data[i]));
447       digit->SetCount(samp, counts);
448       data[i] = 0; // Reset cache 
449     }
450   }
451   return kTRUE;
452 }
453 #else
454 //____________________________________________________________________
455 Bool_t
456 AliFMDRawReader::ReadAdcs(TClonesArray* array) 
457 {
458   // Read raw data into the digits array, using AliFMDAltroReader. 
459   if (!array) {
460     AliError("No TClonesArray passed");
461     return kFALSE;
462   }
463   //  if (!fReader->ReadHeader()) {
464   //    AliError("Couldn't read header");
465   //    return kFALSE;
466   //  }
467   // Get sample rate 
468   AliFMDParameters* pars = AliFMDParameters::Instance();
469
470   // Select FMD DDL's 
471   fReader->Select("FMD");
472
473   UShort_t stripMin = 0;
474   UShort_t stripMax = 127;
475   UShort_t preSamp  = 0;
476   
477   do {
478     UChar_t* cdata;
479     if (!fReader->ReadNextData(cdata)) break;
480     size_t   nchar = fReader->GetDataSize();
481     UShort_t ddl   = fReader->GetDDLID();
482     UShort_t rate  = 0;
483     AliFMDDebug(1, ("Reading %d bytes (%d 10bit words) from %d", 
484                      nchar, nchar * 8 / 10, ddl));
485     // Make a stream to read from 
486     std::string str((char*)(cdata), nchar);
487     std::istringstream s(str);
488     // Prep the reader class.
489     AliFMDAltroReader r(s);
490     // Data array is approx twice the size needed. 
491     UShort_t data[2048], hwaddr, last;
492     while (r.ReadChannel(hwaddr, last, data) > 0) {
493       AliFMDDebug(5, ("Read channel 0x%x of size %d", hwaddr, last));
494       UShort_t det, sec, str;
495       Char_t   ring;
496       if (!pars->Hardware2Detector(ddl, hwaddr, det, ring, sec, str)) {
497         AliError(Form("Failed to detector id from DDL %d "
498                       "and hardware address 0x%x", ddl, hwaddr));
499         continue;
500       }
501       rate     = pars->GetSampleRate(det, ring, sec, str);
502       stripMin = pars->GetMinStrip(det, ring, sec, str);
503       stripMax = pars->GetMaxStrip(det, ring, sec, str);
504       AliFMDDebug(5, ("DDL 0x%04x, address 0x%03x maps to FMD%d%c[%2d,%3d]", 
505                        ddl, hwaddr, det, ring, sec, str));
506
507       // Loop over the `timebins', and make the digits
508       for (size_t i = 0; i < last; i++) {
509         if (i < preSamp) continue;
510         Int_t    n      = array->GetEntries();
511         UShort_t curStr = str + stripMin + i / rate;
512         if ((curStr-str) > stripMax) {
513           AliError(Form("Current strip is %d but DB says max is %d", 
514                         curStr, stripMax));
515         }
516         AliFMDDebug(5, ("making digit for FMD%d%c[%2d,%3d] from sample %4d", 
517                          det, ring, sec, curStr, i));
518         new ((*array)[n]) AliFMDDigit(det, ring, sec, curStr, data[i], 
519                                       (rate >= 2 ? data[i+1] : 0),
520                                       (rate >= 3 ? data[i+2] : 0));
521         if (rate >= 2) i++;
522         if (rate >= 3) i++;
523         }
524         if (r.IsBof()) break;
525     }
526   } while (true);
527   return kTRUE;
528 }
529
530   
531
532 // This is the old method, for comparison.   It's really ugly, and far
533 // too convoluted. 
534 //____________________________________________________________________
535 void
536 AliFMDRawReader::Exec(Option_t*) 
537 {
538   // Read raw data into the digits array
539   //  if (!fReader->ReadHeader()) {
540   //    Error("ReadAdcs", "Couldn't read header");
541   //    return;
542   //  }
543
544   Int_t n = 0;
545   TClonesArray* array = new TClonesArray("AliFMDDigit");
546   fTree->Branch("FMD", &array);
547
548   // Get sample rate 
549   AliFMDParameters* pars = AliFMDParameters::Instance();
550   fSampleRate = pars->GetSampleRate(0);
551
552   // Use AliAltroRawStream to read the ALTRO format.  No need to
553   // reinvent the wheel :-) 
554   AliFMDRawStream input(fReader, fSampleRate);
555   // Select FMD DDL's 
556   fReader->Select("FMD");
557   
558   Int_t    oldDDL      = -1;
559   Int_t    count       = 0;
560   UShort_t detector    = 1; // Must be one here
561   UShort_t oldDetector = 0;
562   Bool_t   next        = kTRUE;
563
564   // local Cache 
565   TArrayI counts(10);
566   counts.Reset(-1);
567   
568   // Loop over data in file 
569   while (next) {
570     next = input.Next();
571
572     count++; 
573     Int_t ddl = fReader->GetDDLID();
574     AliFMDDebug(10, ("Current DDL is %d", ddl));
575     if (ddl != oldDDL || input.IsNewStrip() || !next) {
576       // Make a new digit, if we have some data (oldDetector == 0,
577       // means that we haven't really read anything yet - that is,
578       // it's the first time we get here). 
579       if (oldDetector > 0) {
580         // Got a new strip. 
581         AliFMDDebug(10, ("Add a new strip: FMD%d%c[%2d,%3d] "
582                           "(current: FMD%d%c[%2d,%3d])", 
583                           oldDetector, input.PrevRing(), 
584                           input.PrevSector() , input.PrevStrip(),
585                           detector , input.Ring(), input.Sector(), 
586                           input.Strip()));
587         new ((*array)[n]) AliFMDDigit(oldDetector, 
588                                       input.PrevRing(), 
589                                       input.PrevSector(), 
590                                       input.PrevStrip(), 
591                                       counts[0], counts[1], counts[2]);
592         n++;
593 #if 0
594         AliFMDDigit* digit = 
595           static_cast<AliFMDDigit*>(fFMD->Digits()->
596                                     UncheckedAt(fFMD->GetNdigits()-1));
597 #endif 
598       }
599         
600       if (!next) { 
601         AliFMDDebug(10, ("Read %d channels for FMD%d", 
602                           count + 1, detector));
603         break;
604       }
605     
606     
607       // If we got a new DDL, it means we have a new detector. 
608       if (ddl != oldDDL) {
609         if (detector != 0) 
610           AliFMDDebug(10, ("Read %d channels for FMD%d", count + 1, detector));
611         // Reset counts, and update the DDL cache 
612         count       = 0;
613         oldDDL      = ddl;
614         // Check that we're processing a FMD detector 
615         Int_t detId = fReader->GetDetectorID();
616         if (detId != (AliDAQ::DetectorID("FMD"))) {
617           AliError(Form("Detector ID %d != %d",
618                         detId, (AliDAQ::DetectorID("FMD"))));
619           break;
620         }
621         // Figure out what detector we're deling with 
622         oldDetector = detector;
623         switch (ddl) {
624         case 0: detector = 1; break;
625         case 1: detector = 2; break;
626         case 2: detector = 3; break;
627         default:
628           AliError(Form("Unknown DDL 0x%x for FMD", ddl));
629           return;
630         }
631         AliFMDDebug(10, ("Reading ADCs for 0x%x  - That is FMD%d",
632                           fReader->GetEquipmentId(), detector));
633       }
634       counts.Reset(-1);
635     }
636     
637     counts[input.Sample()] = input.Count();
638     
639     AliFMDDebug(10, ("ADC of FMD%d%c[%2d,%3d] += %d",
640                       detector, input.Ring(), input.Sector(), 
641                       input.Strip(), input.Count()));
642     oldDetector = detector;
643   }
644   fTree->Fill();
645   return;
646
647 }
648 #endif
649 //____________________________________________________________________
650 Bool_t AliFMDRawReader::ReadSODevent(AliFMDCalibSampleRate* sampleRate, 
651                                      AliFMDCalibStripRange* stripRange, 
652                                      TArrayS &pulseSize, 
653                                      TArrayS &pulseLength) 
654 {
655
656   AliFMDDebug(0, ("Start of SOD/EOD"));
657   
658   UInt_t shift_clk[18];
659   UInt_t sample_clk[18];
660   UInt_t strip_low[18];
661   UInt_t strip_high[18];
662   UInt_t pulse_size[18];
663   UInt_t pulse_length[18];  
664   AliFMDParameters*   param = AliFMDParameters::Instance();
665   AliFMDAltroMapping* map   = param->GetAltroMap();
666   
667   while(fReader->ReadNextData(fData)) {
668     
669     Int_t ddl   = fReader->GetDDLID();
670     Int_t detID = fReader->GetDetectorID();
671     AliFMDDebug(0, (" From reader: DDL number is %d , det ID is %d",ddl,detID));
672     
673     fNbytes = fReader->GetDataSize();
674     
675     ULong_t nWords      = GetNwords();
676     UInt_t  trailerLast = Get32bitWord(nWords);
677     AliFMDDebug(20, (" # Bytes: %d, # Words: %d, Last word: 0x%08x", 
678                      fNbytes, nWords, trailerLast));
679     if ((trailerLast & 0xFFFF0000) != 0xAAAA0000) {
680       AliWarning(Form("Last word 0x%08x does not match RCU II trailer",
681                       trailerLast));
682       return kFALSE;
683     }
684     
685     // 7 last bits is size of trailer
686     ULong_t nTrailerWords = trailerLast & 0x7f; 
687     ULong_t nPayloadWords = Get32bitWord(nWords - nTrailerWords+1);
688
689     AliFMDDebug(20, (" # trailer words: %d, # payload words: %d", 
690                      nTrailerWords, nPayloadWords));
691     
692     
693     for (ULong_t i = 1; i <= nPayloadWords ; i++) {
694       UInt_t payloadWord = Get32bitWord(i);
695       
696       // address is only 24 bit
697       UInt_t address       = (0xffffff & payloadWord);
698       UInt_t type          = ((address >> 21) & 0xf);
699       UInt_t error         = ((address >> 20) & 0x1);
700       UInt_t bcast         = ((address >> 18) & 0x1);
701       UInt_t bc_not_altro  = ((address >> 17) & 0x1);
702       UInt_t board         = ((address >> 12) & 0x1f);
703       UInt_t instruction   = 0;
704       UInt_t chip          = 0;
705       UInt_t channel       = 0;
706       if(bc_not_altro)
707         instruction        = address & 0xfff;
708       else {
709         chip               = ((address >> 9) & 0x7);
710         channel            = ((address >> 5) & 0x5);
711         instruction        = (address & 0x1f);
712       }
713         
714       Bool_t readDataWord = kFALSE;
715       switch(type) {
716       case 0x0: // Fec read
717         readDataWord = kTRUE;  
718       case 0x1: // Fec cmd
719       case 0x2: // Fec write
720         i++;  
721         break;
722       case 0x4: // Loop
723       case 0x5: // Wait
724         break;
725       case 0x6: // End sequence
726       case 0x7: // End Mem
727         i = nPayloadWords + 1;
728         break;
729       default:    
730         break;
731       }
732         
733       //Don't read unless we have a FEC_RD
734       if(!readDataWord)  continue;
735
736       UInt_t dataWord      = Get32bitWord(i);
737       UInt_t data          = (0xFFFFF & dataWord) ;
738       //UInt_t data          = (0xFFFF & dataWord) ;
739         
740       if(error) {
741         AliWarning(Form("error bit detected at Word 0x%06x; "
742                         "error % d, type %d, bc_not_altro %d, "
743                         "bcast %d, board 0x%02x, chip 0x%x, "
744                         "channel 0x%02x, instruction 0x%03x",
745                         address, error, type, bc_not_altro, 
746                         bcast,board,chip,channel,instruction));
747         //process error
748         continue;
749       }
750         
751         
752       switch(instruction) {
753           
754       case 0x01: break;  // First ADC T           
755       case 0x02: break; // I  3.3 V              
756       case 0x03: break; // I  2.5 V altro digital
757       case 0x04: break; // I  2.5 V altro analog 
758       case 0x05: break; // I  2.5 V VA           
759       case 0x06: break; // First ADC T           
760       case 0x07: break; // I  3.3 V              
761       case 0x08: break; // I  2.5 V altro digital
762       case 0x09: break; // I  2.5 V altro analog 
763       case 0x0A: break; // I  2.5 V VA           
764       case 0x2D: break; // Second ADC T           
765       case 0x2E: break; // I  1.5 V VA            
766       case 0x2F: break; // I -2.0 V               
767       case 0x30: break; // I -2.0 V VA            
768       case 0x31: break; //    2.5 V Digital driver
769       case 0x32: break; // Second ADC T           
770       case 0x33: break; // I  1.5 V VA            
771       case 0x34: break; // I -2.0 V               
772       case 0x35: break; // I -2.0 V VA            
773       case 0x36: break; //    2.5 V Digital driver
774       case 0x37: break; // Third ADC T             
775       case 0x38: break; // Temperature sens. 1     
776       case 0x39: break; // Temperature sens. 2     
777       case 0x3A: break; // U  2.5 altro digital (m)
778       case 0x3B: break; // U  2.5 altro analog (m) 
779       case 0x3C: break; // Third ADC T             
780       case 0x3D: break; // Temperature sens. 1     
781       case 0x3E: break; // Temperature sens. 2     
782       case 0x3F: break; // U  2.5 altro digital (m)
783       case 0x40: break; // U  2.5 altro analog (m) 
784       case 0x41: break; // Forth ADC T  
785       case 0x42: break; // U  2.5 VA (m)
786       case 0x43: break; // U  1.5 VA (m)
787       case 0x44: break; // U -2.0 VA (m)
788       case 0x45: break; // U -2.0 (m)   
789       case 0x46: break; // Forth ADC T  
790       case 0x47: break; // U  2.5 VA (m)
791       case 0x48: break; // U  1.5 VA (m)
792       case 0x49: break; // U -2.0 VA (m)
793       case 0x4A: break;  // U -2.0 (m)   
794         // Counters 
795       case 0x0B: break; // L1 trigger CouNTer
796       case 0x0C: break; // L2 trigger CouNTer
797       case 0x0D: break; // Sampling CLK CouNTer
798       case 0x0E: break; // DSTB CouNTer
799         // Test mode 
800       case 0x0F: break; // Test mode word
801       case 0x10: break; // Undersampling ratio.
802         // Configuration and status 
803       case 0x11: break; // Config/Status Register 0
804       case 0x12: break; // Config/Status Register 1
805       case 0x13: break; // Config/Status Register 2
806       case 0x14: break; // Config/Status Register 3
807       case 0x15: break; // Free
808         // Comands:
809       case 0x16: break; // Latch L1, L2, SCLK Counters
810       case 0x17: break; // Clear counters
811       case 0x18: break; // Clear CSR1
812       case 0x19: break; // rstb ALTROs
813       case 0x1A: break; // rstb BC
814       case 0x1B: break; // Start conversion
815       case 0x1C: break; // Scan event length
816       case 0x1D: break; // Read event length
817       case 0x1E: break; // Start test mode
818       case 0x1F: break; // Read acquisition memory
819         // FMD
820       case 0x20: break; // FMDD status
821       case 0x21: break; // L0 counters
822       case 0x22: break; // FMD: Wait to hold
823       case 0x23: break; // FMD: L1 timeout
824       case 0x24: break; // FMD: L2 timeout
825       case 0x25: // FMD: Shift clk 
826         shift_clk[board] = ((data >> 8 ) & 0xFF); 
827         break; 
828       case 0x26: // FMD: Strips 
829         strip_low[board]  = ((data >> 0 ) & 0xFF); 
830         strip_high[board] = ((data >> 8 ) & 0xFF);  
831         break; 
832       case 0x27: // FMD: Cal pulse 
833         pulse_size[board]  =  ((data >> 8 ) & 0xFF);
834         break; 
835       case 0x28: break; // FMD: Shape bias
836       case 0x29: break; // FMD: Shape ref
837       case 0x2A: break; // FMD: Preamp ref
838       case 0x2B: // FMD: Sample clk 
839         sample_clk[board] = ((data >> 8 ) & 0xFF); 
840         break; 
841       case 0x2C: break; // FMD: Commands
842       case 0x4B: // FMD: Cal events 
843         pulse_length[board] = ((data >> 0 ) & 0xFF);
844         break; 
845       default: break;
846         
847       }
848       AliFMDDebug(50, ("instruction 0x%x, dataword 0x%x",
849                        instruction,dataWord));
850     }
851     
852     UShort_t det,sector;
853     Short_t strip;
854     Char_t ring;
855    
856     const UInt_t boards[4] = {0,1,16,17};
857     for(Int_t i=0;i<4;i++) {
858       if(ddl==0 && (i==1 || i==3)) continue;
859
860       UInt_t chip =0, channel=0;
861       det = map->DDL2Detector(ddl);
862       map->Channel2StripBase(boards[i], chip, channel, ring, sector, strip);
863      
864       UInt_t samplerate = 1;
865       if(sample_clk[boards[i]] == 0) {
866         if(ddl == 0) {
867           Int_t sample1 = sample_clk[boards[0]];
868           Int_t sample2 = sample_clk[boards[2]];            
869           if(sample1) sample_clk[boards[i]] = sample1;
870           else sample_clk[boards[i]] = sample2;
871         }
872         else {
873           Int_t sample1 = sample_clk[boards[0]];
874           Int_t sample2 = sample_clk[boards[1]];
875           Int_t sample3 = sample_clk[boards[2]];
876           Int_t sample4 = sample_clk[boards[3]];
877           Int_t agreement = 0;
878           if(sample1 == sample2) agreement++;
879           if(sample1 == sample3) agreement++;
880           if(sample1 == sample4) agreement++;
881           if(sample2 == sample3) agreement++;
882           if(sample2 == sample4) agreement++;
883           if(sample3 == sample4) agreement++;
884             
885           Int_t idx = 0;
886           if(i<3) idx = i+1;
887           else  idx = i-1;
888           if(agreement == 3) {
889             sample_clk[boards[i]] = sample_clk[boards[idx]];
890             shift_clk[boards[i]] = shift_clk[boards[idx]];
891             strip_low[boards[i]] = strip_low[boards[idx]];
892             strip_high[boards[i]] = strip_high[boards[idx]];
893             pulse_length[boards[i]] = pulse_length[boards[idx]];
894             pulse_size[boards[i]] = pulse_size[boards[idx]];
895             AliFMDDebug(0, ("Vote taken for ddl %d, board 0x%x",
896                             ddl,boards[i]));
897           }
898         }
899       } 
900       
901       if(sample_clk[boards[i]])
902         samplerate = shift_clk[boards[i]]/sample_clk[boards[i]];
903       sampleRate->Set(det,ring,sector,0,samplerate);
904       stripRange->Set(det,ring,sector,0,
905                       strip_low[boards[i]],strip_high[boards[i]]);
906       
907       AliFMDDebug(20, ("det %d, ring %c, ",det,ring));
908       pulseLength.AddAt(pulse_length[boards[i]],
909                         GetHalfringIndex(det,ring,boards[i]/16));
910       pulseSize.AddAt(pulse_size[boards[i]],
911                       GetHalfringIndex(det,ring,boards[i]/16));
912       
913       
914       
915       AliFMDDebug(20, (": Board: 0x%02x\n"
916                        "\tstrip_low  %3d, strip_high   %3d\n"
917                        "\tshift_clk  %3d, sample_clk   %3d\n"
918                        "\tpulse_size %3d, pulse_length %3d",
919                        boards[i], 
920                        strip_low[boards[i]], strip_high[boards[i]],
921                        shift_clk[boards[i]], sample_clk[boards[i]],
922                        pulse_size[boards[i]],pulse_length[boards[i]]));
923     }
924     
925   }
926   
927   AliFMDParameters::Instance()->SetSampleRate(sampleRate);
928   AliFMDParameters::Instance()->SetStripRange(stripRange);
929   
930   AliFMDDebug(0, ("End of SOD/EOD"));
931   
932   return kTRUE;
933 }
934 //____________________________________________________________________
935
936 UInt_t AliFMDRawReader::Get32bitWord(Int_t idx)
937 {
938   // This method returns the 32 bit word at a given
939   // position inside the raw data payload.
940   // The 'index' points to the beginning of the next word.
941   // The method is supposed to be endian (platform)
942   // independent.
943   if (!fData) {
944     AliFatal("Raw data paylod buffer is not yet initialized !");
945   }
946
947   Int_t index = 4*idx;
948   
949   if (index < 4) {
950     //  fRawReader->AddFatalErrorLog(k32bitWordReadErr,Form("pos = %d",index));
951     //    PrintDebug();
952     AliWarning(Form("Invalid raw data payload position (%d) !",index));
953   }
954
955   UInt_t word = 0;
956    
957   word  = fData[--index] << 24;
958   word |= fData[--index] << 16;
959   word |= fData[--index] << 8;
960   word |= fData[--index] << 0 ;
961
962   return word;
963 }
964 //_____________________________________________________________________ 
965 Int_t AliFMDRawReader::GetHalfringIndex(UShort_t det, Char_t ring, 
966                                         UShort_t board) {
967
968   UShort_t iring  =  (ring == 'I' ? 1 : 0);
969   
970   Int_t index = (((det-1) << 2) | (iring << 1) | (board << 0));
971   
972   return index-2;
973   
974 }
975 //____________________________________________________________________
976 // 
977 // EOF
978 //