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