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