]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDRawReader.cxx
Remove print
[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 "AliFMDUShortMap.h"
62 // #include "AliFMDAltroIO.h"   // ALIFMDALTROIO_H 
63 #include "AliAltroRawStreamV3.h"
64 #include <TArrayS.h>            // ROOT_TArrayS
65 #include <TTree.h>              // ROOT_TTree
66 #include <TClonesArray.h>       // ROOT_TClonesArray
67 #include <TString.h>
68 #include <iostream>
69 #include <climits>
70 // #include <iomanip>
71
72 //____________________________________________________________________
73 ClassImp(AliFMDRawReader)
74 #if 0
75   ; // This is here to keep Emacs for indenting the next line
76 #endif
77
78 //____________________________________________________________________
79 AliFMDRawReader::AliFMDRawReader(AliRawReader* reader, TTree* tree) 
80   : TTask("FMDRawReader", "Reader of Raw ADC values from the FMD"),
81     fTree(tree),
82     fReader(reader), 
83     // fSampleRate(1),
84     fData(0),
85     fNbytes(0), 
86     fMinStrip(0),
87     fMaxStrip(127), 
88     fPreSamp(14+5),
89     fSeen(0)
90 {
91   // Default CTOR
92   for (Int_t i = 0; i < 3; i++) { 
93     fSampleRate[i]   = 0;
94     fZeroSuppress[i] = kFALSE;
95     fNoiseFactor[i]  = 1;
96     fL1Phase[i]      = 0;
97     fNErrors[i]      = 0;
98   }
99 }
100
101 //____________________________________________________________________
102 void
103 AliFMDRawReader::Exec(Option_t*) 
104 {
105   // Read the data 
106   TClonesArray* array = new TClonesArray("AliFMDDigit");
107   if (!fTree) {
108     AliError("No tree");
109     return;
110   }
111   fTree->Branch("FMD", &array);
112   
113   
114   ReadAdcs(array);
115   Int_t nWrite = fTree->Fill();
116   AliDebug(1,Form("Got a grand total of %d digits, wrote %d bytes to tree", 
117                    array->GetEntriesFast(), nWrite));
118   delete array;
119 }
120
121 //____________________________________________________________________
122 Int_t
123 AliFMDRawReader::NewDDL(AliAltroRawStreamV3& input, UShort_t& det)
124 {
125   // Process a new DDL.  Sets the internal data members fZeroSuppress, 
126   // fSampleRate, and fNoiseFactor based on information in the RCU trailer. 
127   // 
128   // Parameters:
129   //   input Input stream
130   //   det   On return, the detector number
131   // 
132   // Return value:
133   //   negative value in case of problems, the DDL number otherwise
134
135   // Get the DDL number
136   UInt_t ddl = input.GetDDLNumber();
137   AliDebug(2,Form("DDL number %d", ddl));
138
139   // Note, previously, the ALTROCFG1 register was interpreted as 
140   // 
141   // Bits    Value    Description
142   //   0- 3      0/1   1st Baseline filter, mode 
143   //   4- 5   Over-1   2nd baseline filter, # of pre-samples
144   //   6- 9   factor   2nd baseline filter, # of post-samples 
145   //  10-          0   2nd baseline filter, enable
146   //  11-12       00   Zero suppression, glitch filter mode
147   //  13-15      001   Zero suppression, # of post samples
148   //  16-17       01   Zero suppression, # of pre  samples
149   //  18         0/1   Zero suppression, enable
150   //
151   // The interpretation used in AliAltroRawStreamerV3 - which
152   // corresponds directly to ALTRO DPCFG register - is 
153   //
154   // Bits    Value  Description
155   //   0- 3    0/1   1st Baseline filter, mode 
156   //   4         0   Polarity (if '1', then "1's inverse")
157   //   5- 6     01   Zero suppression, # of pre samples
158   //   7-10   0001   Zero suppression, # of post samples
159   //  11         0   2nd baseline filter, enable
160   //  12-13     00   Zero suppression, glitch filter mode
161   //  14-16 factor   2nd baseline filter, # of post-samples
162   //  17-18     01   2nd baseline filter, # of pre-samples 
163   //  19       0/1   Zero suppression, enable
164   //
165   //  Writing 'x' for variable values, that means we have the
166   //  following patterns for the 2 cases 
167   //
168   //    bit #  20   16   12    8    4    0
169   //     old    |0x01|0010|00xx|xxxx|xxxx|
170   //     new    |x01x|xx00|0000|1010|xxxx|
171   //
172   //  That means that we can check if bits 10-13 are '1000' or
173   //  '0000', which will tell us if the value was written with the
174   //  new or the old interpretation.    That is, we can check that 
175   //
176   //    if (((altrocfg1 >> 10) & 0x8) == 0x8) { 
177   //      // old interpretation 
178   //    }
179   //    else { 
180   //      // New interpretation 
181   //    }
182   //
183   // That means, that we should never 
184   //
185   //  - change the # of zero suppression post samples 
186   //  - Turn on 2nd baseline correction 
187   //  - Change the zero-suppression glitch filter mode
188   //
189   // This change as introduced in version 1.2 of Rcu++
190   //
191   UInt_t cfg1 = input.GetAltroCFG1();
192   if (((cfg1 >> 10) & 0x8) == 0x8) {
193     UInt_t cfg2 = input.GetAltroCFG2();
194     AliDebug(3,Form("We have data from older MiniConf 0x%x cfg2=0x%08x", 
195                     ((cfg1 >> 10) & 0x8), cfg2));
196     fZeroSuppress[ddl] = (cfg1 >>  0) & 0x1;
197     fNoiseFactor[ddl]  = (cfg1 >>  6) & 0xF;
198     fSampleRate[ddl]   = (cfg2 >> 20) & 0xF;
199   }
200   else {
201     AliDebug(3,Form("We have data from newer MiniConf 0x%x", 
202                     ((cfg1 >> 10) & 0x8)));
203     fZeroSuppress[ddl] = input.GetZeroSupp();
204     // WARNING: We store the noise factor in the 2nd baseline
205     // filters excluded post samples, since we'll never use that
206     // mode. 
207     fNoiseFactor[ddl]  = input.GetNPostsamples();
208     // WARNING: We store the sample rate in the number of pre-trigger
209     // samples, since we'll never use that mode.
210     fSampleRate[ddl]     = input.GetNPretriggerSamples();
211     // 
212   }
213   AliDebug(10,Form("Phase of DDL=%d is %g (%d)", ddl, input.GetL1Phase(),
214                    input.GetAltroCFG2() & 0x1F));
215   fL1Phase[ddl] = input.GetAltroCFG2() & 0x1F; // input.GetL1Phase();
216   AliDebug(3,Form("RCU @ DDL %d zero suppression: %s", 
217                    ddl, (fZeroSuppress[ddl] ? "yes" : "no")));
218   AliDebug(3,Form("RCU @ DDL %d noise factor: %d", ddl,fNoiseFactor[ddl]));    
219   AliDebug(3,Form("RCU @ DDL %d sample rate: %d", ddl,fSampleRate[ddl]));
220
221
222   // Get Errors seen 
223   Int_t nChAddrMismatch = input.GetNChAddrMismatch();
224   Int_t nChLenMismatch  = input.GetNChLengthMismatch();
225   if (nChAddrMismatch != 0) 
226     AliWarning(Form("Got %d channels with address mis-matches for 0x%03x",
227                     nChAddrMismatch, ddl));
228   if (nChLenMismatch != 0) 
229     AliWarning(Form("Got %d channels with length mis-matches for 0x%03x",
230                     nChLenMismatch, ddl));
231
232   // Map DDL number to the detector number 
233   AliFMDParameters*    pars   = AliFMDParameters::Instance();
234   AliFMDAltroMapping*  map    = pars->GetAltroMap();
235   if (map->DDL2Detector(ddl) < 0) return -1;
236   det = map->DDL2Detector(ddl);
237
238   if (AliLog::GetDebugLevel("FMD", 0) > 5) 
239     input.PrintRCUTrailer();
240   return ddl;
241 }
242
243 //____________________________________________________________________
244 Int_t
245 AliFMDRawReader::NewChannel(const AliAltroRawStreamV3& input,  UShort_t det,
246                             Char_t&  ring, UShort_t& sec, Short_t& strbase)
247 {
248   // Processs a new channel.  Sets the internal data members
249   // fMinStrip, fMaxStrip, and fPreSamp. 
250   //
251   // Parameter:
252   //   input   Input stream
253   //   ring    On return, the ring identifier 
254   //   sec     On return, the sector number
255   //   strbase On return, the strip base
256   // 
257   // Return value
258   //   negative value in case of problems, hardware address otherwise
259
260   // Get the hardware address, and map that to detector coordinates 
261   UShort_t board, chip, channel;
262   Int_t    ddl    = input.GetDDLNumber();
263   Int_t    hwaddr = input.GetHWAddress();
264   if (input.IsChannelBad()) { 
265     const char* msg = Form("Ignoring channel %03d/0x%03x with errors", 
266                            ddl, hwaddr); 
267     AliWarning(msg); 
268     if (AliDebugLevel() > 10) input.HexDumpChannel();
269     fReader->AddMinorErrorLog(AliAltroRawStreamV3::kAltroPayloadErr,msg);
270     fNErrors[ddl] += 1;
271     return 0xFFFF;
272   }
273   
274   AliFMDParameters*    pars   = AliFMDParameters::Instance();
275   AliFMDAltroMapping*  map    = pars->GetAltroMap();
276   // Map to hardware stuff 
277   map->ChannelAddress(hwaddr, board, chip, channel);
278   // Then try to map to detector address
279   if (!map->Channel2StripBase(board, chip, channel, ring, sec, strbase)) {
280     AliError(Form("Failed to get detector id from DDL %d, "
281                   "hardware address 0x%03x", ddl, hwaddr));
282     return -1;
283   }
284   AliDebug(4,Form("Board: 0x%02x, Altro: 0x%x, Channel: 0x%x", 
285                   board, chip, channel));
286
287   // Get the 'conditions'
288   fMinStrip = pars->GetMinStrip(det, ring, sec, strbase);
289   fMaxStrip = pars->GetMaxStrip(det, ring, sec, strbase);
290   fPreSamp  = pars->GetPreSamples(det, ring, sec, strbase);
291   if (fSampleRate[ddl] == 0) 
292     fSampleRate[ddl] = pars->GetSampleRate(det, ring, sec, strbase);
293
294   return hwaddr;
295 }
296
297 //____________________________________________________________________
298 Bool_t
299 AliFMDRawReader::NewBunch(const AliAltroRawStreamV3& input, 
300                           UShort_t&  start, UShort_t& length)
301 {
302   // 
303   // Do some checks on the bunch data 
304   // 
305   Int_t    ddl      = input.GetDDLNumber();
306   Int_t    hwaddr   = input.GetHWAddress();  
307   UShort_t nSamples = input.GetNSamplesPerCh() + fPreSamp;
308   UShort_t tstart   = input.GetStartTimeBin();
309   length            = input.GetBunchLength();
310
311   if (tstart >= nSamples) {
312     const char* msg = Form("Bunch in %03d/0x%03x has an start time greater "
313                            "than number of samples: 0x%x >= 0x%x", 
314                            ddl, hwaddr, tstart, nSamples);
315     AliWarning(msg);
316     if (AliDebugLevel() > 10) input.HexDumpChannel();
317     fReader->AddMinorErrorLog(AliAltroRawStreamV3::kAltroPayloadErr,msg);
318     fNErrors[ddl]++;
319     return false;
320   }
321   if ((int(tstart) - length + 1) < 0) { 
322     const char* msg = Form("Bunch in %03d/0x%03x has an invalid length and "
323                            "start time: 0x%x,0x%x (%d-%d+1=%d<0)", 
324                            ddl, hwaddr, length, tstart, tstart, length, 
325                            int(tstart)-length+1);
326     AliWarning(msg);
327     if (AliDebugLevel() > 10) input.HexDumpChannel();
328     fReader->AddMinorErrorLog(AliAltroRawStreamV3::kAltroPayloadErr,msg);
329     fNErrors[ddl]++;                            
330     return false;
331   }
332   if (tstart >= start) { 
333     const char* msg = Form("Bunch in %03d/0x%03x has early start time: "
334                            "0x%x >= 0x%x", ddl, hwaddr, tstart, start);
335     AliWarning(msg);
336     if (AliDebugLevel() > 10) input.HexDumpChannel();
337     fReader->AddMinorErrorLog(AliAltroRawStreamV3::kAltroPayloadErr,msg);
338     fNErrors[ddl]++;
339     return false;
340   }
341   start = tstart;
342   return true;
343 }
344
345 //____________________________________________________________________
346 Int_t
347 AliFMDRawReader::NewSample(const AliAltroRawStreamV3& input, 
348                            Int_t i, UShort_t t, UShort_t sec,
349                            UShort_t  strbase, Short_t&  str, UShort_t& samp)
350 {
351   // Process a new timebin
352   // 
353   // Parameters:
354   //   input   Input stream
355   //   i       Index into bunch data
356   //   t       Time
357   //   strbase Base of strip numbers for this channel
358   //   str     On return, the strip number
359   //   samp    On return, the sample number
360   // 
361   // Return value
362   //   negative value in case of problems, ADC value otherwise
363   if (t < fPreSamp) return -1;
364
365   Int_t           ddl  = input.GetDDLNumber();
366   Int_t           hwa  = input.GetHWAddress();
367   const UShort_t* data = input.GetSignals();
368   Short_t         adc  = data[i];
369   AliDebug(10,Form("0x%04x/0x%03x/%04d %4d", ddl, hwa, t, adc));
370
371   AliFMDParameters*    pars   = AliFMDParameters::Instance();
372   AliFMDAltroMapping*  map    = pars->GetAltroMap();
373
374   samp            = 0;
375   Short_t  stroff = 0;
376   map->Timebin2Strip(sec, t, fPreSamp, fSampleRate[ddl], stroff, samp);
377   str             = strbase + stroff;
378       
379   AliDebug(20,Form("0x%04x/0x%03x/%04d=%4d maps to strip %3d sample %d " 
380                    "(pre: %d, min: %d, max: %d, rate: %d)",
381                   ddl, hwa, t, adc, str, samp, fPreSamp, 
382                   fMinStrip, fMaxStrip, fSampleRate[ddl]));
383   if (str < 0) { 
384     AliDebug(10,Form("Got presamples at timebin %d", i));
385     return -1;
386   }
387           
388   // VA1 Local strip number 
389   Short_t lstrip = (t - fPreSamp) / fSampleRate[ddl] + fMinStrip;
390       
391   AliDebug(15,Form("Checking if strip %d (%d) in range [%d,%d]", 
392                    lstrip, str, fMinStrip, fMaxStrip));
393   if (lstrip < fMinStrip || lstrip > fMaxStrip) {
394     AliDebug(10,Form("Strip %03d-%d (%d,%d) from t=%d out of range (%3d->%3d)", 
395                     str, samp, lstrip, stroff, t, fMinStrip, fMaxStrip));
396     adc = -1;
397   }
398   // Possibly do pedestal subtraction of signal 
399   if (adc > 1023) 
400     AliWarning(Form("ADC value out of range: %4d", adc));
401   return adc;
402 }
403
404 //____________________________________________________________________
405 Int_t
406 AliFMDRawReader::NextSample(UShort_t& det, Char_t&   rng, UShort_t& sec, 
407                             UShort_t& str, UShort_t& sam, UShort_t& rat, 
408                             Short_t&  adc, Bool_t&   zs,  UShort_t& fac)
409 {
410   // Scan current event for next signal.   It returns kFALSE when
411   // there's no more data in the event. 
412   // 
413   // Note, that this member function is in principle very fast, but
414   // contains less error checking.  In particular, channels that have
415   // bad bunches cannot be checked here.  Seeing a bad bunch will only
416   // skip the remainder of the channel and not reset the already read
417   // digits.   This is potentially dangerous. 
418   //
419   // Parameters: 
420   //    det         On return, contain the detector number 
421   //    rng         On return, contain the ring identifier 
422   //    sec         On return, contain the sector number 
423   //    str         On return, contain the strip number 
424   //    sam         On return, contain the sample number 
425   //    rat         On return, contain the sample rate 
426   //    adc         On return, contain the ADC counts 
427   //    zs          On return, contain the zero-supp. flag 
428   //    fac         On return, contain the zero-supp. noise factor 
429   // 
430   // Return values: 
431   //    0    No more data 
432   //    -1   Read sample belongs to a bad bunch 
433   //    >0   Good status - contains bit mask of values 
434   //       Bit 1    New DDL
435   //       Bit 2    New Channel
436   //       Bit 3    New Bunch
437   //       Bit 4    New Sample
438   static AliAltroRawStreamV3 stream(fReader); //    = 0;
439   static Int_t               ddl      = -1;
440   static UShort_t            tdet     = 0;
441   static Char_t              trng     = '\0';
442   static UShort_t            tsec     = 0;
443   static Short_t             tstr     = 0;   
444   static Short_t             bstr     = -1;
445   static UShort_t            tsam     = 0;   
446   static UInt_t              trate    = 0;
447   static Int_t               hwaddr   = -1;
448   static UShort_t            start    = 0;
449   static UShort_t            length   = 0;
450   static Short_t             t        = -1;
451   static Int_t               i        = 0; 
452   // First entry!
453   if (stream.GetDDLNumber() < 0) { 
454     fReader->Reset();
455     fReader->Select("FMD");
456     stream.Reset();
457     stream.SelectRawData("FMD");
458     stream.SetCheckAltroPayload(false);
459     for (Int_t j = 0; j < kNDDL; j++) fNErrors[j] = 0;
460
461     // Reset variables
462     ddl    = -1;  
463     trate  = 0;   
464     tdet   = 0;   
465     trng   = '\0';
466     tsec   = 0;   
467     tstr   = 0;  
468     tsam   = -1;
469     hwaddr = -1;
470   }
471
472   UShort_t ret = 0;
473   do { 
474     AliDebug(15,Form("t=%4d, start=%4d, length=%4d", t, start, length));
475     if (t < start - length + 1) { 
476       AliDebug(10,Form("Time t=%d < start-length+1=%d-%d+1 (%3d/0x%03x)", 
477                        t, start, length, ddl, hwaddr));
478       if (hwaddr > 0xFFF || 
479           hwaddr < 0 || 
480           !stream.NextBunch()) { 
481         if (AliDebugLevel() >= 10 && hwaddr > 0xFFF) {
482           AliDebug(10,"Last channel read was marked bad");
483         }
484         if (AliDebugLevel() >= 10 && hwaddr < 0) {
485           AliDebug(10,"No more channels");
486         }
487         AliDebug(10,"No next bunch, or first entry");
488         if (ddl < 0 || !stream.NextChannel()) { 
489           if (AliDebugLevel() >= 10 && ddl < 0) { 
490             AliDebug(10,"No DDL");
491           }
492           AliDebug(10,"No next channel, or first entry");
493           if (!stream.NextDDL()) {
494             AliDebug(10,"No more DDLs");
495             stream.Reset();
496             return 0;
497           }
498           ddl = NewDDL(stream, tdet);
499           AliDebug(5,Form("New DDL: %d (%d)", ddl, tdet));
500           ret |= 0x1;
501           continue;
502         }
503         hwaddr = NewChannel(stream, tdet, trng, tsec, bstr);
504         if (hwaddr > 0xFFF) fNErrors[ddl] += 1;
505         AliDebug(5,Form("New Channel: %3d/0x%03x", ddl, hwaddr));
506         start  = 1024;
507         ret |= 0x2;
508         continue;
509       }
510       if (!NewBunch(stream, start, length)) { 
511         // AliWarning(Form("Bad bunch in %3d/0x%03x read - "
512         //                 "should progress to next channel "
513         //                 "(t=%4d,start=%4d,length=%4d)", 
514         //                 ddl, hwaddr, t,start, length));
515         hwaddr = 0xFFFF; // Bad channel
516         return -1;
517       }
518       AliDebug(5, Form("New bunch in  %3d/0x%03x: start=0x%03x, length=%4d", 
519                        ddl, hwaddr, start, length));
520       ret |= 0x4;
521       t      = start;
522       i      = 0;
523       AliDebug(10,Form("Got new bunch FMD%d%c[%2d], bunch @ %d, length=%d", 
524                        tdet, trng, tsec, start, length));
525     }
526     Int_t tadc = NewSample(stream, i, t, tsec, bstr, tstr, tsam);
527     AliDebug(10,Form("New sample FMD%d%c[%2d,%3d]-%d = 0x%03x", 
528                  tdet, trng, tsec, tstr, tsam, tadc));
529     ret |= 0x8;
530     if (tadc >= 0) { 
531       det = tdet;
532       rng = trng;
533       sec = tsec;
534       str = tstr;
535       sam = tsam;
536       adc = tadc;
537       rat = fSampleRate[ddl];
538       zs  = fZeroSuppress[ddl];
539       fac = fNoiseFactor[ddl];
540       t--;
541       i++;
542       AliDebug(10,Form("Returning FMD%d%c[%2d,%3d]-%d = 0x%03x (%d,%d,%d)",
543                    det, rng, sec, str, sam, adc, rat, zs, fac));
544       break;
545     }
546     t--;
547     i++;
548   } while (true);
549   AliDebug(5,Form("Returning 0x%02x", ret));
550   return ret;
551 }
552
553
554 //____________________________________________________________________
555 Int_t
556 AliFMDRawReader::NextSignal(UShort_t& det, Char_t&   rng, 
557                             UShort_t& sec, UShort_t& str, 
558                             Short_t&  adc, Bool_t&   zs, 
559                             UShort_t& fac)
560 {
561   // 
562   // Get the next signal
563   // 
564   // Parameters:
565   //    det  On return, the detector
566   //    rng  On return, the ring
567   //    sec  On return, the sector
568   //    str  On return, the strip
569   //    adc  On return, the ADC value
570   //    zs   On return, whether zero-supp. is enabled
571   //    fac  On return, the usd noise factor
572   // 
573   // Return:
574   //    true if valid data is returned
575   //
576   Int_t ret = 0;
577   do { 
578     UShort_t samp, rate;
579     if ((ret = NextSample(det, rng, sec, str, samp, rate, adc, zs, fac)) <= 0)
580       return ret;
581
582     Bool_t take = SelectSample(samp, rate);
583     if (!take) continue;
584     break;
585   } while (true);
586   return ret;
587 }
588
589 //____________________________________________________________________
590 Bool_t
591 AliFMDRawReader::SelectSample(UShort_t samp, UShort_t rate) 
592 {
593   // Check if the passed sample is the one we need
594   Bool_t take = kFALSE;
595   switch (rate) { 
596   case 1:                      take = kTRUE; break;
597   case 2:  if (samp == 1)      take = kTRUE; break;
598   case 3:  if (samp == 1)      take = kTRUE; break; 
599   case 4:  if (samp == 2)      take = kTRUE; break;
600   default: if (samp == rate-2) take = kTRUE; break;
601   }
602   
603   return take;
604 }
605   
606 //____________________________________________________________________
607 Bool_t
608 AliFMDRawReader::ReadAdcs(TClonesArray* array) 
609 {
610   // Read ADC values from raw input into passed TClonesArray of AliFMDDigit
611   // objects. 
612   AliDebug(3,Form("Reading ADC values into a TClonesArray"));
613
614   // Read raw data into the digits array, using AliFMDAltroReader. 
615   if (!array) {
616     AliError("No TClonesArray passed");
617     return kFALSE;
618   }
619   const UShort_t kUShortMax = (1 << 16) - 1;
620   fSeen.Reset(kUShortMax);
621   for (Int_t ddl = 0; ddl < kNDDL; ddl++) fNErrors[ddl] = 0;
622
623   AliAltroRawStreamV3  input(fReader);
624   input.Reset();
625   input.SetCheckAltroPayload(false);
626   input.SelectRawData("FMD");
627   
628   // Loop over input RORCs
629   while (input.NextDDL()) { 
630     UShort_t det = 0;
631     Int_t    ddl = NewDDL(input, det);
632     if (ddl < 0) break;
633     fNErrors[ddl] = 0;
634
635     while (input.NextChannel()) { 
636       // Get the hardware address, and map that to detector coordinates 
637       Char_t   ring;
638       UShort_t sec;
639       Short_t  strbase;
640       Int_t    hwaddr   = NewChannel(input, det, ring, sec, strbase);
641       if (hwaddr < 0) break;
642       if (hwaddr > 0xFFF) continue;  
643
644       UShort_t start    = 0x3FF;
645       Bool_t   errors   = false;
646       Int_t    first    = -1;
647       Int_t    last     = -1;
648       // Loop over bunches 
649       while (input.NextBunch()) { 
650         // Get Lenght of bunch, and pointer to the data 
651         const UShort_t* data   = input.GetSignals();
652         UShort_t        length;
653         if (!NewBunch(input, start, length)) {
654           errors = true;
655           break;
656         }
657
658       
659         // Loop over the data and store it. 
660         for (Int_t i = 0; i < length; i++) { 
661           // Time 
662           Short_t  str;
663           UShort_t samp;
664           Int_t    t    = start - i;
665           Int_t    adc  = NewSample(input, i, t, sec, strbase, str, samp);
666           if (adc < 0) continue;
667           UShort_t counts = adc;
668       
669           AliDebug(10, Form("FMD%d%c[%02d,%03d]-%d: %4d", 
670                             det, ring, sec, str, samp, counts));
671           // Check the cache of indicies
672           Int_t idx = fSeen(det, ring, sec, str);
673           AliFMDDigit* digit = 0;
674           if (idx == kUShortMax) { 
675             // We haven't seen this strip yet. 
676             fSeen(det, ring, sec, str) = idx = array->GetEntriesFast();
677             AliDebug(7,Form("making digit @ %5d for FMD%d%c[%2d,%3d]-%d "
678                             "from %3d/0x%03x/%4d", 
679                             idx, det, ring, sec, str, samp, ddl, hwaddr, t));
680             digit = new ((*array)[idx]) AliFMDDigit(det, ring, sec, str);
681             digit->SetDefaultCounts(fSampleRate[ddl]);
682           }
683           else {
684             digit = static_cast<AliFMDDigit*>(array->At(idx));
685           }
686           if (first < 0) first = idx;
687           last = idx;
688           AliDebug(10, Form("Setting FMD%d%c[%2d,%3d]-%d from timebin "
689                             "%4d=%4d (%4d)", det, ring, sec, str, samp, t, 
690                             counts, data[i]));
691           digit->SetCount(samp, counts);
692         } // for (i)
693       } // while (bunch)
694       if (errors) { 
695         AliWarning(Form("Channel %3d/0x%03x contain errors, "
696                         "resetting index %d to %d", ddl, hwaddr, first, last));
697         if (first >= 0) {
698           for (Int_t i = first; i <= last; i++) { 
699             AliFMDDigit* digit = static_cast<AliFMDDigit*>(array->At(i));
700             for (Int_t j = 0; j < fSampleRate[ddl]; j++) {
701               AliDebug(10,Form("Resetting strip %s=%d",
702                                digit->GetName(),digit->Counts()));
703               digit->SetCount(j, kBadSignal);
704             }
705           }
706         }
707       }
708       // if (errors && (AliDebugLevel() > 0)) input.HexDumpChannel();
709     } // while (channel)
710   } // while (ddl)
711   return kTRUE;
712 }
713 //____________________________________________________________________
714 Bool_t
715 AliFMDRawReader::ReadAdcs(AliFMDUShortMap& map) 
716 {
717   // Read ADC values from raw input into passed TClonesArray of AliFMDDigit
718   // objects. 
719   AliDebug(3,Form("Reading ADC values into a map"));
720
721   // const UShort_t kUShortMax = (1 << 16) - 1;
722   for (Int_t ddl = 0; ddl < kNDDL; ddl++) fNErrors[ddl] = 0;
723
724   AliAltroRawStreamV3  input(fReader);
725   input.Reset();
726   input.SetCheckAltroPayload(false);
727   input.SelectRawData("FMD");
728   
729   // Loop over input RORCs
730   while (input.NextDDL()) { 
731     UShort_t det = 0;
732     Int_t    ddl = NewDDL(input, det);
733     if (ddl < 0) break;
734     fNErrors[ddl] = 0;
735
736     while (input.NextChannel()) { 
737       // Get the hardware address, and map that to detector coordinates 
738       Char_t   ring;
739       UShort_t sec;
740       Short_t  strbase;
741       Int_t    hwaddr   = NewChannel(input, det, ring, sec, strbase);
742       if (hwaddr < 0) break;
743       if (hwaddr > 0xFFF) continue;  
744
745       UShort_t start    = 0x3FF;
746       Bool_t   errors   = false;
747       Int_t    first    = -1;
748       Int_t    last     = -1;
749       // Loop over bunches 
750       while (input.NextBunch()) { 
751         // Get Lenght of bunch, and pointer to the data 
752         // const UShort_t* data   = input.GetSignals();
753         UShort_t        length;
754         if (!NewBunch(input, start, length)) {
755           errors = true;
756           break;
757         }
758
759       
760         // Loop over the data and store it. 
761         for (Int_t i = 0; i < length; i++) { 
762           // Time 
763           Short_t  str;
764           UShort_t samp;
765           Int_t    t    = start - i;
766           Int_t    adc  = NewSample(input, i, t, sec, strbase, str, samp);
767           if (adc < 0) continue;
768           UShort_t counts = adc;
769       
770           AliDebug(10, Form("FMD%d%c[%02d,%03d]-%d: %4d", 
771                             det, ring, sec, str, samp, counts));
772           if (SelectSample(samp, fSampleRate[ddl]))
773             map(det,ring,sec,str) = counts; 
774           if (first < 0) first = str;
775           last = str;
776         } // for (i)
777       } // while (bunch)
778       if (errors) { 
779         AliWarning(Form("Channel %3d/0x%03x contain errors, "
780                         "resetting strips %d to %d", ddl, hwaddr, first, last));
781         if (first >= 0) {
782           Int_t ds = first <= last ? 1 : -1;
783           for (Int_t i = first; i != last+ds; i += ds) { 
784             AliDebug(10, Form("Resetting strip FMD%d%c[%02d,%03d]=%d",
785                               det,ring,sec,i,map(det,ring,sec,i)));
786             map(det,ring,sec,i) = kBadSignal;
787           }
788         }
789       }
790     } // while (channel)
791   } // while (ddl)
792   return kTRUE;
793 }
794
795 //____________________________________________________________________
796 Bool_t AliFMDRawReader::ReadSODevent(AliFMDCalibSampleRate* sampleRate, 
797                                      AliFMDCalibStripRange* stripRange, 
798                                      TArrayS &pulseSize, 
799                                      TArrayS &pulseLength, 
800                                      Bool_t* detectors) 
801 {
802   // 
803   // Read SOD event into passed objects.
804   // 
805   // Parameters:
806   //    samplerate   The sample rate object to fill
807   //    striprange   The strip range object to fill
808   //    pulseSize    The pulse size object to fill
809   //    pulseLength  The pulse length (in events) object to fill
810   // 
811   // Return:
812   //    @c true on success
813   //  
814   AliDebug(0,Form("Start of SOD/EOD"));
815   
816   UInt_t shift_clk[18];
817   UInt_t sample_clk[18];
818   UInt_t strip_low[18];
819   UInt_t strip_high[18];
820   UInt_t pulse_size[18];
821   UInt_t pulse_length[18];  
822   for (size_t i = 0; i < 18; i++) { 
823     shift_clk[i]    = 0;
824     sample_clk[i]   = 0;
825     strip_low[i]    = 0;
826     strip_high[i]   = 0;
827     pulse_size[i]   = 0;
828     pulse_length[i] = 0;
829   }
830   AliFMDParameters*   param = AliFMDParameters::Instance();
831   AliFMDAltroMapping* map   = param->GetAltroMap();
832   
833   AliAltroRawStreamV3  streamer(fReader);
834   streamer.Reset();
835   streamer.SelectRawData("FMD");
836   //fReader->GetDDLID();
837   while (streamer.NextDDL()) {
838     Int_t ddl   = streamer.GetDDLNumber();
839     Int_t detID = fReader->GetDetectorID();
840     if (detectors) detectors[map->DDL2Detector(ddl)-1] = kTRUE;
841     AliDebug(0,Form(" From reader: DDL number is %d , det ID is %d",ddl,detID));
842     
843     ULong_t  nPayloadWords = streamer.GetRCUPayloadSizeInSOD();
844     UChar_t* payloadData   = streamer.GetRCUPayloadInSOD();
845     UInt_t*  payloadWords  = reinterpret_cast<UInt_t*>(payloadData);
846     //UInt_t*   payloadWords  = streamer.GetRCUPayloadInSOD();
847
848     //std::cout<<nPayloadWords<<"    "<<ddl<<std::endl;
849     for (ULong_t i = 1; i <= nPayloadWords ; i++, payloadWords++) {
850       UInt_t payloadWord = *payloadWords; // Get32bitWord(i);
851     
852       //std::cout<<i<<Form("  word: 0x%x",payloadWord)<<std::endl;
853       // address is only 24 bit
854       UInt_t address       = (0xffffff & payloadWord);
855       UInt_t type          = ((address >> 21) & 0xf);
856       UInt_t error         = ((address >> 20) & 0x1);
857       UInt_t bcast         = ((address >> 18) & 0x1);
858       UInt_t bc_not_altro  = ((address >> 17) & 0x1);
859       UInt_t board         = ((address >> 12) & 0x1f);
860       UInt_t instruction   = 0;
861       UInt_t chip          = 0;
862       UInt_t channel       = 0;
863       if(bc_not_altro)
864         instruction        = address & 0xfff;
865       else {
866         chip               = ((address >> 9) & 0x7);
867         channel            = ((address >> 5) & 0x5);
868         instruction        = (address & 0x1f);
869       }
870         
871       Bool_t readDataWord = kFALSE;
872       switch(type) {
873       case 0x0: // Fec read
874         readDataWord = kTRUE;  
875       case 0x1: // Fec cmd
876       case 0x2: // Fec write
877         i++;  
878         payloadWords++;
879         break;
880       case 0x4: // Loop
881       case 0x5: // Wait
882         break;
883       case 0x6: // End sequence
884       case 0x7: // End Mem
885         i = nPayloadWords + 1;
886         break;
887       default:    
888         break;
889       }
890         
891       //Don't read unless we have a FEC_RD
892       if(!readDataWord)  continue;
893
894       UInt_t dataWord      = *payloadWords;//Get32bitWord(i);
895       UInt_t data          = (0xFFFFF & dataWord) ;
896       //UInt_t data          = (0xFFFF & dataWord) ;
897         
898       if(error) {
899         AliWarning(Form("error bit detected at Word 0x%06x; "
900                         "error % d, type %d, bc_not_altro %d, "
901                         "bcast %d, board 0x%02x, chip 0x%x, "
902                         "channel 0x%02x, instruction 0x%03x",
903                         address, error, type, bc_not_altro, 
904                         bcast,board,chip,channel,instruction));
905         //process error
906         continue;
907       }
908         
909         
910       switch(instruction) {
911           
912       case 0x01: break;  // First ADC T           
913       case 0x02: break; // I  3.3 V              
914       case 0x03: break; // I  2.5 V altro digital
915       case 0x04: break; // I  2.5 V altro analog 
916       case 0x05: break; // I  2.5 V VA           
917       case 0x06: break; // First ADC T           
918       case 0x07: break; // I  3.3 V              
919       case 0x08: break; // I  2.5 V altro digital
920       case 0x09: break; // I  2.5 V altro analog 
921       case 0x0A: break; // I  2.5 V VA           
922       case 0x2D: break; // Second ADC T           
923       case 0x2E: break; // I  1.5 V VA            
924       case 0x2F: break; // I -2.0 V               
925       case 0x30: break; // I -2.0 V VA            
926       case 0x31: break; //    2.5 V Digital driver
927       case 0x32: break; // Second ADC T           
928       case 0x33: break; // I  1.5 V VA            
929       case 0x34: break; // I -2.0 V               
930       case 0x35: break; // I -2.0 V VA            
931       case 0x36: break; //    2.5 V Digital driver
932       case 0x37: break; // Third ADC T             
933       case 0x38: break; // Temperature sens. 1     
934       case 0x39: break; // Temperature sens. 2     
935       case 0x3A: break; // U  2.5 altro digital (m)
936       case 0x3B: break; // U  2.5 altro analog (m) 
937       case 0x3C: break; // Third ADC T             
938       case 0x3D: break; // Temperature sens. 1     
939       case 0x3E: break; // Temperature sens. 2     
940       case 0x3F: break; // U  2.5 altro digital (m)
941       case 0x40: break; // U  2.5 altro analog (m) 
942       case 0x41: break; // Forth ADC T  
943       case 0x42: break; // U  2.5 VA (m)
944       case 0x43: break; // U  1.5 VA (m)
945       case 0x44: break; // U -2.0 VA (m)
946       case 0x45: break; // U -2.0 (m)   
947       case 0x46: break; // Forth ADC T  
948       case 0x47: break; // U  2.5 VA (m)
949       case 0x48: break; // U  1.5 VA (m)
950       case 0x49: break; // U -2.0 VA (m)
951       case 0x4A: break;  // U -2.0 (m)   
952         // Counters 
953       case 0x0B: break; // L1 trigger CouNTer
954       case 0x0C: break; // L2 trigger CouNTer
955       case 0x0D: break; // Sampling CLK CouNTer
956       case 0x0E: break; // DSTB CouNTer
957         // Test mode 
958       case 0x0F: break; // Test mode word
959       case 0x10: break; // Undersampling ratio.
960         // Configuration and status 
961       case 0x11: break; // Config/Status Register 0
962       case 0x12: break; // Config/Status Register 1
963       case 0x13: break; // Config/Status Register 2
964       case 0x14: break; // Config/Status Register 3
965       case 0x15: break; // Free
966         // Comands:
967       case 0x16: break; // Latch L1, L2, SCLK Counters
968       case 0x17: break; // Clear counters
969       case 0x18: break; // Clear CSR1
970       case 0x19: break; // rstb ALTROs
971       case 0x1A: break; // rstb BC
972       case 0x1B: break; // Start conversion
973       case 0x1C: break; // Scan event length
974       case 0x1D: break; // Read event length
975       case 0x1E: break; // Start test mode
976       case 0x1F: break; // Read acquisition memory
977         // FMD
978       case 0x20: break; // FMDD status
979       case 0x21: break; // L0 counters
980       case 0x22: break; // FMD: Wait to hold
981       case 0x23: break; // FMD: L1 timeout
982       case 0x24: break; // FMD: L2 timeout
983       case 0x25: // FMD: Shift clk 
984         shift_clk[board] = ((data >> 8 ) & 0xFF); 
985         AliDebug(30,Form("Read shift_clk=%d for board 0x%02x", 
986                          shift_clk[board], board));
987         break; 
988       case 0x26: // FMD: Strips 
989         strip_low[board]  = ((data >> 0 ) & 0xFF); 
990         strip_high[board] = ((data >> 8 ) & 0xFF);  
991         break; 
992       case 0x27: // FMD: Cal pulse 
993         pulse_size[board]  =  ((data >> 8 ) & 0xFF);
994         break; 
995       case 0x28: break; // FMD: Shape bias
996       case 0x29: break; // FMD: Shape ref
997       case 0x2A: break; // FMD: Preamp ref
998       case 0x2B: // FMD: Sample clk 
999         sample_clk[board] = ((data >> 8 ) & 0xFF); 
1000         AliDebug(30,Form("Read sample_clk=%d for board 0x%02x", 
1001                          sample_clk[board], board));
1002         break; 
1003       case 0x2C: break; // FMD: Commands
1004       case 0x4B: // FMD: Cal events 
1005         pulse_length[board] = ((data >> 0 ) & 0xFF);
1006         break; 
1007       default: break;
1008         
1009       }
1010       AliDebug(50,Form("instruction 0x%x, dataword 0x%x",
1011                        instruction,dataWord));
1012     } // End of loop over Result memory event
1013     
1014     UShort_t det    = 0;
1015     UShort_t sector = 0;
1016     Short_t  strip  = -1;
1017     Char_t   ring   = '\0';
1018    
1019     const UInt_t boards[4] = {0,1,16,17};
1020     for(Int_t i=0;i<4;i++) {
1021       if(ddl==0 && (i==1 || i==3)) continue;
1022
1023       UInt_t chip =0, channel=0;
1024       det = map->DDL2Detector(ddl);
1025       map->Channel2StripBase(boards[i], chip, channel, ring, sector, strip);
1026      
1027       UInt_t samplerate = 0;
1028 #if USE_VOTE
1029       if(sample_clk[boards[i]] == 0) {
1030         if(ddl == 0) {
1031           Int_t sample1 = sample_clk[boards[0]];
1032           Int_t sample2 = sample_clk[boards[2]];            
1033           if(sample1) sample_clk[boards[i]] = sample1;
1034           else sample_clk[boards[i]] = sample2;
1035         }
1036         else {
1037           Int_t sample1 = sample_clk[boards[0]];
1038           Int_t sample2 = sample_clk[boards[1]];
1039           Int_t sample3 = sample_clk[boards[2]];
1040           Int_t sample4 = sample_clk[boards[3]];
1041           Int_t agreement = 0;
1042           if(sample1 == sample2) agreement++;
1043           if(sample1 == sample3) agreement++;
1044           if(sample1 == sample4) agreement++;
1045           if(sample2 == sample3) agreement++;
1046           if(sample2 == sample4) agreement++;
1047           if(sample3 == sample4) agreement++;
1048             
1049           Int_t idx = 0;
1050           if(i<3) idx = i+1;
1051           else  idx = i-1;
1052           if(agreement == 3) {
1053             sample_clk[boards[i]] = sample_clk[boards[idx]];
1054             shift_clk[boards[i]] = shift_clk[boards[idx]];
1055             strip_low[boards[i]] = strip_low[boards[idx]];
1056             strip_high[boards[i]] = strip_high[boards[idx]];
1057             pulse_length[boards[i]] = pulse_length[boards[idx]];
1058             pulse_size[boards[i]] = pulse_size[boards[idx]];
1059             AliDebug(3,Form("Vote taken for ddl %d, board 0x%x",
1060                             ddl,boards[i]));
1061           }
1062         }
1063       } 
1064 #endif
1065       
1066       if(sample_clk[boards[i]])
1067         samplerate = shift_clk[boards[i]]/sample_clk[boards[i]];
1068       AliDebug(10,Form("Sample rate for board 0x%02x is %d", 
1069                       boards[i], samplerate));
1070       sampleRate->Set(det,ring,sector,0,samplerate);
1071       stripRange->Set(det,ring,sector,0,
1072                       strip_low[boards[i]],strip_high[boards[i]]);
1073       
1074       AliDebug(20,Form("det %d, ring %c, ",det,ring));
1075       pulseLength.AddAt(pulse_length[boards[i]],
1076                         GetHalfringIndex(det,ring,boards[i]/16));
1077       pulseSize.AddAt(pulse_size[boards[i]],
1078                       GetHalfringIndex(det,ring,boards[i]/16));
1079       
1080       
1081       
1082       AliDebug(20,Form(": Board: 0x%02x\n"
1083                        "\tstrip_low  %3d, strip_high   %3d\n"
1084                        "\tshift_clk  %3d, sample_clk   %3d\n"
1085                        "\tpulse_size %3d, pulse_length %3d",
1086                        boards[i], 
1087                        strip_low[boards[i]], strip_high[boards[i]],
1088                        shift_clk[boards[i]], sample_clk[boards[i]],
1089                        pulse_size[boards[i]],pulse_length[boards[i]]));
1090     }
1091     
1092   }
1093   
1094   AliFMDParameters::Instance()->SetSampleRate(sampleRate);
1095   AliFMDParameters::Instance()->SetStripRange(stripRange);
1096   
1097   AliDebug(0,Form("End of SOD/EOD"));
1098   
1099   return kTRUE;
1100 }
1101 //____________________________________________________________________
1102
1103 UInt_t AliFMDRawReader::Get32bitWord(Int_t idx)
1104 {
1105   // This method returns the 32 bit word at a given
1106   // position inside the raw data payload.
1107   // The 'index' points to the beginning of the next word.
1108   // The method is supposed to be endian (platform)
1109   // independent.
1110   if (!fData) {
1111     AliFatal("Raw data paylod buffer is not yet initialized !");
1112   }
1113
1114   Int_t index = 4*idx;
1115   
1116   if (index < 4) {
1117     //  fRawReader->AddFatalErrorLog(k32bitWordReadErr,Form("pos = %d",index));
1118     //    PrintDebug();
1119     AliWarning(Form("Invalid raw data payload position (%d) !",index));
1120   }
1121
1122   UInt_t word = 0;
1123    
1124   word  = fData[--index] << 24;
1125   word |= fData[--index] << 16;
1126   word |= fData[--index] << 8;
1127   word |= fData[--index] << 0 ;
1128
1129   return word;
1130 }
1131 //_____________________________________________________________________ 
1132 Int_t AliFMDRawReader::GetHalfringIndex(UShort_t det, Char_t ring, 
1133                                         UShort_t board) const
1134 {
1135   // 
1136   // Get short index for a given half-ring
1137   // 
1138   // Parameters:
1139   //    det   Detector number
1140   //    ring  Ring identifer
1141   //    board Board number
1142   // 
1143   // Return:
1144   //    
1145   //
1146   UShort_t iring  =  (ring == 'I' ? 1 : 0);
1147   
1148   Int_t index = (((det-1) << 2) | (iring << 1) | (board << 0));
1149   
1150   return index-2;
1151   
1152 }
1153
1154 //____________________________________________________________________
1155 // 
1156 // EOF
1157 //
1158