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;
101 //____________________________________________________________________
103 AliFMDRawReader::Exec(Option_t*)
106 TClonesArray* array = new TClonesArray("AliFMDDigit");
111 fTree->Branch("FMD", &array);
115 Int_t nWrite = fTree->Fill();
116 AliDebug(1,Form("Got a grand total of %d digits, wrote %d bytes to tree",
117 array->GetEntriesFast(), nWrite));
121 //____________________________________________________________________
123 AliFMDRawReader::NewDDL(AliAltroRawStreamV3& input, UShort_t& det)
125 // Process a new DDL. Sets the internal data members fZeroSuppress,
126 // fSampleRate, and fNoiseFactor based on information in the RCU trailer.
129 // input Input stream
130 // det On return, the detector number
133 // negative value in case of problems, the DDL number otherwise
135 // Get the DDL number
136 UInt_t ddl = input.GetDDLNumber();
137 AliDebug(2,Form("DDL number %d", ddl));
139 // Note, previously, the ALTROCFG1 register was interpreted as
141 // Bits Value Description
142 // 0- 3 0/1 1st Baseline filter, mode
143 // 4- 5 Over-1 2nd baseline filter, # of pre-samples
144 // 6- 9 factor 2nd baseline filter, # of post-samples
145 // 10- 0 2nd baseline filter, enable
146 // 11-12 00 Zero suppression, glitch filter mode
147 // 13-15 001 Zero suppression, # of post samples
148 // 16-17 01 Zero suppression, # of pre samples
149 // 18 0/1 Zero suppression, enable
151 // The interpretation used in AliAltroRawStreamerV3 - which
152 // corresponds directly to ALTRO DPCFG register - is
154 // Bits Value Description
155 // 0- 3 0/1 1st Baseline filter, mode
156 // 4 0 Polarity (if '1', then "1's inverse")
157 // 5- 6 01 Zero suppression, # of pre samples
158 // 7-10 0001 Zero suppression, # of post samples
159 // 11 0 2nd baseline filter, enable
160 // 12-13 00 Zero suppression, glitch filter mode
161 // 14-16 factor 2nd baseline filter, # of post-samples
162 // 17-18 01 2nd baseline filter, # of pre-samples
163 // 19 0/1 Zero suppression, enable
165 // Writing 'x' for variable values, that means we have the
166 // following patterns for the 2 cases
168 // bit # 20 16 12 8 4 0
169 // old |0x01|0010|00xx|xxxx|xxxx|
170 // new |x01x|xx00|0000|1010|xxxx|
172 // That means that we can check if bits 10-13 are '1000' or
173 // '0000', which will tell us if the value was written with the
174 // new or the old interpretation. That is, we can check that
176 // if (((altrocfg1 >> 10) & 0x8) == 0x8) {
177 // // old interpretation
180 // // New interpretation
183 // That means, that we should never
185 // - change the # of zero suppression post samples
186 // - Turn on 2nd baseline correction
187 // - Change the zero-suppression glitch filter mode
189 // This change as introduced in version 1.2 of Rcu++
191 UInt_t cfg1 = input.GetAltroCFG1();
192 if (((cfg1 >> 10) & 0x8) == 0x8) {
193 UInt_t cfg2 = input.GetAltroCFG2();
194 AliDebug(3,Form("We have data from older MiniConf 0x%x cfg2=0x%08x",
195 ((cfg1 >> 10) & 0x8), cfg2));
196 fZeroSuppress[ddl] = (cfg1 >> 0) & 0x1;
197 fNoiseFactor[ddl] = (cfg1 >> 6) & 0xF;
198 fSampleRate[ddl] = (cfg2 >> 20) & 0xF;
201 AliDebug(3,Form("We have data from newer MiniConf 0x%x",
202 ((cfg1 >> 10) & 0x8)));
203 fZeroSuppress[ddl] = input.GetZeroSupp();
204 // WARNING: We store the noise factor in the 2nd baseline
205 // filters excluded post samples, since we'll never use that
207 fNoiseFactor[ddl] = input.GetNPostsamples();
208 // WARNING: We store the sample rate in the number of pre-trigger
209 // samples, since we'll never use that mode.
210 fSampleRate[ddl] = input.GetNPretriggerSamples();
213 AliDebug(10,Form("Phase of DDL=%d is %g (%d)", ddl, input.GetL1Phase(),
214 input.GetAltroCFG2() & 0x1F));
215 fL1Phase[ddl] = input.GetAltroCFG2() & 0x1F; // input.GetL1Phase();
216 AliDebug(3,Form("RCU @ DDL %d zero suppression: %s",
217 ddl, (fZeroSuppress[ddl] ? "yes" : "no")));
218 AliDebug(3,Form("RCU @ DDL %d noise factor: %d", ddl,fNoiseFactor[ddl]));
219 AliDebug(3,Form("RCU @ DDL %d sample rate: %d", ddl,fSampleRate[ddl]));
223 Int_t nChAddrMismatch = input.GetNChAddrMismatch();
224 Int_t nChLenMismatch = input.GetNChLengthMismatch();
225 if (nChAddrMismatch != 0)
226 AliWarning(Form("Got %d channels with address mis-matches for 0x%03x",
227 nChAddrMismatch, ddl));
228 if (nChLenMismatch != 0)
229 AliWarning(Form("Got %d channels with length mis-matches for 0x%03x",
230 nChLenMismatch, ddl));
232 // Map DDL number to the detector number
233 AliFMDParameters* pars = AliFMDParameters::Instance();
234 AliFMDAltroMapping* map = pars->GetAltroMap();
235 if (map->DDL2Detector(ddl) < 0) return -1;
236 det = map->DDL2Detector(ddl);
238 if (AliLog::GetDebugLevel("FMD", 0) > 5)
239 input.PrintRCUTrailer();
243 //____________________________________________________________________
245 AliFMDRawReader::NewChannel(const AliAltroRawStreamV3& input, UShort_t det,
246 Char_t& ring, UShort_t& sec, Short_t& strbase)
248 // Processs a new channel. Sets the internal data members
249 // fMinStrip, fMaxStrip, and fPreSamp.
252 // input Input stream
253 // ring On return, the ring identifier
254 // sec On return, the sector number
255 // strbase On return, the strip base
258 // negative value in case of problems, hardware address otherwise
260 // Get the hardware address, and map that to detector coordinates
261 UShort_t board, chip, channel;
262 Int_t ddl = input.GetDDLNumber();
263 Int_t hwaddr = input.GetHWAddress();
264 if (input.IsChannelBad()) {
265 const char* msg = Form("Ignoring channel %03d/0x%03x with errors",
268 if (AliDebugLevel() > 10) input.HexDumpChannel();
269 fReader->AddMinorErrorLog(AliAltroRawStreamV3::kAltroPayloadErr,msg);
274 AliFMDParameters* pars = AliFMDParameters::Instance();
275 AliFMDAltroMapping* map = pars->GetAltroMap();
276 // Map to hardware stuff
277 map->ChannelAddress(hwaddr, board, chip, channel);
278 // Then try to map to detector address
279 if (!map->Channel2StripBase(board, chip, channel, ring, sec, strbase)) {
280 AliError(Form("Failed to get detector id from DDL %d, "
281 "hardware address 0x%03x", ddl, hwaddr));
284 AliDebug(4,Form("Board: 0x%02x, Altro: 0x%x, Channel: 0x%x",
285 board, chip, channel));
287 // Get the 'conditions'
288 fMinStrip = pars->GetMinStrip(det, ring, sec, strbase);
289 fMaxStrip = pars->GetMaxStrip(det, ring, sec, strbase);
290 fPreSamp = pars->GetPreSamples(det, ring, sec, strbase);
291 if (fSampleRate[ddl] == 0) {
292 AliDebug(3,Form("Get sample rate for RCU @ DDL %d from OCDB", ddl));
293 fSampleRate[ddl] = pars->GetSampleRate(det, ring, sec, strbase);
295 AliDebug(3,Form("RCU @ DDL %d sample rate: %d", ddl,fSampleRate[ddl]));
300 //____________________________________________________________________
302 AliFMDRawReader::NewBunch(const AliAltroRawStreamV3& input,
303 UShort_t& start, UShort_t& length)
306 // Do some checks on the bunch data
308 Int_t ddl = input.GetDDLNumber();
309 Int_t hwaddr = input.GetHWAddress();
310 UShort_t nSamples = input.GetNSamplesPerCh() + fPreSamp;
311 UShort_t tstart = input.GetStartTimeBin();
312 length = input.GetBunchLength();
314 if (tstart >= nSamples) {
315 const char* msg = Form("Bunch in %03d/0x%03x has an start time greater "
316 "than number of samples: 0x%x >= 0x%x",
317 ddl, hwaddr, tstart, nSamples);
319 if (AliDebugLevel() > 10) input.HexDumpChannel();
320 fReader->AddMinorErrorLog(AliAltroRawStreamV3::kAltroPayloadErr,msg);
324 if ((int(tstart) - length + 1) < 0) {
325 const char* msg = Form("Bunch in %03d/0x%03x has an invalid length and "
326 "start time: 0x%x,0x%x (%d-%d+1=%d<0)",
327 ddl, hwaddr, length, tstart, tstart, length,
328 int(tstart)-length+1);
330 if (AliDebugLevel() > 10) input.HexDumpChannel();
331 fReader->AddMinorErrorLog(AliAltroRawStreamV3::kAltroPayloadErr,msg);
335 if (tstart >= start) {
336 const char* msg = Form("Bunch in %03d/0x%03x has early start time: "
337 "0x%x >= 0x%x", ddl, hwaddr, tstart, start);
339 if (AliDebugLevel() > 10) input.HexDumpChannel();
340 fReader->AddMinorErrorLog(AliAltroRawStreamV3::kAltroPayloadErr,msg);
348 //____________________________________________________________________
350 AliFMDRawReader::NewSample(const AliAltroRawStreamV3& input,
351 Int_t i, UShort_t t, UShort_t sec,
352 UShort_t strbase, Short_t& str, UShort_t& samp)
354 // Process a new timebin
357 // input Input stream
358 // i Index into bunch data
360 // strbase Base of strip numbers for this channel
361 // str On return, the strip number
362 // samp On return, the sample number
365 // negative value in case of problems, ADC value otherwise
366 if (t < fPreSamp) return -1;
368 Int_t ddl = input.GetDDLNumber();
369 Int_t hwa = input.GetHWAddress();
370 const UShort_t* data = input.GetSignals();
371 Short_t adc = data[i];
372 AliDebug(10,Form("0x%04x/0x%03x/%04d %4d", ddl, hwa, t, adc));
374 AliFMDParameters* pars = AliFMDParameters::Instance();
375 AliFMDAltroMapping* map = pars->GetAltroMap();
379 map->Timebin2Strip(sec, t, fPreSamp, fSampleRate[ddl], stroff, samp);
380 str = strbase + stroff;
382 AliDebug(20,Form("0x%04x/0x%03x/%04d=%4d maps to strip %3d sample %d "
383 "(pre: %d, min: %d, max: %d, rate: %d)",
384 ddl, hwa, t, adc, str, samp, fPreSamp,
385 fMinStrip, fMaxStrip, fSampleRate[ddl]));
387 AliDebug(10,Form("Got presamples at timebin %d", i));
391 // VA1 Local strip number
392 Short_t lstrip = (t - fPreSamp) / fSampleRate[ddl] + fMinStrip;
394 AliDebug(15,Form("Checking if strip %d (%d) in range [%d,%d]",
395 lstrip, str, fMinStrip, fMaxStrip));
396 if (lstrip < fMinStrip || lstrip > fMaxStrip) {
397 AliDebug(10,Form("Strip %03d-%d (%d,%d) from t=%d out of range (%3d->%3d)",
398 str, samp, lstrip, stroff, t, fMinStrip, fMaxStrip));
401 // Possibly do pedestal subtraction of signal
403 AliWarning(Form("ADC value out of range: %4d", adc));
407 //____________________________________________________________________
409 AliFMDRawReader::NextSample(UShort_t& det, Char_t& rng, UShort_t& sec,
410 UShort_t& str, UShort_t& sam, UShort_t& rat,
411 Short_t& adc, Bool_t& zs, UShort_t& fac)
413 // Scan current event for next signal. It returns kFALSE when
414 // there's no more data in the event.
416 // Note, that this member function is in principle very fast, but
417 // contains less error checking. In particular, channels that have
418 // bad bunches cannot be checked here. Seeing a bad bunch will only
419 // skip the remainder of the channel and not reset the already read
420 // digits. This is potentially dangerous.
423 // det On return, contain the detector number
424 // rng On return, contain the ring identifier
425 // sec On return, contain the sector number
426 // str On return, contain the strip number
427 // sam On return, contain the sample number
428 // rat On return, contain the sample rate
429 // adc On return, contain the ADC counts
430 // zs On return, contain the zero-supp. flag
431 // fac On return, contain the zero-supp. noise factor
435 // -1 Read sample belongs to a bad bunch
436 // >0 Good status - contains bit mask of values
441 static AliAltroRawStreamV3 stream(fReader); // = 0;
442 static Int_t ddl = -1;
443 static UShort_t tdet = 0;
444 static Char_t trng = '\0';
445 static UShort_t tsec = 0;
446 static Short_t tstr = 0;
447 static Short_t bstr = -1;
448 static UShort_t tsam = 0;
449 static UInt_t trate = 0;
450 static Int_t hwaddr = -1;
451 static UShort_t start = 0;
452 static UShort_t length = 0;
453 static Short_t t = -1;
456 if (stream.GetDDLNumber() < 0) {
458 fReader->Select("FMD");
460 stream.SelectRawData("FMD");
461 stream.SetCheckAltroPayload(false);
462 for (Int_t j = 0; j < kNDDL; j++) fNErrors[j] = 0;
477 AliDebug(15,Form("t=%4d, start=%4d, length=%4d", t, start, length));
478 if (t < start - length + 1) {
479 AliDebug(10,Form("Time t=%d < start-length+1=%d-%d+1 (%3d/0x%03x)",
480 t, start, length, ddl, hwaddr));
481 if (hwaddr > 0xFFF ||
483 !stream.NextBunch()) {
484 if (AliDebugLevel() >= 10 && hwaddr > 0xFFF) {
485 AliDebug(10,"Last channel read was marked bad");
487 if (AliDebugLevel() >= 10 && hwaddr < 0) {
488 AliDebug(10,"No more channels");
490 AliDebug(10,"No next bunch, or first entry");
491 if (ddl < 0 || !stream.NextChannel()) {
492 if (AliDebugLevel() >= 10 && ddl < 0) {
493 AliDebug(10,"No DDL");
495 AliDebug(10,"No next channel, or first entry");
496 if (!stream.NextDDL()) {
497 AliDebug(10,"No more DDLs");
501 ddl = NewDDL(stream, tdet);
502 AliDebug(5,Form("New DDL: %d (%d)", ddl, tdet));
506 hwaddr = NewChannel(stream, tdet, trng, tsec, bstr);
507 if (hwaddr > 0xFFF) fNErrors[ddl] += 1;
508 AliDebug(5,Form("New Channel: %3d/0x%03x", ddl, hwaddr));
513 if (!NewBunch(stream, start, length)) {
514 // AliWarning(Form("Bad bunch in %3d/0x%03x read - "
515 // "should progress to next channel "
516 // "(t=%4d,start=%4d,length=%4d)",
517 // ddl, hwaddr, t,start, length));
518 hwaddr = 0xFFFF; // Bad channel
521 AliDebug(5, Form("New bunch in %3d/0x%03x: start=0x%03x, length=%4d",
522 ddl, hwaddr, start, length));
526 AliDebug(10,Form("Got new bunch FMD%d%c[%2d], bunch @ %d, length=%d",
527 tdet, trng, tsec, start, length));
529 Int_t tadc = NewSample(stream, i, t, tsec, bstr, tstr, tsam);
530 AliDebug(10,Form("New sample FMD%d%c[%2d,%3d]-%d = 0x%03x",
531 tdet, trng, tsec, tstr, tsam, tadc));
540 rat = fSampleRate[ddl];
541 zs = fZeroSuppress[ddl];
542 fac = fNoiseFactor[ddl];
545 AliDebug(10,Form("Returning FMD%d%c[%2d,%3d]-%d = 0x%03x (%d,%d,%d)",
546 det, rng, sec, str, sam, adc, rat, zs, fac));
552 AliDebug(5,Form("Returning 0x%02x", ret));
557 //____________________________________________________________________
559 AliFMDRawReader::NextSignal(UShort_t& det, Char_t& rng,
560 UShort_t& sec, UShort_t& str,
561 Short_t& adc, Bool_t& zs,
565 // Get the next signal
568 // det On return, the detector
569 // rng On return, the ring
570 // sec On return, the sector
571 // str On return, the strip
572 // adc On return, the ADC value
573 // zs On return, whether zero-supp. is enabled
574 // fac On return, the usd noise factor
577 // true if valid data is returned
582 if ((ret = NextSample(det, rng, sec, str, samp, rate, adc, zs, fac)) <= 0)
585 Bool_t take = SelectSample(samp, rate);
592 //____________________________________________________________________
594 AliFMDRawReader::SelectSample(UShort_t samp, UShort_t rate)
596 // Check if the passed sample is the one we need
597 Bool_t take = kFALSE;
599 case 1: take = kTRUE; break;
600 case 2: if (samp == 1) take = kTRUE; break;
601 case 3: if (samp == 1) take = kTRUE; break;
602 case 4: if (samp == 2) take = kTRUE; break;
603 default: if (samp == rate-2) take = kTRUE; break;
609 //____________________________________________________________________
611 AliFMDRawReader::ReadAdcs(TClonesArray* array)
613 // Read ADC values from raw input into passed TClonesArray of AliFMDDigit
615 AliDebug(3,Form("Reading ADC values into a TClonesArray"));
617 // Read raw data into the digits array, using AliFMDAltroReader.
619 AliError("No TClonesArray passed");
622 const UShort_t kUShortMax = (1 << 16) - 1;
623 fSeen.Reset(kUShortMax);
624 for (Int_t ddl = 0; ddl < kNDDL; ddl++) fNErrors[ddl] = 0;
626 AliAltroRawStreamV3 input(fReader);
628 input.SetCheckAltroPayload(false);
629 input.SelectRawData("FMD");
631 // Loop over input RORCs
632 while (input.NextDDL()) {
634 Int_t ddl = NewDDL(input, det);
638 while (input.NextChannel()) {
639 // Get the hardware address, and map that to detector coordinates
643 Int_t hwaddr = NewChannel(input, det, ring, sec, strbase);
644 if (hwaddr < 0) break;
645 if (hwaddr > 0xFFF) continue;
647 UShort_t start = 0x3FF;
648 Bool_t errors = false;
652 while (input.NextBunch()) {
653 // Get Lenght of bunch, and pointer to the data
654 const UShort_t* data = input.GetSignals();
656 if (!NewBunch(input, start, length)) {
662 // Loop over the data and store it.
663 for (Int_t i = 0; i < length; i++) {
668 Int_t adc = NewSample(input, i, t, sec, strbase, str, samp);
669 if (adc < 0) continue;
670 UShort_t counts = adc;
672 AliDebug(10, Form("FMD%d%c[%02d,%03d]-%d: %4d",
673 det, ring, sec, str, samp, counts));
674 // Check the cache of indicies
675 Int_t idx = fSeen(det, ring, sec, str);
676 AliFMDDigit* digit = 0;
677 if (idx == kUShortMax) {
678 // We haven't seen this strip yet.
679 fSeen(det, ring, sec, str) = idx = array->GetEntriesFast();
680 AliDebug(7,Form("making digit @ %5d for FMD%d%c[%2d,%3d]-%d "
681 "from %3d/0x%03x/%4d",
682 idx, det, ring, sec, str, samp, ddl, hwaddr, t));
683 digit = new ((*array)[idx]) AliFMDDigit(det, ring, sec, str);
684 digit->SetDefaultCounts(fSampleRate[ddl]);
687 digit = static_cast<AliFMDDigit*>(array->At(idx));
689 if (first < 0) first = idx;
691 AliDebug(10, Form("Setting FMD%d%c[%2d,%3d]-%d from timebin "
692 "%4d=%4d (%4d)", det, ring, sec, str, samp, t,
694 digit->SetCount(samp, counts);
698 AliWarning(Form("Channel %3d/0x%03x contain errors, "
699 "resetting index %d to %d", ddl, hwaddr, first, last));
701 for (Int_t i = first; i <= last; i++) {
702 AliFMDDigit* digit = static_cast<AliFMDDigit*>(array->At(i));
703 for (Int_t j = 0; j < fSampleRate[ddl]; j++) {
704 AliDebug(10,Form("Resetting strip %s=%d",
705 digit->GetName(),digit->Counts()));
706 digit->SetCount(j, kBadSignal);
711 // if (errors && (AliDebugLevel() > 0)) input.HexDumpChannel();
716 //____________________________________________________________________
718 AliFMDRawReader::ReadAdcs(AliFMDUShortMap& map)
720 // Read ADC values from raw input into passed TClonesArray of AliFMDDigit
722 AliDebug(3,Form("Reading ADC values into a map"));
724 // const UShort_t kUShortMax = (1 << 16) - 1;
725 for (Int_t ddl = 0; ddl < kNDDL; ddl++) fNErrors[ddl] = 0;
727 AliAltroRawStreamV3 input(fReader);
729 input.SetCheckAltroPayload(false);
730 input.SelectRawData("FMD");
732 // Loop over input RORCs
733 while (input.NextDDL()) {
735 Int_t ddl = NewDDL(input, det);
739 while (input.NextChannel()) {
740 // Get the hardware address, and map that to detector coordinates
744 Int_t hwaddr = NewChannel(input, det, ring, sec, strbase);
745 if (hwaddr < 0) break;
746 if (hwaddr > 0xFFF) continue;
748 UShort_t start = 0x3FF;
749 Bool_t errors = false;
753 while (input.NextBunch()) {
754 // Get Lenght of bunch, and pointer to the data
755 // const UShort_t* data = input.GetSignals();
757 if (!NewBunch(input, start, length)) {
763 // Loop over the data and store it.
764 for (Int_t i = 0; i < length; i++) {
769 Int_t adc = NewSample(input, i, t, sec, strbase, str, samp);
770 if (adc < 0) continue;
771 UShort_t counts = adc;
773 AliDebug(10, Form("FMD%d%c[%02d,%03d]-%d: %4d",
774 det, ring, sec, str, samp, counts));
775 if (SelectSample(samp, fSampleRate[ddl]))
776 map(det,ring,sec,str) = counts;
777 if (first < 0) first = str;
782 AliWarning(Form("Channel %3d/0x%03x contain errors, "
783 "resetting strips %d to %d", ddl, hwaddr, first, last));
785 Int_t ds = first <= last ? 1 : -1;
786 for (Int_t i = first; i != last+ds; i += ds) {
787 AliDebug(10, Form("Resetting strip FMD%d%c[%02d,%03d]=%d",
788 det,ring,sec,i,map(det,ring,sec,i)));
789 map(det,ring,sec,i) = kBadSignal;
798 //____________________________________________________________________
799 Bool_t AliFMDRawReader::ReadSODevent(AliFMDCalibSampleRate* sampleRate,
800 AliFMDCalibStripRange* stripRange,
802 TArrayS &pulseLength,
806 // Read SOD event into passed objects.
809 // samplerate The sample rate object to fill
810 // striprange The strip range object to fill
811 // pulseSize The pulse size object to fill
812 // pulseLength The pulse length (in events) object to fill
815 // @c true on success
817 AliDebug(0,Form("Start of SOD/EOD"));
819 UInt_t shift_clk[18];
820 UInt_t sample_clk[18];
821 UInt_t strip_low[18];
822 UInt_t strip_high[18];
823 UInt_t pulse_size[18];
824 UInt_t pulse_length[18];
825 for (size_t i = 0; i < 18; i++) {
833 AliFMDParameters* param = AliFMDParameters::Instance();
834 AliFMDAltroMapping* map = param->GetAltroMap();
836 AliAltroRawStreamV3 streamer(fReader);
838 streamer.SelectRawData("FMD");
839 //fReader->GetDDLID();
840 while (streamer.NextDDL()) {
841 Int_t ddl = streamer.GetDDLNumber();
842 Int_t detID = fReader->GetDetectorID();
843 if (detectors) detectors[map->DDL2Detector(ddl)-1] = kTRUE;
844 AliDebug(0,Form(" From reader: DDL number is %d , det ID is %d",ddl,detID));
846 ULong_t nPayloadWords = streamer.GetRCUPayloadSizeInSOD();
847 UChar_t* payloadData = streamer.GetRCUPayloadInSOD();
848 UInt_t* payloadWords = reinterpret_cast<UInt_t*>(payloadData);
849 //UInt_t* payloadWords = streamer.GetRCUPayloadInSOD();
851 //std::cout<<nPayloadWords<<" "<<ddl<<std::endl;
852 for (ULong_t i = 1; i <= nPayloadWords ; i++, payloadWords++) {
853 UInt_t payloadWord = *payloadWords; // Get32bitWord(i);
855 //std::cout<<i<<Form(" word: 0x%x",payloadWord)<<std::endl;
856 // address is only 24 bit
857 UInt_t address = (0xffffff & payloadWord);
858 UInt_t type = ((address >> 21) & 0xf);
859 UInt_t error = ((address >> 20) & 0x1);
860 UInt_t bcast = ((address >> 18) & 0x1);
861 UInt_t bc_not_altro = ((address >> 17) & 0x1);
862 UInt_t board = ((address >> 12) & 0x1f);
863 UInt_t instruction = 0;
867 instruction = address & 0xfff;
869 chip = ((address >> 9) & 0x7);
870 channel = ((address >> 5) & 0x5);
871 instruction = (address & 0x1f);
874 Bool_t readDataWord = kFALSE;
876 case 0x0: // Fec read
877 readDataWord = kTRUE;
879 case 0x2: // Fec write
886 case 0x6: // End sequence
888 i = nPayloadWords + 1;
894 //Don't read unless we have a FEC_RD
895 if(!readDataWord) continue;
897 UInt_t dataWord = *payloadWords;//Get32bitWord(i);
898 UInt_t data = (0xFFFFF & dataWord) ;
899 //UInt_t data = (0xFFFF & dataWord) ;
902 AliWarning(Form("error bit detected at Word 0x%06x; "
903 "error % d, type %d, bc_not_altro %d, "
904 "bcast %d, board 0x%02x, chip 0x%x, "
905 "channel 0x%02x, instruction 0x%03x",
906 address, error, type, bc_not_altro,
907 bcast,board,chip,channel,instruction));
913 switch(instruction) {
915 case 0x01: break; // First ADC T
916 case 0x02: break; // I 3.3 V
917 case 0x03: break; // I 2.5 V altro digital
918 case 0x04: break; // I 2.5 V altro analog
919 case 0x05: break; // I 2.5 V VA
920 case 0x06: break; // First ADC T
921 case 0x07: break; // I 3.3 V
922 case 0x08: break; // I 2.5 V altro digital
923 case 0x09: break; // I 2.5 V altro analog
924 case 0x0A: break; // I 2.5 V VA
925 case 0x2D: break; // Second ADC T
926 case 0x2E: break; // I 1.5 V VA
927 case 0x2F: break; // I -2.0 V
928 case 0x30: break; // I -2.0 V VA
929 case 0x31: break; // 2.5 V Digital driver
930 case 0x32: break; // Second ADC T
931 case 0x33: break; // I 1.5 V VA
932 case 0x34: break; // I -2.0 V
933 case 0x35: break; // I -2.0 V VA
934 case 0x36: break; // 2.5 V Digital driver
935 case 0x37: break; // Third ADC T
936 case 0x38: break; // Temperature sens. 1
937 case 0x39: break; // Temperature sens. 2
938 case 0x3A: break; // U 2.5 altro digital (m)
939 case 0x3B: break; // U 2.5 altro analog (m)
940 case 0x3C: break; // Third ADC T
941 case 0x3D: break; // Temperature sens. 1
942 case 0x3E: break; // Temperature sens. 2
943 case 0x3F: break; // U 2.5 altro digital (m)
944 case 0x40: break; // U 2.5 altro analog (m)
945 case 0x41: break; // Forth ADC T
946 case 0x42: break; // U 2.5 VA (m)
947 case 0x43: break; // U 1.5 VA (m)
948 case 0x44: break; // U -2.0 VA (m)
949 case 0x45: break; // U -2.0 (m)
950 case 0x46: break; // Forth ADC T
951 case 0x47: break; // U 2.5 VA (m)
952 case 0x48: break; // U 1.5 VA (m)
953 case 0x49: break; // U -2.0 VA (m)
954 case 0x4A: break; // U -2.0 (m)
956 case 0x0B: break; // L1 trigger CouNTer
957 case 0x0C: break; // L2 trigger CouNTer
958 case 0x0D: break; // Sampling CLK CouNTer
959 case 0x0E: break; // DSTB CouNTer
961 case 0x0F: break; // Test mode word
962 case 0x10: break; // Undersampling ratio.
963 // Configuration and status
964 case 0x11: break; // Config/Status Register 0
965 case 0x12: break; // Config/Status Register 1
966 case 0x13: break; // Config/Status Register 2
967 case 0x14: break; // Config/Status Register 3
968 case 0x15: break; // Free
970 case 0x16: break; // Latch L1, L2, SCLK Counters
971 case 0x17: break; // Clear counters
972 case 0x18: break; // Clear CSR1
973 case 0x19: break; // rstb ALTROs
974 case 0x1A: break; // rstb BC
975 case 0x1B: break; // Start conversion
976 case 0x1C: break; // Scan event length
977 case 0x1D: break; // Read event length
978 case 0x1E: break; // Start test mode
979 case 0x1F: break; // Read acquisition memory
981 case 0x20: break; // FMDD status
982 case 0x21: break; // L0 counters
983 case 0x22: break; // FMD: Wait to hold
984 case 0x23: break; // FMD: L1 timeout
985 case 0x24: break; // FMD: L2 timeout
986 case 0x25: // FMD: Shift clk
987 shift_clk[board] = ((data >> 8 ) & 0xFF);
988 AliDebug(30,Form("Read shift_clk=%d for board 0x%02x",
989 shift_clk[board], board));
991 case 0x26: // FMD: Strips
992 strip_low[board] = ((data >> 0 ) & 0xFF);
993 strip_high[board] = ((data >> 8 ) & 0xFF);
995 case 0x27: // FMD: Cal pulse
996 pulse_size[board] = ((data >> 8 ) & 0xFF);
998 case 0x28: break; // FMD: Shape bias
999 case 0x29: break; // FMD: Shape ref
1000 case 0x2A: break; // FMD: Preamp ref
1001 case 0x2B: // FMD: Sample clk
1002 sample_clk[board] = ((data >> 8 ) & 0xFF);
1003 AliDebug(30,Form("Read sample_clk=%d for board 0x%02x",
1004 sample_clk[board], board));
1006 case 0x2C: break; // FMD: Commands
1007 case 0x4B: // FMD: Cal events
1008 pulse_length[board] = ((data >> 0 ) & 0xFF);
1013 AliDebug(50,Form("instruction 0x%x, dataword 0x%x",
1014 instruction,dataWord));
1015 } // End of loop over Result memory event
1018 UShort_t sector = 0;
1022 const UInt_t boards[4] = {0,1,16,17};
1023 for(Int_t i=0;i<4;i++) {
1024 if(ddl==0 && (i==1 || i==3)) continue;
1026 UInt_t chip =0, channel=0;
1027 det = map->DDL2Detector(ddl);
1028 map->Channel2StripBase(boards[i], chip, channel, ring, sector, strip);
1030 UInt_t samplerate = 0;
1032 if(sample_clk[boards[i]] == 0) {
1034 Int_t sample1 = sample_clk[boards[0]];
1035 Int_t sample2 = sample_clk[boards[2]];
1036 if(sample1) sample_clk[boards[i]] = sample1;
1037 else sample_clk[boards[i]] = sample2;
1040 Int_t sample1 = sample_clk[boards[0]];
1041 Int_t sample2 = sample_clk[boards[1]];
1042 Int_t sample3 = sample_clk[boards[2]];
1043 Int_t sample4 = sample_clk[boards[3]];
1044 Int_t agreement = 0;
1045 if(sample1 == sample2) agreement++;
1046 if(sample1 == sample3) agreement++;
1047 if(sample1 == sample4) agreement++;
1048 if(sample2 == sample3) agreement++;
1049 if(sample2 == sample4) agreement++;
1050 if(sample3 == sample4) agreement++;
1055 if(agreement == 3) {
1056 sample_clk[boards[i]] = sample_clk[boards[idx]];
1057 shift_clk[boards[i]] = shift_clk[boards[idx]];
1058 strip_low[boards[i]] = strip_low[boards[idx]];
1059 strip_high[boards[i]] = strip_high[boards[idx]];
1060 pulse_length[boards[i]] = pulse_length[boards[idx]];
1061 pulse_size[boards[i]] = pulse_size[boards[idx]];
1062 AliDebug(3,Form("Vote taken for ddl %d, board 0x%x",
1069 if(sample_clk[boards[i]])
1070 samplerate = shift_clk[boards[i]]/sample_clk[boards[i]];
1071 AliDebug(10,Form("Sample rate for board 0x%02x is %d",
1072 boards[i], samplerate));
1073 sampleRate->Set(det,ring,sector,0,samplerate);
1074 stripRange->Set(det,ring,sector,0,
1075 strip_low[boards[i]],strip_high[boards[i]]);
1077 AliDebug(20,Form("det %d, ring %c, ",det,ring));
1078 pulseLength.AddAt(pulse_length[boards[i]],
1079 GetHalfringIndex(det,ring,boards[i]/16));
1080 pulseSize.AddAt(pulse_size[boards[i]],
1081 GetHalfringIndex(det,ring,boards[i]/16));
1085 AliDebug(20,Form(": Board: 0x%02x\n"
1086 "\tstrip_low %3d, strip_high %3d\n"
1087 "\tshift_clk %3d, sample_clk %3d\n"
1088 "\tpulse_size %3d, pulse_length %3d",
1090 strip_low[boards[i]], strip_high[boards[i]],
1091 shift_clk[boards[i]], sample_clk[boards[i]],
1092 pulse_size[boards[i]],pulse_length[boards[i]]));
1097 AliFMDParameters::Instance()->SetSampleRate(sampleRate);
1098 AliFMDParameters::Instance()->SetStripRange(stripRange);
1100 AliDebug(0,Form("End of SOD/EOD"));
1104 //____________________________________________________________________
1106 UInt_t AliFMDRawReader::Get32bitWord(Int_t idx)
1108 // This method returns the 32 bit word at a given
1109 // position inside the raw data payload.
1110 // The 'index' points to the beginning of the next word.
1111 // The method is supposed to be endian (platform)
1114 AliFatal("Raw data paylod buffer is not yet initialized !");
1117 Int_t index = 4*idx;
1120 // fRawReader->AddFatalErrorLog(k32bitWordReadErr,Form("pos = %d",index));
1122 AliWarning(Form("Invalid raw data payload position (%d) !",index));
1127 word = fData[--index] << 24;
1128 word |= fData[--index] << 16;
1129 word |= fData[--index] << 8;
1130 word |= fData[--index] << 0 ;
1134 //_____________________________________________________________________
1135 Int_t AliFMDRawReader::GetHalfringIndex(UShort_t det, Char_t ring,
1136 UShort_t board) const
1139 // Get short index for a given half-ring
1142 // det Detector number
1143 // ring Ring identifer
1144 // board Board number
1149 UShort_t iring = (ring == 'I' ? 1 : 0);
1151 Int_t index = (((det-1) << 2) | (iring << 1) | (board << 0));
1157 //____________________________________________________________________