From 20daa34d89d94f7b5c29b9b8ff9df29376036f57 Mon Sep 17 00:00:00 2001 From: cvetan Date: Fri, 17 Mar 2006 15:45:41 +0000 Subject: [PATCH] New methods to read and write trailers. New method to read whole altro channel --- RAW/AliAltroBuffer.cxx | 138 +++++++++++++++++++++++++++++++++++++---- RAW/AliAltroBuffer.h | 14 ++++- 2 files changed, 139 insertions(+), 13 deletions(-) diff --git a/RAW/AliAltroBuffer.cxx b/RAW/AliAltroBuffer.cxx index d9e6350a236..0b6cfba54ce 100644 --- a/RAW/AliAltroBuffer.cxx +++ b/RAW/AliAltroBuffer.cxx @@ -292,6 +292,17 @@ void AliAltroBuffer::WriteTrailer(Int_t wordsNumber, Int_t padNumber, rowNumber,secNumber); } + Short_t hwAdress = fMapping->GetHWAdress(rowNumber,padNumber,secNumber); + if (hwAdress == -1) + AliFatal(Form("No hardware (ALTRO) adress found for these pad-row (%d) and pad (%d) indeces !",rowNumber,padNumber)); + WriteTrailer(wordsNumber,hwAdress); +} + +//_____________________________________________________________________________ +void AliAltroBuffer::WriteTrailer(Int_t wordsNumber, Short_t hwAdress) +{ +//Writes a trailer of 40 bits using +//a given hardware adress Int_t num = fFreeCellBuffer % 4; for(Int_t i = 0; i < num; i++) { FillBuffer(0x2AA); @@ -305,10 +316,6 @@ void AliAltroBuffer::WriteTrailer(Int_t wordsNumber, Int_t padNumber, temp = (wordsNumber << 6) & 0x3FF; temp |= (0xA << 2); - Short_t hwAdress = fMapping->GetHWAdress(rowNumber,padNumber,secNumber); - if (hwAdress == -1) - AliFatal(Form("No hardware (ALTRO) adress found for these pad-row (%d) and pad (%d) indeces !",rowNumber,padNumber)); - temp |= (hwAdress >> 10) & 0x3; FillBuffer(temp); temp = hwAdress & 0x3FF; @@ -343,6 +350,20 @@ Bool_t AliAltroBuffer::ReadTrailer(Int_t& wordsNumber, Int_t& padNumber, rowNumber,secNumber); } + Short_t hwAdress; + if (!ReadTrailer(wordsNumber,hwAdress)) return kFALSE; + rowNumber = fMapping->GetPadRow(hwAdress); + padNumber = fMapping->GetPad(hwAdress); + secNumber = fMapping->GetSector(hwAdress); + + return kTRUE; +} + +//_____________________________________________________________________________ +Bool_t AliAltroBuffer::ReadTrailer(Int_t& wordsNumber, Short_t& hwAdress) +{ +//Read a trailer of 40 bits in the forward reading mode + Int_t temp = GetNext(); if (temp != 0x2AA) AliFatal(Form("Incorrect trailer found ! Expecting 0x2AA but found %x !",temp)); @@ -354,23 +375,19 @@ Bool_t AliAltroBuffer::ReadTrailer(Int_t& wordsNumber, Int_t& padNumber, temp = GetNext(); wordsNumber |= (temp >> 6); - if ((temp & 0xF) != 0xA) + if (((temp >> 2) & 0xF) != 0xA) AliFatal(Form("Incorrect trailer found ! Expecting second 0xA but found %x !",temp >> 6)); - Int_t hwAdress = (temp & 0x3) << 10; + hwAdress = (temp & 0x3) << 10; temp = GetNext(); hwAdress |= temp; - rowNumber = fMapping->GetPadRow(hwAdress); - padNumber = fMapping->GetPad(hwAdress); - secNumber = fMapping->GetSector(hwAdress); return kTRUE; } //_____________________________________________________________________________ -Bool_t AliAltroBuffer::ReadTrailerBackward(Int_t& wordsNumber, - Int_t& padNumber, - Int_t& rowNumber, Int_t& secNumber) +Bool_t AliAltroBuffer::ReadDummyTrailerBackward(Int_t& wordsNumber, Int_t& padNumber, + Int_t& rowNumber, Int_t& secNumber) { //Read a trailer of 40 bits in the backward reading mode @@ -392,6 +409,59 @@ Bool_t AliAltroBuffer::ReadTrailerBackward(Int_t& wordsNumber, return kTRUE; } +//_____________________________________________________________________________ +Bool_t AliAltroBuffer::ReadTrailerBackward(Int_t& wordsNumber, Int_t& padNumber, + Int_t& rowNumber, Int_t& secNumber) +{ +//Read a trailer of 40 bits in the backward reading mode + if (!fMapping) { + AliError("No ALTRO mapping information is loaded! Reading a dummy trailer!"); + return ReadDummyTrailerBackward(wordsNumber,padNumber, + rowNumber,secNumber); + } + + Short_t hwAdress; + if (!ReadTrailerBackward(wordsNumber,hwAdress)) return kFALSE; + rowNumber = fMapping->GetPadRow(hwAdress); + padNumber = fMapping->GetPad(hwAdress); + secNumber = fMapping->GetSector(hwAdress); + + return kTRUE; +} + +//_____________________________________________________________________________ +Bool_t AliAltroBuffer::ReadTrailerBackward(Int_t& wordsNumber, Short_t& hwAdress) +{ +//Read a trailer of 40 bits in the backward reading mode + + Int_t temp; + fEndingFillWords = 0; + do { + temp = GetNextBackWord(); + fEndingFillWords++; + if (temp == -1) return kFALSE; + } while (temp == 0x2AA); + fEndingFillWords--; + + hwAdress = temp; + + temp = GetNextBackWord(); + hwAdress |= (temp & 0x3) << 10; + if (((temp >> 2) & 0xF) != 0xA) + AliFatal(Form("Incorrect trailer found ! Expecting second 0xA but found %x !",temp >> 6)); + wordsNumber = (temp >> 6); + + temp = GetNextBackWord(); + wordsNumber |= (temp << 4) & 0x3FF; + if ((temp >> 6) != 0xA) + AliFatal(Form("Incorrect trailer found ! Expecting 0xA but found %x !",temp >> 6)); + + temp = GetNextBackWord(); + if (temp != 0x2AA) + AliFatal(Form("Incorrect trailer found ! Expecting 0x2AA but found %x !",temp)); + + return kTRUE; +} //_____________________________________________________________________________ void AliAltroBuffer::WriteChannel(Int_t padNumber, Int_t rowNumber, @@ -432,6 +502,50 @@ void AliAltroBuffer::WriteChannel(Int_t padNumber, Int_t rowNumber, WriteTrailer(nWords, padNumber, rowNumber, secNumber); } +//_____________________________________________________________________________ +void AliAltroBuffer::ReadChannel(Int_t padNumber, Int_t rowNumber, + Int_t secNumber, + Int_t& nTimeBins, Int_t* adcValues) +{ +//Read all ADC values and the trailer of a channel + + Int_t wordsNumber; + if (!ReadTrailer(wordsNumber,padNumber, + rowNumber,secNumber)) return; + + if (wordsNumber < 0) return; + // Number of fill words + Int_t nFillWords; + if ((wordsNumber % 4) == 0) + nFillWords = 0; + else + nFillWords = 4 - wordsNumber % 4; + // Read the fill words + for (Int_t i = 0; i < nFillWords; i++) { + Int_t temp = GetNext(); + if (temp != 0x2AA) + AliFatal(Form("Invalid fill word, expected 0x2AA, but got %X", temp)); + } + + // Decoding + Int_t lastWord = wordsNumber; + nTimeBins = -1; + while (lastWord > 0) { + Int_t l = GetNext(); + if (l < 0) AliFatal(Form("Bad bunch length (%d) !", l)); + Int_t t = GetNext(); + if (t < 0) AliFatal(Form("Bad bunch time (%d) !", t)); + lastWord -= 2; + if (nTimeBins == -1) nTimeBins = t + 1; + for (Int_t i = 2; i < l; i++) { + Int_t amp = GetNext(); + if (amp < 0) AliFatal(Form("Bad adc value (%X) !", amp)); + adcValues[t - (i-2)] = amp; + lastWord--; + } + } + +} //_____________________________________________________________________________ void AliAltroBuffer::WriteDataHeader(Bool_t dummy, Bool_t compressed) diff --git a/RAW/AliAltroBuffer.h b/RAW/AliAltroBuffer.h index 9cb46b62922..ce27658c029 100644 --- a/RAW/AliAltroBuffer.h +++ b/RAW/AliAltroBuffer.h @@ -43,23 +43,35 @@ class AliAltroBuffer: public TObject { void WriteTrailer(Int_t wordsNumber, Int_t padNumber, Int_t rowNumber, Int_t secNumber); //this method is used to write the trailer + void WriteTrailer(Int_t wordsNumber, Short_t hwAdress); + //this method is used to write the trailer void WriteDummyTrailer(Int_t wordsNumber, Int_t padNumber, Int_t rowNumber, Int_t secNumber); //this method is used to write a dummy trailer Bool_t ReadTrailer(Int_t& wordsNumber, Int_t& padNumber, Int_t& rowNumber, Int_t &secNumber); //this method is used to read the trailer when the file is read forward + Bool_t ReadTrailer(Int_t& wordsNumber, Short_t& hwAdress); + //this method is used to read the trailer when the file is read forward Bool_t ReadDummyTrailer(Int_t& wordsNumber, Int_t& padNumber, Int_t& rowNumber, Int_t &secNumber); //this method is used to read the trailer when the file is read forward Bool_t ReadTrailerBackward(Int_t& wordsNumber, Int_t& padNumber, Int_t& rowNumber, Int_t& secNumber); //this method is used to read the trailer when the file is read backward + Bool_t ReadTrailerBackward(Int_t& wordsNumber, Short_t& hwAdress); + //this method is used to read the trailer when the file is read backward + Bool_t ReadDummyTrailerBackward(Int_t& wordsNumber, Int_t& padNumber, + Int_t& rowNumber, Int_t& secNumber); + //this method is used to read the trailer when the file is read backward void WriteChannel(Int_t padNumber, Int_t rowNumber, Int_t secNumber, Int_t nTimeBins, const Int_t* adcValues, Int_t threshold = 0); //this method is used to write all ADC values and the trailer of a channel + void ReadChannel(Int_t padNumber, Int_t rowNumber, Int_t secNumber, + Int_t& nTimeBins, Int_t* adcValues); + //this method is used to read all ADC values and the trailer of a channel void WriteDataHeader(Bool_t dummy, Bool_t compressed); //this method is used to write the data header @@ -75,7 +87,7 @@ class AliAltroBuffer: public TObject { void SetMapping(AliAltroMapping *mapping) { fMapping = mapping; } - private: + protected: AliAltroBuffer(const AliAltroBuffer& source); AliAltroBuffer& operator = (const AliAltroBuffer& source); -- 2.39.3