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