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