/* ======================================================================== *\ ! ! * ! * This file is part of MARS, the MAGIC Analysis and Reconstruction ! * Software. It is distributed to you in the hope that it can be a useful ! * and timesaving tool in analysing Data of imaging Cerenkov telescopes. ! * It is distributed WITHOUT ANY WARRANTY. ! * ! * Permission to use, copy, modify and distribute this software and its ! * documentation for any purpose is hereby granted without fee, ! * provided that the above copyright notice appear in all copies and ! * that both that copyright notice and this permission notice appear ! * in supporting documentation. It is provided "as is" without express ! * or implied warranty. ! * ! ! ! Author(s): Thomas Bretz, 12/2000 ! ! Copyright: MAGIC Software Development, 2000-2004 ! ! \* ======================================================================== */ ///////////////////////////////////////////////////////////////////////////// // // MRawEvtData // // Storage container to store the raw FADC values. // // MArrayS fHiGainPixId // --------------------- // Array of Pixel Numbers for their high voltage channel in the order the // FADC values are stored in fHiGainFadcSamples // // MArrayB fHiGainFadcSaples // ------------------------- // FADC samples (hi gain) of all pixels // // MArrayS fLoGainPixId // -------------------- // see fHiGainPixId // // MArrayB fLoGainFadcSamples // -------------------------- // see fHiGainFadcSaples // // // Version 5 (0.8.5): // ------------------ // - Changed type of ABFlags from TArrayC to MArrayB // // Version 4 (0.8.4): // ------------------ // - Changed derivement from MCamEvent to MParContainer and MCamEvent // // Version 3 (0.8.4): // ------------------ // - Added fABFlags // // Version 2: // ---------- // - Derives from MCamEvent now // // Version 1: // ---------- // - First implementation // ///////////////////////////////////////////////////////////////////////////// #include "MRawEvtData.h" #include #include #include #include #include #include "MLog.h" #include "MLogManip.h" #include "MArrayS.h" #include "MArrayB.h" #include "MGeomCam.h" #include "MRawCrateArray.h" #include "MRawCrateData.h" #include "MRawRunHeader.h" #include "MRawEvtPixelIter.h" ClassImp(MRawEvtData); using namespace std; // -------------------------------------------------------------------------- // // Default constructor. It initializes all arrays with zero size. // MRawEvtData::MRawEvtData(const char *name, const char *title) : fRunHeader(0) { fName = name ? name : "MRawEvtData"; fTitle = title ? title : "Raw Event Data Information"; InitArrays(); } // -------------------------------------------------------------------------- // // Destructor. Deletes all the arrays. // MRawEvtData::~MRawEvtData() { DeleteArrays(); } // -------------------------------------------------------------------------- // // reset all arrays // void MRawEvtData::Clear(Option_t *) { /* FIXME: Is Reset (set all entries to zero) what you want to do or Set(0) (delete the array) */ fHiGainPixId->Reset(); fLoGainPixId->Reset(); fHiGainFadcSamples->Reset(); fLoGainFadcSamples->Reset(); fABFlags->Reset(); } // -------------------------------------------------------------------------- // // return the number of hi gain samples per pixel // Byte_t MRawEvtData::GetNumHiGainSamples() const { return fHiGainPixId->GetSize() ? fHiGainFadcSamples->GetSize()/fHiGainPixId->GetSize() : 0; } // -------------------------------------------------------------------------- // // return the number of lo gain samples per pixel // Byte_t MRawEvtData::GetNumLoGainSamples() const { return fLoGainPixId->GetSize() ? fLoGainFadcSamples->GetSize()/fLoGainPixId->GetSize() : 0; } // -------------------------------------------------------------------------- // // return the number of stored pixel // UShort_t MRawEvtData::GetNumPixels() const { return fHiGainPixId->GetSize(); } // -------------------------------------------------------------------------- // // Print out the onformation to *fLog. // Options: // "hex" Prints the time slices hexadecimal (default) // "dec" Prints the time slices decimal // void MRawEvtData::Print(Option_t *opt) const { // // print fadc inforation to screen // Possible Options: // - DEC: Print values decimal instead of hexadecimal (default) // const Byte_t nHiSamp = GetNumHiGainSamples(); const Byte_t nLoSamp = GetNumLoGainSamples(); const UShort_t nHiPix = fHiGainPixId->GetSize();; const UShort_t nLoPix = fLoGainPixId->GetSize();; fLog->unsetf(ios::showbase); *fLog << dec << all; *fLog << GetDescriptor() << ": " << endl; *fLog << "HiGain: " << nHiPix << " Pixels with " << (Int_t)nHiSamp << " Samples "; *fLog << "LoGain: " << nLoPix << " Pixels with " << (Int_t)nLoSamp << " Samples"; TString str(opt); Int_t manip = str.Contains("dec", TString::kIgnoreCase); Int_t l=0; for (int i=0; iGetSize()) { *fLog << "<" << hex << setfill('0') << setw(2); *fLog << ((Int_t)(*fABFlags)[idx/8]&0xff) << "> "; } *fLog << (manip?dec:hex) << (manip ? setfill(' ') : setfill('0')); for (int j=0; j The pixel with the given index is drawn // void MRawEvtData::Draw(Option_t *opt) { if (GetNumPixels()==0) { *fLog << warn << "Sorry, no pixel to draw!" << endl; return; } TString str(opt); str.ToLower(); UInt_t id = 0; if (str.BeginsWith("graph")) if (str.Length()>5) sscanf(&str[5], "%d", &id); if (str.BeginsWith("hist")) if (str.Length()>4) sscanf(&str[4], "%d", &id); MRawEvtPixelIter pix(this); if (!pix.Jump(id)) { *fLog << warn << dec << "Pixel Idx #" << id << " doesn't exist!" << endl; return; } const Byte_t *higains = pix.GetHiGainSamples(); const Byte_t *logains = pix.GetLoGainSamples(); const Int_t nh = GetNumHiGainSamples(); const Int_t nl = GetNumLoGainSamples(); TString name = "Pixel Idx."; name += pix.GetPixelId(); Bool_t same = str.Contains("same"); if (str.BeginsWith("graph")) { *fLog << inf << "Drawing Graph: Pixel Idx #" << dec << pix.GetPixelId(); *fLog << " of " << (int)GetNumPixels() << "Pixels" << endl; TGraph *graphhi = new TGraph; for (int i=0; iSetPoint(graphhi->GetN(), i, higains[i]); graphhi->SetMaximum(256); graphhi->SetMinimum(0); graphhi->SetBit(kCanDelete); graphhi->Draw(same ? "C*" : "AC*"); TH1F *histhi = graphhi->GetHistogram(); histhi->SetMinimum(0); histhi->SetMaximum(255); histhi->SetXTitle("Time/FADC Slices"); histhi->SetYTitle("Signal/FADC Units"); if (nl>0) { TGraph *graphlo = new TGraph; for (int i=0; iSetPoint(graphlo->GetN(), i, logains[i]); graphlo->SetMaximum(256); graphlo->SetMinimum(0); graphlo->SetLineColor(kBlue); graphlo->SetBit(kCanDelete); graphlo->Draw("C*"); TH1F *histlo = graphlo->GetHistogram(); histlo->SetXTitle("Time/FADC Slices"); histlo->SetYTitle("Signal/FADC Units"); } return; } if (str.BeginsWith("hist")) { // FIXME: Add Legend *fLog << inf << "Drawing Histogram of Pixel with Idx #" << dec << pix.GetPixelId() << " to " << gPad << endl; TH1F *histh = new TH1F(name, "FADC Samples", nh, -0.5, nh-.5); histh->SetMinimum(0); histh->SetMaximum(255); histh->SetXTitle("Time [FADC Slices]"); histh->SetYTitle("Signal [FADC Units]"); histh->SetDirectory(NULL); for (int i=0; iFill(i, higains[i]); histh->SetBit(kCanDelete); histh->Draw(same ? "same" : ""); if (nl>0) { TH1F *histl = new TH1F(name+";2", "FADC Samples", nl, -0.5, nl-.5); histl->SetLineColor(kBlue); histl->SetDirectory(NULL); for (int i=0; iFill(i, logains[i]); histl->SetBit(kCanDelete); histl->Draw("same"); } return; } *fLog << warn << dbginf << "Warning - You must specify either 'GRAPH' or 'HIST'" << endl; } // -------------------------------------------------------------------------- // // Deletes all the arrays // The flag is for future usage. // void MRawEvtData::InitArrays(UShort_t numconnected, UShort_t maxid) { const Int_t numlo = fRunHeader ? fRunHeader->GetNumSamplesHiGain() : 0; const Int_t numhi = fRunHeader ? fRunHeader->GetNumSamplesHiGain() : 0; fHiGainPixId = new MArrayS(numconnected); fLoGainPixId = new MArrayS(numconnected); fHiGainFadcSamples = new MArrayB(numconnected*numhi); fLoGainFadcSamples = new MArrayB(numconnected*numlo); fABFlags = new MArrayB(maxid/8+1); fConnectedPixels = 0; } // -------------------------------------------------------------------------- // // Deletes all the arrays // void MRawEvtData::DeleteArrays() { delete fHiGainPixId; delete fLoGainPixId; delete fHiGainFadcSamples; delete fLoGainFadcSamples; delete fABFlags; } // -------------------------------------------------------------------------- // // Deletes all arrays describing the pixel Id and Samples in pixels. // The flag is for future usage. // void MRawEvtData::ResetPixels(UShort_t numconnected, UShort_t maxid) { //const int npix = fRunHeader->GetNumCrates()*fRunHeader->GetNumPixInCrate(); if (fHiGainPixId && fHiGainPixId->GetSize()==numconnected && (UShort_t)fABFlags->GetSize()==(maxid/8+1)) { fConnectedPixels = 0; return; } DeleteArrays(); InitArrays(numconnected, maxid); } // -------------------------------------------------------------------------- // // This is to fill the data of one pixel to the MRawEvtHeader Class. // The parameters are the pixelnumber and the FADC_SLICES values of ADCs // Add to lo gains if lflag = 1 // void MRawEvtData::AddPixel(UShort_t nOfPixel, TArrayC *data, Bool_t lflag) { MArrayS *arrpix = lflag ? fLoGainPixId : fHiGainPixId; MArrayB *arrsam = lflag ? fLoGainFadcSamples : fHiGainFadcSamples; // // check whether we got the right number of new samples // if there are no samples already stored: this is the new number of samples // const Byte_t ns = data->GetSize(); const Byte_t nSamp = lflag ? GetNumLoGainSamples() : GetNumHiGainSamples(); if (nSamp && ns!=nSamp) { *fLog << err << "RawEvtData::AddPixel: Error, number of samples in "; *fLog << "TArrayC " << ns << " doesn't match current number " << nSamp << endl; return; } // // enhance pixel array by one // arrpix->Set(arrpix->GetSize()+1); // // add the number of the new pixel to the array as last entry // arrpix->AddAt(nOfPixel, arrpix->GetSize()-1); // // enhance the array by the number of new samples // arrsam->Set(arrsam->GetSize()+ns); // // add the new slices as last entries to array // arrsam->AddAt((Byte_t*)data->GetArray(), arrsam->GetSize()-ns, ns); } void MRawEvtData::ReadPixel(istream &fin, Int_t npix, Bool_t ab) { Byte_t *poshi = fHiGainFadcSamples->GetArray() + fConnectedPixels*fRunHeader->GetNumSamplesHiGain(); Byte_t *poslo = fLoGainFadcSamples->GetArray() + fConnectedPixels*fRunHeader->GetNumSamplesLoGain(); // Read data for one pixel fin.read((char*)poshi, fRunHeader->GetNumSamplesHiGain()); fin.read((char*)poslo, fRunHeader->GetNumSamplesLoGain()); // // This is to fill the data of one pixel to the MRawEvtHeader Class. // The parameters are the pixelnumber and the FADC_SLICES values of ADCs // Add to lo gains if lflag = 1 // fHiGainPixId->AddAt(npix, fConnectedPixels); // FIXME: Not implemented in the raw files yet //if (IsLoGainOn(i, j)) //{ fLoGainPixId->AddAt(npix, fConnectedPixels); //} if (ab) SETBIT((*fABFlags)[npix/8], npix%8); else CLRBIT((*fABFlags)[npix/8], npix%8); fConnectedPixels++; } // -------------------------------------------------------------------------- // // Fills members with information from a magic binary file. // WARNING: you have to use Init() before you can do this // /* void MRawEvtData::ReadEvt(istream &fin, Int_t posinarray) { const UShort_t npic = fRunHeader->GetNumPixInCrate(); const UShort_t npos = npic*posinarray; //const Byte_t ab = fCrateArray->GetEntry(posinarray)->GetABFlags(); for (int i=npos; iGetPixAssignment(ipos); // Check whether the pixel is connected or not if (hwid==0) { const UShort_t n = fRunHeader->GetNumSamplesLoGain()+fRunHeader->GetNumSamplesHiGain(); fin.seekg(n, ios::cur); return; } // -1 converts the hardware pixel Id into the software pixel index const Int_t npix = (Int_t)hwid-1; const Byte_t ab = fCrateArray->GetEntry(posinarray)->GetABFlags(); AddPixel(fin, npix, TESTBIT(ab, i-npos)); } } */ // -------------------------------------------------------------------------- // // Return the size in bytes of one event data block // Int_t MRawEvtData::GetNumBytes() const { const UShort_t nlo = fRunHeader->GetNumSamplesLoGain(); const UShort_t nhi = fRunHeader->GetNumSamplesHiGain(); const UShort_t npic = fRunHeader->GetNumPixInCrate(); return (nhi+nlo)*npic; } // -------------------------------------------------------------------------- // // Make sure, that you skip the whole event. This function only skips a part // of the event - see MRawRead::SkipEvent // void MRawEvtData::SkipEvt(istream &fin) { fin.seekg(GetNumBytes(), ios::cur); } Bool_t MRawEvtData::GetPixelContent(Double_t &val, Int_t idx, const MGeomCam &cam, Int_t type) const { MRawEvtPixelIter Next(const_cast(this)); if (!Next.Jump(idx)) return kFALSE; switch (type) { case 0: val = Next.GetSumHiGainSamples()-(float)GetNumHiGainSamples()*fHiGainFadcSamples->GetArray()[0]; val*= cam.GetPixRatio(idx); break; case 1: val = Next.GetMaxHiGainSample(); break; case 2: val = Next.GetMaxLoGainSample(); break; case 3: val = Next.GetIdxMaxHiGainSample(); break; case 4: val = Next.GetIdxMaxLoGainSample(); break; case 5: val = Next.GetIdxMaxHiLoGainSample(); return val >= 0; } return kTRUE; } void MRawEvtData::Copy(TObject &named) #if ROOT_VERSION_CODE > ROOT_VERSION(3,04,01) const #endif { MRawEvtData &evt = (MRawEvtData &)named; *evt.fHiGainPixId = *fHiGainPixId; *evt.fLoGainPixId = *fLoGainPixId; *evt.fHiGainFadcSamples = *fHiGainFadcSamples; *evt.fLoGainFadcSamples = *fLoGainFadcSamples; *evt.fABFlags = *fABFlags; evt.fConnectedPixels = fConnectedPixels; }