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