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