]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDRawReader.cxx
small corrections to previous checkin
[u/mrichter/AliRoot.git] / FMD / AliFMDRawReader.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15 /* $Id$ */
16 /** @file    AliFMDRawReader.cxx
17     @author  Christian Holm Christensen <cholm@nbi.dk>
18     @date    Mon Mar 27 12:45:23 2006
19     @brief   Class to read raw data 
20     @ingroup FMD_rec
21 */
22 //____________________________________________________________________
23 //
24 // Class to read ADC values from a AliRawReader object. 
25 //
26 // This class uses the AliFMDRawStreamer class to read the ALTRO
27 // formatted data. 
28 // 
29 //          +-------+
30 //          | TTask |
31 //          +-------+
32 //              ^
33 //              |
34 //      +-----------------+  <<references>>  +--------------+
35 //      | AliFMDRawReader |<>----------------| AliRawReader |
36 //      +-----------------+                  +--------------+
37 //              |                                  ^
38 //              | <<uses>>                         |
39 //              V                                  |
40 //      +-----------------+      <<uses>>          |
41 //      | AliFMDRawStream |------------------------+
42 //      +-----------------+
43 //              |
44 //              V
45 //      +----------------+
46 //      | AliAltroStream |
47 //      +----------------+
48 //
49 // #include <AliLog.h>          // ALILOG_H
50 #include "AliFMDDebug.h" // Better debug macros
51 #include "AliFMDParameters.h"   // ALIFMDPARAMETERS_H
52 #include "AliFMDDigit.h"        // ALIFMDDIGIT_H
53 #include "AliFMDSDigit.h"       // ALIFMDSDIGIT_H
54 // #include "AliFMDRawStream.h" // ALIFMDRAWSTREAM_H 
55 #include "AliRawReader.h"       // ALIRAWREADER_H 
56 #include "AliFMDRawReader.h"    // ALIFMDRAWREADER_H 
57 #include "AliFMDDebug.h"
58 #include "AliFMDCalibSampleRate.h"
59 #include "AliFMDCalibStripRange.h"
60 #include "AliFMDAltroMapping.h"
61 // #include "AliFMDAltroIO.h"   // ALIFMDALTROIO_H 
62 #include "AliAltroRawStreamV3.h"
63 #include <TArrayS.h>            // ROOT_TArrayS
64 #include <TTree.h>              // ROOT_TTree
65 #include <TClonesArray.h>       // ROOT_TClonesArray
66 #include <TString.h>
67 #include <iostream>
68 #include <climits>
69 // #include <iomanip>
70
71 //____________________________________________________________________
72 ClassImp(AliFMDRawReader)
73 #if 0
74   ; // This is here to keep Emacs for indenting the next line
75 #endif
76
77 //____________________________________________________________________
78 AliFMDRawReader::AliFMDRawReader(AliRawReader* reader, TTree* tree) 
79   : TTask("FMDRawReader", "Reader of Raw ADC values from the FMD"),
80     fTree(tree),
81     fReader(reader), 
82     // fSampleRate(1),
83     fData(0),
84     fNbytes(0), 
85     fMinStrip(0),
86     fMaxStrip(127), 
87     fPreSamp(14+5),
88     fSeen(0)
89 {
90   // Default CTOR
91   for (Int_t i = 0; i < 3; i++) { 
92     fSampleRate[i]   = 0;
93     fZeroSuppress[i] = kFALSE;
94     fNoiseFactor[i]  = 1;
95   }
96 }
97
98 //____________________________________________________________________
99 void
100 AliFMDRawReader::Exec(Option_t*) 
101 {
102   // Read the data 
103   TClonesArray* array = new TClonesArray("AliFMDDigit");
104   if (!fTree) {
105     AliError("No tree");
106     return;
107   }
108   fTree->Branch("FMD", &array);
109   
110   
111   ReadAdcs(array);
112   Int_t nWrite = fTree->Fill();
113   AliFMDDebug(1, ("Got a grand total of %d digits, wrote %d bytes to tree", 
114                    array->GetEntriesFast(), nWrite));
115 }
116
117 //____________________________________________________________________
118 Int_t
119 AliFMDRawReader::NewDDL(AliAltroRawStreamV3& input, UShort_t& det)
120 {
121   // Process a new DDL.  Sets the internal data members fZeroSuppress, 
122   // fSampleRate, and fNoiseFactor based on information in the RCU trailer. 
123   // 
124   // Parameters:
125   //   input Input stream
126   //   det   On return, the detector number
127   // 
128   // Return value:
129   //   negative value in case of problems, the DDL number otherwise
130
131   // Get the DDL number
132   UInt_t ddl = input.GetDDLNumber();
133   AliFMDDebug(2, ("DDL number %d", ddl));
134
135   // Get zero suppression flag
136   fZeroSuppress[ddl] = input.GetZeroSupp();
137   AliFMDDebug(3, ("RCU @ DDL %d zero suppression: %s", 
138                    ddl, (fZeroSuppress[ddl] ? "yes" : "no")));
139
140   // WARNING: We store the noise factor in the 2nd baseline
141   // filters excluded post samples, since we'll never use that
142   // mode. 
143   fNoiseFactor[ddl]  = input.GetNPostsamples();
144   AliFMDDebug(3, ("RCU @ DDL %d noise factor: %d", ddl,fNoiseFactor[ddl]));
145     
146   // WARNING: We store the sample rate in the number of pre-trigger
147   // samples, since we'll never use that mode.
148   fSampleRate[ddl]     = input.GetNPretriggerSamples();
149   AliFMDDebug(3, ("RCU @ DDL %d sample rate: %d", ddl,fSampleRate[ddl]));
150
151   // Get Errors seen 
152   Int_t nChAddrMismatch = input.GetNChAddrMismatch();
153   Int_t nChLenMismatch  = input.GetNChLengthMismatch();
154   if (nChAddrMismatch != 0) 
155     AliWarning(Form("Got %d channels with address mis-matches for 0x%03x",
156                     nChAddrMismatch, ddl));
157   if (nChLenMismatch != 0) 
158     AliWarning(Form("Got %d channels with length mis-matches for 0x%03x",
159                     nChLenMismatch, ddl));
160
161   // Map DDL number to the detector number 
162   AliFMDParameters*    pars   = AliFMDParameters::Instance();
163   AliFMDAltroMapping*  map    = pars->GetAltroMap();
164   if (map->DDL2Detector(ddl) < 0) return -1;
165   det = map->DDL2Detector(ddl);
166
167   if (AliLog::GetDebugLevel("FMD", 0) > 5) 
168     input.PrintRCUTrailer();
169   return ddl;
170 }
171
172 //____________________________________________________________________
173 Int_t
174 AliFMDRawReader::NewChannel(AliAltroRawStreamV3& input,  UShort_t det,
175                             Char_t&  ring, UShort_t& sec, Short_t& strbase)
176 {
177   // Processs a new channel.  Sets the internal data members
178   // fMinStrip, fMaxStrip, and fPreSamp. 
179   //
180   // Parameter:
181   //   input   Input stream
182   //   ring    On return, the ring identifier 
183   //   sec     On return, the sector number
184   //   strbase On return, the strip base
185   // 
186   // Return value
187   //   negative value in case of problems, hardware address otherwise
188
189   // Get the hardware address, and map that to detector coordinates 
190   UShort_t board, chip, channel;
191   Int_t    ddl    = input.GetDDLNumber();
192   Int_t    hwaddr = input.GetHWAddress();
193   if (input.IsChannelBad()) { 
194     AliError(Form("Channel 0x%03x is marked as bad!", hwaddr));
195   }
196   
197   AliFMDParameters*    pars   = AliFMDParameters::Instance();
198   AliFMDAltroMapping*  map    = pars->GetAltroMap();
199   // Map to hardware stuff 
200   map->ChannelAddress(hwaddr, board, chip, channel);
201   // Then try to map to detector address
202   if (!map->Channel2StripBase(board, chip, channel, ring, sec, strbase)) {
203     AliError(Form("Failed to get detector id from DDL %d, "
204                   "hardware address 0x%03x", ddl, hwaddr));
205     return -1;
206   }
207   AliFMDDebug(3, ("Board: 0x%02x, Altro: 0x%x, Channel: 0x%x", 
208                   board, chip, channel));
209
210   // Get the 'conditions'
211   fMinStrip = pars->GetMinStrip(det, ring, sec, strbase);
212   fMaxStrip = pars->GetMaxStrip(det, ring, sec, strbase);
213   fPreSamp  = pars->GetPreSamples(det, ring, sec, strbase);
214   if (fSampleRate[ddl] == 0) 
215     fSampleRate[ddl] = pars->GetSampleRate(det, ring, sec, strbase);
216
217   return hwaddr;
218 }
219
220 //____________________________________________________________________
221 Int_t
222 AliFMDRawReader::NewSample(AliAltroRawStreamV3& input, 
223                            Int_t i, UShort_t t, UShort_t sec,
224                            UShort_t  strbase, Short_t&  str, UShort_t& samp)
225 {
226   // Process a new timebin
227   // 
228   // Parameters:
229   //   input   Input stream
230   //   i       Index into bunch data
231   //   t       Time
232   //   strbase Base of strip numbers for this channel
233   //   str     On return, the strip number
234   //   samp    On return, the sample number
235   // 
236   // Return value
237   //   negative value in case of problems, ADC value otherwise
238   if (t < fPreSamp) return -1;
239
240   Int_t           ddl  = input.GetDDLNumber();
241   Int_t           hwa  = input.GetHWAddress();
242   const UShort_t* data = input.GetSignals();
243   Short_t         adc  = data[i];
244   AliFMDDebug(10, ("0x%04x/0x%03x/%04d %4d", ddl, hwa, t, adc));
245
246   AliFMDParameters*    pars   = AliFMDParameters::Instance();
247   AliFMDAltroMapping*  map    = pars->GetAltroMap();
248
249   samp            = 0;
250   Short_t  stroff = 0;
251   map->Timebin2Strip(sec, t, fPreSamp, fSampleRate[ddl], stroff, samp);
252   str             = strbase + stroff;
253       
254   AliFMDDebug(20, ("0x%04x/0x%03x/%04d=%4d maps to strip %3d sample %d " 
255                    "(pre: %d, min: %d, max: %d, rate: %d)",
256                   ddl, hwa, t, adc, str, samp, fPreSamp, 
257                   fMinStrip, fMaxStrip, fSampleRate[ddl]));
258   if (str < 0) { 
259     AliFMDDebug(10, ("Got presamples at timebin %d", i));
260     return -1;
261   }
262           
263   // VA1 Local strip number 
264   Short_t lstrip = (t - fPreSamp) / fSampleRate[ddl] + fMinStrip;
265       
266   AliFMDDebug(15, ("Checking if strip %d (%d) in range [%d,%d]", 
267                    lstrip, str, fMinStrip, fMaxStrip));
268   if (lstrip < fMinStrip || lstrip > fMaxStrip) {
269     AliFMDDebug(10, ("Strip %03d-%d (%d,%d) from t=%d out of range (%3d->%3d)", 
270                     str, samp, lstrip, stroff, t, fMinStrip, fMaxStrip));
271     adc = -1;
272   }
273   // Possibly do pedestal subtraction of signal 
274   if (adc > 1023) 
275     AliWarning(Form("ADC value out of range: %4d", adc));
276   return adc;
277 }
278
279 //____________________________________________________________________
280 Bool_t
281 AliFMDRawReader::NextSample(UShort_t& det, Char_t&   rng, UShort_t& sec, 
282                             UShort_t& str, UShort_t& sam, UShort_t& rat, 
283                             Short_t&  adc, Bool_t&   zs,  UShort_t& fac)
284 {
285   // Scan current event for next signal.   It returns kFALSE when
286   // there's no more data in the event. 
287   static AliAltroRawStreamV3 stream(fReader); //    = 0;
288   static Int_t               ddl      = -1;
289   static UShort_t            tdet     = 0;
290   static Char_t              trng     = '\0';
291   static UShort_t            tsec     = 0;
292   static Short_t             tstr     = 0;   
293   static Short_t             bstr     = -1;
294   static UShort_t            tsam     = 0;   
295   static UInt_t              trate    = 0;
296   static Int_t               hwaddr   = -1;
297   static UShort_t            start    = 0;
298   static UShort_t            length   = 0;
299   static Short_t             t        = -1;
300   static Int_t               i        = 0; 
301   // First entry!
302   if (stream.GetDDLNumber() < 0) { 
303     fReader->Select("FMD");
304
305     // Reset variables
306     ddl    = -1;  
307     trate  = 0;   
308     tdet   = 0;   
309     trng   = '\0';
310     tsec   = 0;   
311     tstr   = 0;  
312     tsam   = -1;
313     hwaddr = -1;
314   }
315
316   do { 
317     if (t < start - length + 1) { 
318       if (!stream.NextBunch()) { 
319         if (!stream.NextChannel()) { 
320           if (!stream.NextDDL()) {
321             stream.Reset();
322             return kFALSE;
323           }
324           ddl = NewDDL(stream, tdet);
325         }
326         hwaddr = NewChannel(stream, tdet, trng, tsec, bstr);
327       }
328       start  = stream.GetStartTimeBin();
329       length = stream.GetBunchLength();
330       t      = start;
331       i      = 0;
332     }
333     Int_t tadc = NewSample(stream, i, t, tsec, bstr, tstr, tsam);
334     if (tadc >= 0) { 
335       det = tdet;
336       rng = trng;
337       sec = tsec;
338       str = tstr;
339       sam = tsam;
340       adc = tadc;
341       rat = fSampleRate[ddl];
342       zs  = fZeroSuppress[ddl];
343       fac = fNoiseFactor[ddl];
344       t--;
345       i++;
346       break;
347     }
348   } while (true);
349
350   return kTRUE;
351 }
352
353
354 //____________________________________________________________________
355 Bool_t
356 AliFMDRawReader::NextSignal(UShort_t& det, Char_t&   rng, 
357                             UShort_t& sec, UShort_t& str, 
358                             Short_t&  adc, Bool_t&   zs, 
359                             UShort_t& fac)
360 {
361   
362   do { 
363     UShort_t samp, rate;
364     if (!NextSample(det, rng, sec, str, samp, rate, adc, zs, fac)) 
365       return kFALSE;
366
367     Bool_t take = kFALSE;
368     switch (rate) { 
369     case 1:                      take = kTRUE; break;
370     case 2:  if (samp == 1)      take = kTRUE; break;
371     case 3:  if (samp == 1)      take = kTRUE; break; 
372     case 4:  if (samp == 2)      take = kTRUE; break;
373     default: if (samp == rate-2) take = kTRUE; break;
374     }
375     if (!take) continue;
376     break;
377   } while (true);
378   return kTRUE;
379 }
380
381 //____________________________________________________________________
382 Bool_t
383 AliFMDRawReader::SelectSample(UShort_t samp, UShort_t rate) 
384 {
385   // Check if the passed sample is the one we need
386   Bool_t take = kFALSE;
387   switch (rate) { 
388   case 1:                      take = kTRUE; break;
389   case 2:  if (samp == 1)      take = kTRUE; break;
390   case 3:  if (samp == 1)      take = kTRUE; break; 
391   case 4:  if (samp == 2)      take = kTRUE; break;
392   default: if (samp == rate-2) take = kTRUE; break;
393   }
394   
395   return take;
396 }
397   
398 //____________________________________________________________________
399 Bool_t
400 AliFMDRawReader::ReadAdcs(TClonesArray* array) 
401 {
402   // Read ADC values from raw input into passed TClonesArray of AliFMDDigit
403   // objects. 
404   AliFMDDebug(2, ("Reading ADC values into a TClonesArray"));
405
406   // Read raw data into the digits array, using AliFMDAltroReader. 
407   if (!array) {
408     AliError("No TClonesArray passed");
409     return kFALSE;
410   }
411   const UShort_t kUShortMax = (1 << 16) - 1;
412   fSeen.Reset(kUShortMax);
413
414   AliAltroRawStreamV3  input(fReader);
415   input.Reset();
416   input.SelectRawData("FMD");
417   
418   // Loop over input RORCs
419   while (input.NextDDL()) { 
420     UShort_t det = 0;
421     Int_t    ddl = NewDDL(input, det);
422     if (ddl < 0) break;
423
424     while (input.NextChannel()) { 
425       // Get the hardware address, and map that to detector coordinates 
426       Char_t   ring;
427       UShort_t sec;
428       Short_t  strbase;
429       Int_t    hwaddr = NewChannel(input, det, ring, sec, strbase);
430       if (hwaddr < 0) break;
431
432       // Loop over bunches 
433       while (input.NextBunch()) { 
434         // Get Lenght of bunch, and pointer to the data 
435         const UShort_t* data   = input.GetSignals();
436         UShort_t        start  = input.GetStartTimeBin();
437         UShort_t        length = input.GetBunchLength();
438       
439         // Loop over the data and store it. 
440         for (Int_t i = 0; i < length; i++) { 
441           // Time 
442           Short_t  str;
443           UShort_t samp;
444           Int_t    t    = start - i;
445           Int_t    adc  = NewSample(input, i, t, sec, strbase, str, samp);
446           if (adc < 0) continue;
447           UShort_t counts = adc;
448       
449           AliFMDDebug(10, ("FMD%d%c[%02d,%03d]-%d: %4d", 
450                            det, ring, sec, str, samp, counts));
451           // Check the cache of indicies
452           Int_t idx = fSeen(det, ring, sec, str);
453           AliFMDDigit* digit = 0;
454           if (idx == kUShortMax) { 
455             // We haven't seen this strip yet. 
456             fSeen(det, ring, sec, str) = idx = array->GetEntriesFast();
457             AliFMDDebug(7,("making digit for FMD%d%c[%2d,%3d]-%d "
458                            "from timebin %4d", 
459                            det, ring, sec, str, samp, t));
460             digit = new ((*array)[idx]) AliFMDDigit(det, ring, sec, str);
461             digit->SetDefaultCounts(fSampleRate[ddl]);
462           }
463           else 
464             digit = static_cast<AliFMDDigit*>(array->At(idx));
465           AliFMDDebug(10, ("Setting FMD%d%c[%2d,%3d]-%d "
466                            "from timebin %4d=%4d (%4d)", 
467                            det, ring, sec, str, samp, t, counts, data[i]));
468           digit->SetCount(samp, counts);
469         } // for (i)
470       } // while (bunch)
471     } // while (channel)
472   } // while (ddl)
473   return kTRUE;
474 }
475
476 //____________________________________________________________________
477 Bool_t AliFMDRawReader::ReadSODevent(AliFMDCalibSampleRate* sampleRate, 
478                                      AliFMDCalibStripRange* stripRange, 
479                                      TArrayS &pulseSize, 
480                                      TArrayS &pulseLength, 
481                                      Bool_t* detectors) 
482 {
483
484   AliFMDDebug(0, ("Start of SOD/EOD"));
485   
486   UInt_t shift_clk[18];
487   UInt_t sample_clk[18];
488   UInt_t strip_low[18];
489   UInt_t strip_high[18];
490   UInt_t pulse_size[18];
491   UInt_t pulse_length[18];  
492   AliFMDParameters*   param = AliFMDParameters::Instance();
493   AliFMDAltroMapping* map   = param->GetAltroMap();
494   
495   AliAltroRawStreamV3  streamer(fReader);
496   streamer.Reset();
497   streamer.SelectRawData("FMD");
498   //fReader->GetDDLID();
499   while (streamer.NextDDL()) {
500     Int_t ddl   = streamer.GetDDLNumber();
501     Int_t detID = fReader->GetDetectorID();
502     if (detectors) detectors[map->DDL2Detector(ddl)-1] = kTRUE;
503     AliFMDDebug(0, (" From reader: DDL number is %d , det ID is %d",ddl,detID));
504
505     ULong_t  nPayloadWords = streamer.GetRCUPayloadSizeInSOD();
506     UChar_t* payloadData   = streamer.GetRCUPayloadInSOD();
507     UInt_t*  payloadWords  = reinterpret_cast<UInt_t*>(payloadData);
508     //UInt_t*   payloadWords  = streamer.GetRCUPayloadInSOD();
509
510     //std::cout<<nPayloadWords<<"    "<<ddl<<std::endl;
511     for (ULong_t i = 1; i <= nPayloadWords ; i++, payloadWords++) {
512       UInt_t payloadWord = *payloadWords; // Get32bitWord(i);
513     
514       //std::cout<<i<<Form("  word: 0x%x",payloadWord)<<std::endl;
515       // address is only 24 bit
516       UInt_t address       = (0xffffff & payloadWord);
517       UInt_t type          = ((address >> 21) & 0xf);
518       UInt_t error         = ((address >> 20) & 0x1);
519       UInt_t bcast         = ((address >> 18) & 0x1);
520       UInt_t bc_not_altro  = ((address >> 17) & 0x1);
521       UInt_t board         = ((address >> 12) & 0x1f);
522       UInt_t instruction   = 0;
523       UInt_t chip          = 0;
524       UInt_t channel       = 0;
525       if(bc_not_altro)
526         instruction        = address & 0xfff;
527       else {
528         chip               = ((address >> 9) & 0x7);
529         channel            = ((address >> 5) & 0x5);
530         instruction        = (address & 0x1f);
531       }
532         
533       Bool_t readDataWord = kFALSE;
534       switch(type) {
535       case 0x0: // Fec read
536         readDataWord = kTRUE;  
537       case 0x1: // Fec cmd
538       case 0x2: // Fec write
539         i++;  
540         payloadWords++;
541         break;
542       case 0x4: // Loop
543       case 0x5: // Wait
544         break;
545       case 0x6: // End sequence
546       case 0x7: // End Mem
547         i = nPayloadWords + 1;
548         break;
549       default:    
550         break;
551       }
552         
553       //Don't read unless we have a FEC_RD
554       if(!readDataWord)  continue;
555
556       UInt_t dataWord      = *payloadWords;//Get32bitWord(i);
557       UInt_t data          = (0xFFFFF & dataWord) ;
558       //UInt_t data          = (0xFFFF & dataWord) ;
559         
560       if(error) {
561         AliWarning(Form("error bit detected at Word 0x%06x; "
562                         "error % d, type %d, bc_not_altro %d, "
563                         "bcast %d, board 0x%02x, chip 0x%x, "
564                         "channel 0x%02x, instruction 0x%03x",
565                         address, error, type, bc_not_altro, 
566                         bcast,board,chip,channel,instruction));
567         //process error
568         continue;
569       }
570         
571         
572       switch(instruction) {
573           
574       case 0x01: break;  // First ADC T           
575       case 0x02: break; // I  3.3 V              
576       case 0x03: break; // I  2.5 V altro digital
577       case 0x04: break; // I  2.5 V altro analog 
578       case 0x05: break; // I  2.5 V VA           
579       case 0x06: break; // First ADC T           
580       case 0x07: break; // I  3.3 V              
581       case 0x08: break; // I  2.5 V altro digital
582       case 0x09: break; // I  2.5 V altro analog 
583       case 0x0A: break; // I  2.5 V VA           
584       case 0x2D: break; // Second ADC T           
585       case 0x2E: break; // I  1.5 V VA            
586       case 0x2F: break; // I -2.0 V               
587       case 0x30: break; // I -2.0 V VA            
588       case 0x31: break; //    2.5 V Digital driver
589       case 0x32: break; // Second ADC T           
590       case 0x33: break; // I  1.5 V VA            
591       case 0x34: break; // I -2.0 V               
592       case 0x35: break; // I -2.0 V VA            
593       case 0x36: break; //    2.5 V Digital driver
594       case 0x37: break; // Third ADC T             
595       case 0x38: break; // Temperature sens. 1     
596       case 0x39: break; // Temperature sens. 2     
597       case 0x3A: break; // U  2.5 altro digital (m)
598       case 0x3B: break; // U  2.5 altro analog (m) 
599       case 0x3C: break; // Third ADC T             
600       case 0x3D: break; // Temperature sens. 1     
601       case 0x3E: break; // Temperature sens. 2     
602       case 0x3F: break; // U  2.5 altro digital (m)
603       case 0x40: break; // U  2.5 altro analog (m) 
604       case 0x41: break; // Forth ADC T  
605       case 0x42: break; // U  2.5 VA (m)
606       case 0x43: break; // U  1.5 VA (m)
607       case 0x44: break; // U -2.0 VA (m)
608       case 0x45: break; // U -2.0 (m)   
609       case 0x46: break; // Forth ADC T  
610       case 0x47: break; // U  2.5 VA (m)
611       case 0x48: break; // U  1.5 VA (m)
612       case 0x49: break; // U -2.0 VA (m)
613       case 0x4A: break;  // U -2.0 (m)   
614         // Counters 
615       case 0x0B: break; // L1 trigger CouNTer
616       case 0x0C: break; // L2 trigger CouNTer
617       case 0x0D: break; // Sampling CLK CouNTer
618       case 0x0E: break; // DSTB CouNTer
619         // Test mode 
620       case 0x0F: break; // Test mode word
621       case 0x10: break; // Undersampling ratio.
622         // Configuration and status 
623       case 0x11: break; // Config/Status Register 0
624       case 0x12: break; // Config/Status Register 1
625       case 0x13: break; // Config/Status Register 2
626       case 0x14: break; // Config/Status Register 3
627       case 0x15: break; // Free
628         // Comands:
629       case 0x16: break; // Latch L1, L2, SCLK Counters
630       case 0x17: break; // Clear counters
631       case 0x18: break; // Clear CSR1
632       case 0x19: break; // rstb ALTROs
633       case 0x1A: break; // rstb BC
634       case 0x1B: break; // Start conversion
635       case 0x1C: break; // Scan event length
636       case 0x1D: break; // Read event length
637       case 0x1E: break; // Start test mode
638       case 0x1F: break; // Read acquisition memory
639         // FMD
640       case 0x20: break; // FMDD status
641       case 0x21: break; // L0 counters
642       case 0x22: break; // FMD: Wait to hold
643       case 0x23: break; // FMD: L1 timeout
644       case 0x24: break; // FMD: L2 timeout
645       case 0x25: // FMD: Shift clk 
646         shift_clk[board] = ((data >> 8 ) & 0xFF); 
647         break; 
648       case 0x26: // FMD: Strips 
649         strip_low[board]  = ((data >> 0 ) & 0xFF); 
650         strip_high[board] = ((data >> 8 ) & 0xFF);  
651         break; 
652       case 0x27: // FMD: Cal pulse 
653         pulse_size[board]  =  ((data >> 8 ) & 0xFF);
654         break; 
655       case 0x28: break; // FMD: Shape bias
656       case 0x29: break; // FMD: Shape ref
657       case 0x2A: break; // FMD: Preamp ref
658       case 0x2B: // FMD: Sample clk 
659         sample_clk[board] = ((data >> 8 ) & 0xFF); 
660         break; 
661       case 0x2C: break; // FMD: Commands
662       case 0x4B: // FMD: Cal events 
663         pulse_length[board] = ((data >> 0 ) & 0xFF);
664         break; 
665       default: break;
666         
667       }
668       AliFMDDebug(50, ("instruction 0x%x, dataword 0x%x",
669                        instruction,dataWord));
670     } // End of loop over Result memory event
671     
672     UShort_t det,sector;
673     Short_t strip;
674     Char_t ring;
675    
676     const UInt_t boards[4] = {0,1,16,17};
677     for(Int_t i=0;i<4;i++) {
678       if(ddl==0 && (i==1 || i==3)) continue;
679
680       UInt_t chip =0, channel=0;
681       det = map->DDL2Detector(ddl);
682       map->Channel2StripBase(boards[i], chip, channel, ring, sector, strip);
683      
684       UInt_t samplerate = 1;
685       if(sample_clk[boards[i]] == 0) {
686         if(ddl == 0) {
687           Int_t sample1 = sample_clk[boards[0]];
688           Int_t sample2 = sample_clk[boards[2]];            
689           if(sample1) sample_clk[boards[i]] = sample1;
690           else sample_clk[boards[i]] = sample2;
691         }
692         else {
693           Int_t sample1 = sample_clk[boards[0]];
694           Int_t sample2 = sample_clk[boards[1]];
695           Int_t sample3 = sample_clk[boards[2]];
696           Int_t sample4 = sample_clk[boards[3]];
697           Int_t agreement = 0;
698           if(sample1 == sample2) agreement++;
699           if(sample1 == sample3) agreement++;
700           if(sample1 == sample4) agreement++;
701           if(sample2 == sample3) agreement++;
702           if(sample2 == sample4) agreement++;
703           if(sample3 == sample4) agreement++;
704             
705           Int_t idx = 0;
706           if(i<3) idx = i+1;
707           else  idx = i-1;
708           if(agreement == 3) {
709             sample_clk[boards[i]] = sample_clk[boards[idx]];
710             shift_clk[boards[i]] = shift_clk[boards[idx]];
711             strip_low[boards[i]] = strip_low[boards[idx]];
712             strip_high[boards[i]] = strip_high[boards[idx]];
713             pulse_length[boards[i]] = pulse_length[boards[idx]];
714             pulse_size[boards[i]] = pulse_size[boards[idx]];
715             AliFMDDebug(0, ("Vote taken for ddl %d, board 0x%x",
716                             ddl,boards[i]));
717           }
718         }
719       } 
720       
721       if(sample_clk[boards[i]])
722         samplerate = shift_clk[boards[i]]/sample_clk[boards[i]];
723       sampleRate->Set(det,ring,sector,0,samplerate);
724       stripRange->Set(det,ring,sector,0,
725                       strip_low[boards[i]],strip_high[boards[i]]);
726       
727       AliFMDDebug(20, ("det %d, ring %c, ",det,ring));
728       pulseLength.AddAt(pulse_length[boards[i]],
729                         GetHalfringIndex(det,ring,boards[i]/16));
730       pulseSize.AddAt(pulse_size[boards[i]],
731                       GetHalfringIndex(det,ring,boards[i]/16));
732       
733       
734       
735       AliFMDDebug(20, (": Board: 0x%02x\n"
736                        "\tstrip_low  %3d, strip_high   %3d\n"
737                        "\tshift_clk  %3d, sample_clk   %3d\n"
738                        "\tpulse_size %3d, pulse_length %3d",
739                        boards[i], 
740                        strip_low[boards[i]], strip_high[boards[i]],
741                        shift_clk[boards[i]], sample_clk[boards[i]],
742                        pulse_size[boards[i]],pulse_length[boards[i]]));
743       }
744     
745     }
746   
747   AliFMDParameters::Instance()->SetSampleRate(sampleRate);
748   AliFMDParameters::Instance()->SetStripRange(stripRange);
749   
750   AliFMDDebug(0, ("End of SOD/EOD"));
751   
752   return kTRUE;
753 }
754 //____________________________________________________________________
755
756 UInt_t AliFMDRawReader::Get32bitWord(Int_t idx)
757 {
758   // This method returns the 32 bit word at a given
759   // position inside the raw data payload.
760   // The 'index' points to the beginning of the next word.
761   // The method is supposed to be endian (platform)
762   // independent.
763   if (!fData) {
764     AliFatal("Raw data paylod buffer is not yet initialized !");
765   }
766
767   Int_t index = 4*idx;
768   
769   if (index < 4) {
770     //  fRawReader->AddFatalErrorLog(k32bitWordReadErr,Form("pos = %d",index));
771     //    PrintDebug();
772     AliWarning(Form("Invalid raw data payload position (%d) !",index));
773   }
774
775   UInt_t word = 0;
776    
777   word  = fData[--index] << 24;
778   word |= fData[--index] << 16;
779   word |= fData[--index] << 8;
780   word |= fData[--index] << 0 ;
781
782   return word;
783 }
784 //_____________________________________________________________________ 
785 Int_t AliFMDRawReader::GetHalfringIndex(UShort_t det, Char_t ring, 
786                                         UShort_t board) {
787
788   UShort_t iring  =  (ring == 'I' ? 1 : 0);
789   
790   Int_t index = (((det-1) << 2) | (iring << 1) | (board << 0));
791   
792   return index-2;
793   
794 }
795
796 //____________________________________________________________________
797 // 
798 // EOF
799 //
800