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