1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
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
22 //____________________________________________________________________
24 // Class to read ADC values from a AliRawReader object.
26 // This class uses the AliFMDRawStreamer class to read the ALTRO
34 // +-----------------+ <<references>> +--------------+
35 // | AliFMDRawReader |<>----------------| AliRawReader |
36 // +-----------------+ +--------------+
40 // +-----------------+ <<uses>> |
41 // | AliFMDRawStream |------------------------+
42 // +-----------------+
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 "AliFMDSDigit.h" // ALIFMDSDIGIT_H
54 // #include "AliFMDRawStream.h" // ALIFMDRAWSTREAM_H
55 #include "AliRawReader.h" // ALIRAWREADER_H
56 #include "AliFMDRawReader.h" // ALIFMDRAWREADER_H
57 #include "AliFMDDebug.h"
58 #include "AliFMDCalibSampleRate.h"
59 #include "AliFMDCalibStripRange.h"
60 #include "AliFMDAltroMapping.h"
61 #include "AliFMDUShortMap.h"
62 // #include "AliFMDAltroIO.h" // ALIFMDALTROIO_H
63 #include "AliAltroRawStreamV3.h"
64 #include <TArrayS.h> // ROOT_TArrayS
65 #include <TTree.h> // ROOT_TTree
66 #include <TClonesArray.h> // ROOT_TClonesArray
72 //____________________________________________________________________
73 ClassImp(AliFMDRawReader)
75 ; // This is here to keep Emacs for indenting the next line
78 //____________________________________________________________________
79 AliFMDRawReader::AliFMDRawReader(AliRawReader* reader, TTree* tree)
80 : TTask("FMDRawReader", "Reader of Raw ADC values from the FMD"),
92 for (Int_t i = 0; i < 3; i++) {
94 fZeroSuppress[i] = kFALSE;
99 //____________________________________________________________________
101 AliFMDRawReader::Exec(Option_t*)
104 TClonesArray* array = new TClonesArray("AliFMDDigit");
109 fTree->Branch("FMD", &array);
113 Int_t nWrite = fTree->Fill();
114 AliDebug(1,Form("Got a grand total of %d digits, wrote %d bytes to tree",
115 array->GetEntriesFast(), nWrite));
119 //____________________________________________________________________
121 AliFMDRawReader::NewDDL(AliAltroRawStreamV3& input, UShort_t& det)
123 // Process a new DDL. Sets the internal data members fZeroSuppress,
124 // fSampleRate, and fNoiseFactor based on information in the RCU trailer.
127 // input Input stream
128 // det On return, the detector number
131 // negative value in case of problems, the DDL number otherwise
133 // Get the DDL number
134 UInt_t ddl = input.GetDDLNumber();
135 AliDebug(2,Form("DDL number %d", ddl));
137 // Note, previously, the ALTROCFG1 register was interpreted as
139 // Bits Value Description
140 // 0- 3 0/1 1st Baseline filter, mode
141 // 4- 5 Over-1 2nd baseline filter, # of pre-samples
142 // 6- 9 factor 2nd baseline filter, # of post-samples
143 // 10- 0 2nd baseline filter, enable
144 // 11-12 00 Zero suppression, glitch filter mode
145 // 13-15 001 Zero suppression, # of post samples
146 // 16-17 01 Zero suppression, # of pre samples
147 // 18 0/1 Zero suppression, enable
149 // The interpretation used in AliAltroRawStreamerV3 - which
150 // corresponds directly to ALTRO DPCFG register - is
152 // Bits Value Description
153 // 0- 3 0/1 1st Baseline filter, mode
154 // 4 0 Polarity (if '1', then "1's inverse")
155 // 5- 6 01 Zero suppression, # of pre samples
156 // 7-10 0001 Zero suppression, # of post samples
157 // 11 0 2nd baseline filter, enable
158 // 12-13 00 Zero suppression, glitch filter mode
159 // 14-16 factor 2nd baseline filter, # of post-samples
160 // 17-18 01 2nd baseline filter, # of pre-samples
161 // 19 0/1 Zero suppression, enable
163 // Writing 'x' for variable values, that means we have the
164 // following patterns for the 2 cases
166 // bit # 20 16 12 8 4 0
167 // old |0x01|0010|00xx|xxxx|xxxx|
168 // new |x01x|xx00|0000|1010|xxxx|
170 // That means that we can check if bits 10-13 are '1000' or
171 // '0000', which will tell us if the value was written with the
172 // new or the old interpretation. That is, we can check that
174 // if (((altrocfg1 >> 10) & 0x8) == 0x8) {
175 // // old interpretation
178 // // New interpretation
181 // That means, that we should never
183 // - change the # of zero suppression post samples
184 // - Turn on 2nd baseline correction
185 // - Change the zero-suppression glitch filter mode
187 // This change as introduced in version 1.2 of Rcu++
189 UInt_t cfg1 = input.GetAltroCFG1();
190 if (((cfg1 >> 10) & 0x8) == 0x8) {
191 UInt_t cfg2 = input.GetAltroCFG2();
192 AliDebug(3,Form("We have data from older MiniConf 0x%x cfg2=0x%08x",
193 ((cfg1 >> 10) & 0x8), cfg2));
194 fZeroSuppress[ddl] = (cfg1 >> 0) & 0x1;
195 fNoiseFactor[ddl] = (cfg1 >> 6) & 0xF;
196 fSampleRate[ddl] = (cfg2 >> 20) & 0xF;
199 AliDebug(3,Form("We have data from newer MiniConf 0x%x",
200 ((cfg1 >> 10) & 0x8)));
201 fZeroSuppress[ddl] = input.GetZeroSupp();
202 // WARNING: We store the noise factor in the 2nd baseline
203 // filters excluded post samples, since we'll never use that
205 fNoiseFactor[ddl] = input.GetNPostsamples();
206 // WARNING: We store the sample rate in the number of pre-trigger
207 // samples, since we'll never use that mode.
208 fSampleRate[ddl] = input.GetNPretriggerSamples();
211 AliDebug(10,Form("Phase of DDL=%d is %g (%d)", ddl, input.GetL1Phase(),
212 input.GetAltroCFG2() & 0x1F));
213 fL1Phase[ddl] = input.GetAltroCFG2() & 0x1F; // input.GetL1Phase();
214 AliDebug(3,Form("RCU @ DDL %d zero suppression: %s",
215 ddl, (fZeroSuppress[ddl] ? "yes" : "no")));
216 AliDebug(3,Form("RCU @ DDL %d noise factor: %d", ddl,fNoiseFactor[ddl]));
217 AliDebug(3,Form("RCU @ DDL %d sample rate: %d", ddl,fSampleRate[ddl]));
221 Int_t nChAddrMismatch = input.GetNChAddrMismatch();
222 Int_t nChLenMismatch = input.GetNChLengthMismatch();
223 if (nChAddrMismatch != 0)
224 AliWarning(Form("Got %d channels with address mis-matches for 0x%03x",
225 nChAddrMismatch, ddl));
226 if (nChLenMismatch != 0)
227 AliWarning(Form("Got %d channels with length mis-matches for 0x%03x",
228 nChLenMismatch, ddl));
230 // Map DDL number to the detector number
231 AliFMDParameters* pars = AliFMDParameters::Instance();
232 AliFMDAltroMapping* map = pars->GetAltroMap();
233 if (map->DDL2Detector(ddl) < 0) return -1;
234 det = map->DDL2Detector(ddl);
236 if (AliLog::GetDebugLevel("FMD", 0) > 5)
237 input.PrintRCUTrailer();
241 //____________________________________________________________________
243 AliFMDRawReader::NewChannel(const AliAltroRawStreamV3& input, UShort_t det,
244 Char_t& ring, UShort_t& sec, Short_t& strbase)
246 // Processs a new channel. Sets the internal data members
247 // fMinStrip, fMaxStrip, and fPreSamp.
250 // input Input stream
251 // ring On return, the ring identifier
252 // sec On return, the sector number
253 // strbase On return, the strip base
256 // negative value in case of problems, hardware address otherwise
258 // Get the hardware address, and map that to detector coordinates
259 UShort_t board, chip, channel;
260 Int_t ddl = input.GetDDLNumber();
261 Int_t hwaddr = input.GetHWAddress();
262 if (input.IsChannelBad()) {
263 const char* msg = Form("Ignoring channel %03d/0x%03x with errors",
266 if (AliDebugLevel() > 10) input.HexDumpChannel();
267 fReader->AddMinorErrorLog(AliAltroRawStreamV3::kAltroPayloadErr,msg);
272 AliFMDParameters* pars = AliFMDParameters::Instance();
273 AliFMDAltroMapping* map = pars->GetAltroMap();
274 // Map to hardware stuff
275 map->ChannelAddress(hwaddr, board, chip, channel);
276 // Then try to map to detector address
277 if (!map->Channel2StripBase(board, chip, channel, ring, sec, strbase)) {
278 AliError(Form("Failed to get detector id from DDL %d, "
279 "hardware address 0x%03x", ddl, hwaddr));
282 AliDebug(4,Form("Board: 0x%02x, Altro: 0x%x, Channel: 0x%x",
283 board, chip, channel));
285 // Get the 'conditions'
286 fMinStrip = pars->GetMinStrip(det, ring, sec, strbase);
287 fMaxStrip = pars->GetMaxStrip(det, ring, sec, strbase);
288 fPreSamp = pars->GetPreSamples(det, ring, sec, strbase);
289 if (fSampleRate[ddl] == 0)
290 fSampleRate[ddl] = pars->GetSampleRate(det, ring, sec, strbase);
295 //____________________________________________________________________
297 AliFMDRawReader::NewBunch(const AliAltroRawStreamV3& input,
298 UShort_t& start, UShort_t& length)
301 // Do some checks on the bunch data
303 Int_t ddl = input.GetDDLNumber();
304 Int_t hwaddr = input.GetHWAddress();
305 UShort_t nSamples = input.GetNSamplesPerCh() + fPreSamp;
306 UShort_t tstart = input.GetStartTimeBin();
307 length = input.GetBunchLength();
309 if (tstart >= nSamples) {
310 const char* msg = Form("Bunch in %03d/0x%03x has an start time greater "
311 "than number of samples: 0x%x >= 0x%x",
312 ddl, hwaddr, tstart, nSamples);
314 if (AliDebugLevel() > 10) input.HexDumpChannel();
315 fReader->AddMinorErrorLog(AliAltroRawStreamV3::kAltroPayloadErr,msg);
319 if ((int(tstart) - length + 1) < 0) {
320 const char* msg = Form("Bunch in %03d/0x%03x has an invalid length and "
321 "start time: 0x%x,0x%x (%d-%d+1=%d<0)",
322 ddl, hwaddr, length, tstart, tstart, length,
323 int(tstart)-length+1);
325 if (AliDebugLevel() > 10) input.HexDumpChannel();
326 fReader->AddMinorErrorLog(AliAltroRawStreamV3::kAltroPayloadErr,msg);
330 if (tstart >= start) {
331 const char* msg = Form("Bunch in %03d/0x%03x has early start time: "
332 "0x%x >= 0x%x", ddl, hwaddr, tstart, start);
334 if (AliDebugLevel() > 10) input.HexDumpChannel();
335 fReader->AddMinorErrorLog(AliAltroRawStreamV3::kAltroPayloadErr,msg);
343 //____________________________________________________________________
345 AliFMDRawReader::NewSample(const AliAltroRawStreamV3& input,
346 Int_t i, UShort_t t, UShort_t sec,
347 UShort_t strbase, Short_t& str, UShort_t& samp)
349 // Process a new timebin
352 // input Input stream
353 // i Index into bunch data
355 // strbase Base of strip numbers for this channel
356 // str On return, the strip number
357 // samp On return, the sample number
360 // negative value in case of problems, ADC value otherwise
361 if (t < fPreSamp) return -1;
363 Int_t ddl = input.GetDDLNumber();
364 Int_t hwa = input.GetHWAddress();
365 const UShort_t* data = input.GetSignals();
366 Short_t adc = data[i];
367 AliDebug(10,Form("0x%04x/0x%03x/%04d %4d", ddl, hwa, t, adc));
369 AliFMDParameters* pars = AliFMDParameters::Instance();
370 AliFMDAltroMapping* map = pars->GetAltroMap();
374 map->Timebin2Strip(sec, t, fPreSamp, fSampleRate[ddl], stroff, samp);
375 str = strbase + stroff;
377 AliDebug(20,Form("0x%04x/0x%03x/%04d=%4d maps to strip %3d sample %d "
378 "(pre: %d, min: %d, max: %d, rate: %d)",
379 ddl, hwa, t, adc, str, samp, fPreSamp,
380 fMinStrip, fMaxStrip, fSampleRate[ddl]));
382 AliDebug(10,Form("Got presamples at timebin %d", i));
386 // VA1 Local strip number
387 Short_t lstrip = (t - fPreSamp) / fSampleRate[ddl] + fMinStrip;
389 AliDebug(15,Form("Checking if strip %d (%d) in range [%d,%d]",
390 lstrip, str, fMinStrip, fMaxStrip));
391 if (lstrip < fMinStrip || lstrip > fMaxStrip) {
392 AliDebug(10,Form("Strip %03d-%d (%d,%d) from t=%d out of range (%3d->%3d)",
393 str, samp, lstrip, stroff, t, fMinStrip, fMaxStrip));
396 // Possibly do pedestal subtraction of signal
398 AliWarning(Form("ADC value out of range: %4d", adc));
402 //____________________________________________________________________
404 AliFMDRawReader::NextSample(UShort_t& det, Char_t& rng, UShort_t& sec,
405 UShort_t& str, UShort_t& sam, UShort_t& rat,
406 Short_t& adc, Bool_t& zs, UShort_t& fac)
408 // Scan current event for next signal. It returns kFALSE when
409 // there's no more data in the event.
411 // Note, that this member function is in principle very fast, but
412 // contains less error checking. In particular, channels that have
413 // bad bunches cannot be checked here. Seeing a bad bunch will only
414 // skip the remainder of the channel and not reset the already read
415 // digits. This is potentially dangerous.
418 // det On return, contain the detector number
419 // rng On return, contain the ring identifier
420 // sec On return, contain the sector number
421 // str On return, contain the strip number
422 // sam On return, contain the sample number
423 // rat On return, contain the sample rate
424 // adc On return, contain the ADC counts
425 // zs On return, contain the zero-supp. flag
426 // fac On return, contain the zero-supp. noise factor
430 // -1 Read sample belongs to a bad bunch
431 // >0 Good status - contains bit mask of values
436 static AliAltroRawStreamV3 stream(fReader); // = 0;
437 static Int_t ddl = -1;
438 static UShort_t tdet = 0;
439 static Char_t trng = '\0';
440 static UShort_t tsec = 0;
441 static Short_t tstr = 0;
442 static Short_t bstr = -1;
443 static UShort_t tsam = 0;
444 static UInt_t trate = 0;
445 static Int_t hwaddr = -1;
446 static UShort_t start = 0;
447 static UShort_t length = 0;
448 static Short_t t = -1;
451 if (stream.GetDDLNumber() < 0) {
453 fReader->Select("FMD");
455 stream.SelectRawData("FMD");
456 stream.SetCheckAltroPayload(false);
457 for (Int_t j = 0; j < kNDDL; j++) fNErrors[j] = 0;
472 AliDebug(15,Form("t=%4d, start=%4d, length=%4d", t, start, length));
473 if (t < start - length + 1) {
474 AliDebug(10,Form("Time t=%d < start-length+1=%d-%d+1 (%3d/0x%03x)",
475 t, start, length, ddl, hwaddr));
476 if (hwaddr > 0xFFF ||
478 !stream.NextBunch()) {
479 if (AliDebugLevel() >= 10 && hwaddr > 0xFFF) {
480 AliDebug(10,"Last channel read was marked bad");
482 if (AliDebugLevel() >= 10 && hwaddr < 0) {
483 AliDebug(10,"No more channels");
485 AliDebug(10,"No next bunch, or first entry");
486 if (ddl < 0 || !stream.NextChannel()) {
487 if (AliDebugLevel() >= 10 && ddl < 0) {
488 AliDebug(10,"No DDL");
490 AliDebug(10,"No next channel, or first entry");
491 if (!stream.NextDDL()) {
492 AliDebug(10,"No more DDLs");
496 ddl = NewDDL(stream, tdet);
497 AliDebug(5,Form("New DDL: %d (%d)", ddl, tdet));
501 hwaddr = NewChannel(stream, tdet, trng, tsec, bstr);
502 if (hwaddr > 0xFFF) fNErrors[ddl] += 1;
503 AliDebug(5,Form("New Channel: %3d/0x%03x", ddl, hwaddr));
508 if (!NewBunch(stream, start, length)) {
509 // AliWarning(Form("Bad bunch in %3d/0x%03x read - "
510 // "should progress to next channel "
511 // "(t=%4d,start=%4d,length=%4d)",
512 // ddl, hwaddr, t,start, length));
513 hwaddr = 0xFFFF; // Bad channel
516 AliDebug(5, Form("New bunch in %3d/0x%03x: start=0x%03x, length=%4d",
517 ddl, hwaddr, start, length));
521 AliDebug(10,Form("Got new bunch FMD%d%c[%2d], bunch @ %d, length=%d",
522 tdet, trng, tsec, start, length));
524 Int_t tadc = NewSample(stream, i, t, tsec, bstr, tstr, tsam);
525 AliDebug(10,Form("New sample FMD%d%c[%2d,%3d]-%d = 0x%03x",
526 tdet, trng, tsec, tstr, tsam, tadc));
535 rat = fSampleRate[ddl];
536 zs = fZeroSuppress[ddl];
537 fac = fNoiseFactor[ddl];
540 AliDebug(10,Form("Returning FMD%d%c[%2d,%3d]-%d = 0x%03x (%d,%d,%d)",
541 det, rng, sec, str, sam, adc, rat, zs, fac));
547 AliDebug(5,Form("Returning 0x%02x", ret));
552 //____________________________________________________________________
554 AliFMDRawReader::NextSignal(UShort_t& det, Char_t& rng,
555 UShort_t& sec, UShort_t& str,
556 Short_t& adc, Bool_t& zs,
560 // Get the next signal
563 // det On return, the detector
564 // rng On return, the ring
565 // sec On return, the sector
566 // str On return, the strip
567 // adc On return, the ADC value
568 // zs On return, whether zero-supp. is enabled
569 // fac On return, the usd noise factor
572 // true if valid data is returned
577 if ((ret = NextSample(det, rng, sec, str, samp, rate, adc, zs, fac)) <= 0)
580 Bool_t take = SelectSample(samp, rate);
587 //____________________________________________________________________
589 AliFMDRawReader::SelectSample(UShort_t samp, UShort_t rate)
591 // Check if the passed sample is the one we need
592 Bool_t take = kFALSE;
594 case 1: take = kTRUE; break;
595 case 2: if (samp == 1) take = kTRUE; break;
596 case 3: if (samp == 1) take = kTRUE; break;
597 case 4: if (samp == 2) take = kTRUE; break;
598 default: if (samp == rate-2) take = kTRUE; break;
604 //____________________________________________________________________
606 AliFMDRawReader::ReadAdcs(TClonesArray* array)
608 // Read ADC values from raw input into passed TClonesArray of AliFMDDigit
610 AliDebug(3,Form("Reading ADC values into a TClonesArray"));
612 // Read raw data into the digits array, using AliFMDAltroReader.
614 AliError("No TClonesArray passed");
617 const UShort_t kUShortMax = (1 << 16) - 1;
618 fSeen.Reset(kUShortMax);
619 for (Int_t ddl = 0; ddl < kNDDL; ddl++) fNErrors[ddl] = 0;
621 AliAltroRawStreamV3 input(fReader);
623 input.SetCheckAltroPayload(false);
624 input.SelectRawData("FMD");
626 // Loop over input RORCs
627 while (input.NextDDL()) {
629 Int_t ddl = NewDDL(input, det);
633 while (input.NextChannel()) {
634 // Get the hardware address, and map that to detector coordinates
638 Int_t hwaddr = NewChannel(input, det, ring, sec, strbase);
639 if (hwaddr < 0) break;
640 if (hwaddr > 0xFFF) continue;
642 UShort_t start = 0x3FF;
643 Bool_t errors = false;
647 while (input.NextBunch()) {
648 // Get Lenght of bunch, and pointer to the data
649 const UShort_t* data = input.GetSignals();
651 if (!NewBunch(input, start, length)) {
657 // Loop over the data and store it.
658 for (Int_t i = 0; i < length; i++) {
663 Int_t adc = NewSample(input, i, t, sec, strbase, str, samp);
664 if (adc < 0) continue;
665 UShort_t counts = adc;
667 AliDebug(10, Form("FMD%d%c[%02d,%03d]-%d: %4d",
668 det, ring, sec, str, samp, counts));
669 // Check the cache of indicies
670 Int_t idx = fSeen(det, ring, sec, str);
671 AliFMDDigit* digit = 0;
672 if (idx == kUShortMax) {
673 // We haven't seen this strip yet.
674 fSeen(det, ring, sec, str) = idx = array->GetEntriesFast();
675 AliDebug(7,Form("making digit @ %5d for FMD%d%c[%2d,%3d]-%d "
676 "from %3d/0x%03x/%4d",
677 idx, det, ring, sec, str, samp, ddl, hwaddr, t));
678 digit = new ((*array)[idx]) AliFMDDigit(det, ring, sec, str);
679 digit->SetDefaultCounts(fSampleRate[ddl]);
682 digit = static_cast<AliFMDDigit*>(array->At(idx));
684 if (first < 0) first = idx;
686 AliDebug(10, Form("Setting FMD%d%c[%2d,%3d]-%d from timebin "
687 "%4d=%4d (%4d)", det, ring, sec, str, samp, t,
689 digit->SetCount(samp, counts);
693 AliWarning(Form("Channel %3d/0x%03x contain errors, "
694 "resetting index %d to %d", ddl, hwaddr, first, last));
696 for (Int_t i = first; i <= last; i++) {
697 AliFMDDigit* digit = static_cast<AliFMDDigit*>(array->At(i));
698 for (Int_t j = 0; j < fSampleRate[ddl]; j++) {
699 AliDebug(10,Form("Resetting strip %s=%d",
700 digit->GetName(),digit->Counts()));
701 digit->SetCount(j, kBadSignal);
706 // if (errors && (AliDebugLevel() > 0)) input.HexDumpChannel();
711 //____________________________________________________________________
713 AliFMDRawReader::ReadAdcs(AliFMDUShortMap& map)
715 // Read ADC values from raw input into passed TClonesArray of AliFMDDigit
717 AliDebug(3,Form("Reading ADC values into a map"));
719 // const UShort_t kUShortMax = (1 << 16) - 1;
720 for (Int_t ddl = 0; ddl < kNDDL; ddl++) fNErrors[ddl] = 0;
722 AliAltroRawStreamV3 input(fReader);
724 input.SetCheckAltroPayload(false);
725 input.SelectRawData("FMD");
727 // Loop over input RORCs
728 while (input.NextDDL()) {
730 Int_t ddl = NewDDL(input, det);
734 while (input.NextChannel()) {
735 // Get the hardware address, and map that to detector coordinates
739 Int_t hwaddr = NewChannel(input, det, ring, sec, strbase);
740 if (hwaddr < 0) break;
741 if (hwaddr > 0xFFF) continue;
743 UShort_t start = 0x3FF;
744 Bool_t errors = false;
748 while (input.NextBunch()) {
749 // Get Lenght of bunch, and pointer to the data
750 // const UShort_t* data = input.GetSignals();
752 if (!NewBunch(input, start, length)) {
758 // Loop over the data and store it.
759 for (Int_t i = 0; i < length; i++) {
764 Int_t adc = NewSample(input, i, t, sec, strbase, str, samp);
765 if (adc < 0) continue;
766 UShort_t counts = adc;
768 AliDebug(10, Form("FMD%d%c[%02d,%03d]-%d: %4d",
769 det, ring, sec, str, samp, counts));
770 if (SelectSample(samp, fSampleRate[ddl]))
771 map(det,ring,sec,str) = counts;
772 if (first < 0) first = str;
777 AliWarning(Form("Channel %3d/0x%03x contain errors, "
778 "resetting strips %d to %d", ddl, hwaddr, first, last));
780 Int_t ds = first <= last ? 1 : -1;
781 for (Int_t i = first; i != last+ds; i += ds) {
782 AliDebug(10, Form("Resetting strip FMD%d%c[%02d,%03d]=%d",
783 det,ring,sec,i,map(det,ring,sec,i)));
784 map(det,ring,sec,i) = kBadSignal;
793 //____________________________________________________________________
794 Bool_t AliFMDRawReader::ReadSODevent(AliFMDCalibSampleRate* sampleRate,
795 AliFMDCalibStripRange* stripRange,
797 TArrayS &pulseLength,
801 // Read SOD event into passed objects.
804 // samplerate The sample rate object to fill
805 // striprange The strip range object to fill
806 // pulseSize The pulse size object to fill
807 // pulseLength The pulse length (in events) object to fill
810 // @c true on success
812 AliDebug(0,Form("Start of SOD/EOD"));
814 UInt_t shift_clk[18];
815 UInt_t sample_clk[18];
816 UInt_t strip_low[18];
817 UInt_t strip_high[18];
818 UInt_t pulse_size[18];
819 UInt_t pulse_length[18];
820 for (size_t i = 0; i < 18; i++) {
828 AliFMDParameters* param = AliFMDParameters::Instance();
829 AliFMDAltroMapping* map = param->GetAltroMap();
831 AliAltroRawStreamV3 streamer(fReader);
833 streamer.SelectRawData("FMD");
834 //fReader->GetDDLID();
835 while (streamer.NextDDL()) {
836 Int_t ddl = streamer.GetDDLNumber();
837 Int_t detID = fReader->GetDetectorID();
838 if (detectors) detectors[map->DDL2Detector(ddl)-1] = kTRUE;
839 AliDebug(0,Form(" From reader: DDL number is %d , det ID is %d",ddl,detID));
841 ULong_t nPayloadWords = streamer.GetRCUPayloadSizeInSOD();
842 UChar_t* payloadData = streamer.GetRCUPayloadInSOD();
843 UInt_t* payloadWords = reinterpret_cast<UInt_t*>(payloadData);
844 //UInt_t* payloadWords = streamer.GetRCUPayloadInSOD();
846 //std::cout<<nPayloadWords<<" "<<ddl<<std::endl;
847 for (ULong_t i = 1; i <= nPayloadWords ; i++, payloadWords++) {
848 UInt_t payloadWord = *payloadWords; // Get32bitWord(i);
850 //std::cout<<i<<Form(" word: 0x%x",payloadWord)<<std::endl;
851 // address is only 24 bit
852 UInt_t address = (0xffffff & payloadWord);
853 UInt_t type = ((address >> 21) & 0xf);
854 UInt_t error = ((address >> 20) & 0x1);
855 UInt_t bcast = ((address >> 18) & 0x1);
856 UInt_t bc_not_altro = ((address >> 17) & 0x1);
857 UInt_t board = ((address >> 12) & 0x1f);
858 UInt_t instruction = 0;
862 instruction = address & 0xfff;
864 chip = ((address >> 9) & 0x7);
865 channel = ((address >> 5) & 0x5);
866 instruction = (address & 0x1f);
869 Bool_t readDataWord = kFALSE;
871 case 0x0: // Fec read
872 readDataWord = kTRUE;
874 case 0x2: // Fec write
881 case 0x6: // End sequence
883 i = nPayloadWords + 1;
889 //Don't read unless we have a FEC_RD
890 if(!readDataWord) continue;
892 UInt_t dataWord = *payloadWords;//Get32bitWord(i);
893 UInt_t data = (0xFFFFF & dataWord) ;
894 //UInt_t data = (0xFFFF & dataWord) ;
897 AliWarning(Form("error bit detected at Word 0x%06x; "
898 "error % d, type %d, bc_not_altro %d, "
899 "bcast %d, board 0x%02x, chip 0x%x, "
900 "channel 0x%02x, instruction 0x%03x",
901 address, error, type, bc_not_altro,
902 bcast,board,chip,channel,instruction));
908 switch(instruction) {
910 case 0x01: break; // First ADC T
911 case 0x02: break; // I 3.3 V
912 case 0x03: break; // I 2.5 V altro digital
913 case 0x04: break; // I 2.5 V altro analog
914 case 0x05: break; // I 2.5 V VA
915 case 0x06: break; // First ADC T
916 case 0x07: break; // I 3.3 V
917 case 0x08: break; // I 2.5 V altro digital
918 case 0x09: break; // I 2.5 V altro analog
919 case 0x0A: break; // I 2.5 V VA
920 case 0x2D: break; // Second ADC T
921 case 0x2E: break; // I 1.5 V VA
922 case 0x2F: break; // I -2.0 V
923 case 0x30: break; // I -2.0 V VA
924 case 0x31: break; // 2.5 V Digital driver
925 case 0x32: break; // Second ADC T
926 case 0x33: break; // I 1.5 V VA
927 case 0x34: break; // I -2.0 V
928 case 0x35: break; // I -2.0 V VA
929 case 0x36: break; // 2.5 V Digital driver
930 case 0x37: break; // Third ADC T
931 case 0x38: break; // Temperature sens. 1
932 case 0x39: break; // Temperature sens. 2
933 case 0x3A: break; // U 2.5 altro digital (m)
934 case 0x3B: break; // U 2.5 altro analog (m)
935 case 0x3C: break; // Third ADC T
936 case 0x3D: break; // Temperature sens. 1
937 case 0x3E: break; // Temperature sens. 2
938 case 0x3F: break; // U 2.5 altro digital (m)
939 case 0x40: break; // U 2.5 altro analog (m)
940 case 0x41: break; // Forth ADC T
941 case 0x42: break; // U 2.5 VA (m)
942 case 0x43: break; // U 1.5 VA (m)
943 case 0x44: break; // U -2.0 VA (m)
944 case 0x45: break; // U -2.0 (m)
945 case 0x46: break; // Forth ADC T
946 case 0x47: break; // U 2.5 VA (m)
947 case 0x48: break; // U 1.5 VA (m)
948 case 0x49: break; // U -2.0 VA (m)
949 case 0x4A: break; // U -2.0 (m)
951 case 0x0B: break; // L1 trigger CouNTer
952 case 0x0C: break; // L2 trigger CouNTer
953 case 0x0D: break; // Sampling CLK CouNTer
954 case 0x0E: break; // DSTB CouNTer
956 case 0x0F: break; // Test mode word
957 case 0x10: break; // Undersampling ratio.
958 // Configuration and status
959 case 0x11: break; // Config/Status Register 0
960 case 0x12: break; // Config/Status Register 1
961 case 0x13: break; // Config/Status Register 2
962 case 0x14: break; // Config/Status Register 3
963 case 0x15: break; // Free
965 case 0x16: break; // Latch L1, L2, SCLK Counters
966 case 0x17: break; // Clear counters
967 case 0x18: break; // Clear CSR1
968 case 0x19: break; // rstb ALTROs
969 case 0x1A: break; // rstb BC
970 case 0x1B: break; // Start conversion
971 case 0x1C: break; // Scan event length
972 case 0x1D: break; // Read event length
973 case 0x1E: break; // Start test mode
974 case 0x1F: break; // Read acquisition memory
976 case 0x20: break; // FMDD status
977 case 0x21: break; // L0 counters
978 case 0x22: break; // FMD: Wait to hold
979 case 0x23: break; // FMD: L1 timeout
980 case 0x24: break; // FMD: L2 timeout
981 case 0x25: // FMD: Shift clk
982 shift_clk[board] = ((data >> 8 ) & 0xFF);
983 AliDebug(30,Form("Read shift_clk=%d for board 0x%02x",
984 shift_clk[board], board));
986 case 0x26: // FMD: Strips
987 strip_low[board] = ((data >> 0 ) & 0xFF);
988 strip_high[board] = ((data >> 8 ) & 0xFF);
990 case 0x27: // FMD: Cal pulse
991 pulse_size[board] = ((data >> 8 ) & 0xFF);
993 case 0x28: break; // FMD: Shape bias
994 case 0x29: break; // FMD: Shape ref
995 case 0x2A: break; // FMD: Preamp ref
996 case 0x2B: // FMD: Sample clk
997 sample_clk[board] = ((data >> 8 ) & 0xFF);
998 AliDebug(30,Form("Read sample_clk=%d for board 0x%02x",
999 sample_clk[board], board));
1001 case 0x2C: break; // FMD: Commands
1002 case 0x4B: // FMD: Cal events
1003 pulse_length[board] = ((data >> 0 ) & 0xFF);
1008 AliDebug(50,Form("instruction 0x%x, dataword 0x%x",
1009 instruction,dataWord));
1010 } // End of loop over Result memory event
1013 UShort_t sector = 0;
1017 const UInt_t boards[4] = {0,1,16,17};
1018 for(Int_t i=0;i<4;i++) {
1019 if(ddl==0 && (i==1 || i==3)) continue;
1021 UInt_t chip =0, channel=0;
1022 det = map->DDL2Detector(ddl);
1023 map->Channel2StripBase(boards[i], chip, channel, ring, sector, strip);
1025 UInt_t samplerate = 0;
1027 if(sample_clk[boards[i]] == 0) {
1029 Int_t sample1 = sample_clk[boards[0]];
1030 Int_t sample2 = sample_clk[boards[2]];
1031 if(sample1) sample_clk[boards[i]] = sample1;
1032 else sample_clk[boards[i]] = sample2;
1035 Int_t sample1 = sample_clk[boards[0]];
1036 Int_t sample2 = sample_clk[boards[1]];
1037 Int_t sample3 = sample_clk[boards[2]];
1038 Int_t sample4 = sample_clk[boards[3]];
1039 Int_t agreement = 0;
1040 if(sample1 == sample2) agreement++;
1041 if(sample1 == sample3) agreement++;
1042 if(sample1 == sample4) agreement++;
1043 if(sample2 == sample3) agreement++;
1044 if(sample2 == sample4) agreement++;
1045 if(sample3 == sample4) agreement++;
1050 if(agreement == 3) {
1051 sample_clk[boards[i]] = sample_clk[boards[idx]];
1052 shift_clk[boards[i]] = shift_clk[boards[idx]];
1053 strip_low[boards[i]] = strip_low[boards[idx]];
1054 strip_high[boards[i]] = strip_high[boards[idx]];
1055 pulse_length[boards[i]] = pulse_length[boards[idx]];
1056 pulse_size[boards[i]] = pulse_size[boards[idx]];
1057 AliDebug(3,Form("Vote taken for ddl %d, board 0x%x",
1064 if(sample_clk[boards[i]])
1065 samplerate = shift_clk[boards[i]]/sample_clk[boards[i]];
1066 AliDebug(10,Form("Sample rate for board 0x%02x is %d",
1067 boards[i], samplerate));
1068 sampleRate->Set(det,ring,sector,0,samplerate);
1069 stripRange->Set(det,ring,sector,0,
1070 strip_low[boards[i]],strip_high[boards[i]]);
1072 AliDebug(20,Form("det %d, ring %c, ",det,ring));
1073 pulseLength.AddAt(pulse_length[boards[i]],
1074 GetHalfringIndex(det,ring,boards[i]/16));
1075 pulseSize.AddAt(pulse_size[boards[i]],
1076 GetHalfringIndex(det,ring,boards[i]/16));
1080 AliDebug(20,Form(": Board: 0x%02x\n"
1081 "\tstrip_low %3d, strip_high %3d\n"
1082 "\tshift_clk %3d, sample_clk %3d\n"
1083 "\tpulse_size %3d, pulse_length %3d",
1085 strip_low[boards[i]], strip_high[boards[i]],
1086 shift_clk[boards[i]], sample_clk[boards[i]],
1087 pulse_size[boards[i]],pulse_length[boards[i]]));
1092 AliFMDParameters::Instance()->SetSampleRate(sampleRate);
1093 AliFMDParameters::Instance()->SetStripRange(stripRange);
1095 AliDebug(0,Form("End of SOD/EOD"));
1099 //____________________________________________________________________
1101 UInt_t AliFMDRawReader::Get32bitWord(Int_t idx)
1103 // This method returns the 32 bit word at a given
1104 // position inside the raw data payload.
1105 // The 'index' points to the beginning of the next word.
1106 // The method is supposed to be endian (platform)
1109 AliFatal("Raw data paylod buffer is not yet initialized !");
1112 Int_t index = 4*idx;
1115 // fRawReader->AddFatalErrorLog(k32bitWordReadErr,Form("pos = %d",index));
1117 AliWarning(Form("Invalid raw data payload position (%d) !",index));
1122 word = fData[--index] << 24;
1123 word |= fData[--index] << 16;
1124 word |= fData[--index] << 8;
1125 word |= fData[--index] << 0 ;
1129 //_____________________________________________________________________
1130 Int_t AliFMDRawReader::GetHalfringIndex(UShort_t det, Char_t ring,
1131 UShort_t board) const
1134 // Get short index for a given half-ring
1137 // det Detector number
1138 // ring Ring identifer
1139 // board Board number
1144 UShort_t iring = (ring == 'I' ? 1 : 0);
1146 Int_t index = (((det-1) << 2) | (iring << 1) | (board << 0));
1152 //____________________________________________________________________