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