Public Member Functions | Private Member Functions | Private Attributes | List of all members
VDColdboxDataInterface Class Reference

#include <VDColdboxDataInterface.h>

Inheritance diagram for VDColdboxDataInterface:
PDSPTPCDataInterfaceParent

Public Member Functions

 VDColdboxDataInterface (fhicl::ParameterSet const &ps)
 
 ~VDColdboxDataInterface ()
 
int retrieveData (art::Event &evt, std::string inputlabel, std::vector< raw::RawDigit > &raw_digits, std::vector< raw::RDTimeStamp > &rd_timestamps, std::vector< raw::RDStatus > &rdstatuses)
 
void readFragmentsForEvent (art::Event &evt)
 
int retrieveDataAPAListWithLabels (art::Event &evt, std::string inputlabel, std::vector< raw::RawDigit > &raw_digits, std::vector< raw::RDTimeStamp > &rd_timestamps, std::vector< raw::RDStatus > &rdstatuses, std::vector< int > &apalist)
 
int retrieveDataForSpecifiedAPAs (art::Event &evt, std::vector< raw::RawDigit > &raw_digits, std::vector< raw::RDTimeStamp > &rd_timestamps, std::vector< raw::RDStatus > &rdstatuses, std::vector< int > &apalist)
 
- Public Member Functions inherited from PDSPTPCDataInterfaceParent
virtual ~PDSPTPCDataInterfaceParent () noexcept=default
 

Private Member Functions

void _collectRDStatus (std::vector< raw::RDStatus > &rdstatuses)
 
void getFragmentsForEvent (hid_t the_group, RawDigits &raw_digits, RDTimeStamps &timestamps, int apano, int maxchan)
 
void getMedianSigma (const raw::RawDigit::ADCvector_t &v_adc, float &median, float &sigma)
 

Private Attributes

std::map< int, std::vector< std::string > > _input_labels_by_apa
 
std::string logname = "VDColdboxDataInterface"
 
hid_t fPrevStoredHandle = -1
 
hid_t fHDFFile = -1
 
bool fForceOpen
 
std::string fFileInfoLabel
 
int fMaxChan = 1000000
 

Detailed Description

Definition at line 22 of file VDColdboxDataInterface.h.

Constructor & Destructor Documentation

VDColdboxDataInterface::VDColdboxDataInterface ( fhicl::ParameterSet const &  ps)

Definition at line 160 of file VDColdboxDataInterface_tool.cc.

161  : fForceOpen(p.get<bool>("ForceOpen", false)),
162  fFileInfoLabel(p.get<std::string>("FileInfoLabel", "daq")),
163  fMaxChan(p.get<int>("MaxChan",1000000)) {
164 }
std::string string
Definition: nybbler.cc:12
p
Definition: test.py:223
VDColdboxDataInterface::~VDColdboxDataInterface ( )
inline

Definition at line 27 of file VDColdboxDataInterface.h.

27  {
28  if (fForceOpen) {
29  H5Fclose(fHDFFile);
30  }
31  };

Member Function Documentation

void VDColdboxDataInterface::_collectRDStatus ( std::vector< raw::RDStatus > &  rdstatuses)
inlineprivate

Definition at line 57 of file VDColdboxDataInterface.h.

57 {};
void VDColdboxDataInterface::getFragmentsForEvent ( hid_t  the_group,
RawDigits raw_digits,
RDTimeStamps timestamps,
int  apano,
int  maxchan 
)
private

Definition at line 243 of file VDColdboxDataInterface_tool.cc.

245  {
246 
247  using namespace dune::HDF5Utils;
248  using dunedaq::detdataformats::wib::WIBFrame;
249  using dunedaq::detdataformats::wib::WIBHeader;
250 
252 
253  std::deque<std::string> det_types
254  = getMidLevelGroupNames(the_group);
255 
256  for (const auto & det : det_types)
257  {
258  if (det != "TPC") continue;
259  //std::cout << " Detector type: " << det << std::endl;
260  hid_t geoGroup = getGroupFromPath(the_group, det);
261  std::deque<std::string> apaNames
262  = getMidLevelGroupNames(geoGroup);
263 
264  std::cout << "Size of apaNames: " << apaNames.size() << std::endl;
265  std::cout << "apaNames[apano]: " << apaNames[apano-1] << std::endl;
266 
267  // apaNames is a vector whose elements start at [0].
268  hid_t linkGroup = getGroupFromPath(geoGroup, apaNames[apano-1]);
269  std::deque<std::string> linkNames
270  = getMidLevelGroupNames(linkGroup);
271  for (const auto & t : linkNames)
272  {
273  hid_t dataset = H5Dopen(linkGroup, t.data(), H5P_DEFAULT);
274  hsize_t ds_size = H5Dget_storage_size(dataset);
275  if (ds_size <= sizeof(FragmentHeader)) continue; //Too small
276 
277  std::vector<char> ds_data(ds_size);
278  H5Dread(dataset, H5T_STD_I8LE, H5S_ALL, H5S_ALL, H5P_DEFAULT,
279  ds_data.data());
280  H5Dclose(dataset);
281 
282  //Each fragment is a collection of WIB Frames
283  Fragment frag(&ds_data[0], Fragment::BufferAdoptionMode::kReadOnlyMode);
284  size_t n_frames = (ds_size - sizeof(FragmentHeader))/sizeof(WIBFrame);
285  std::vector<raw::RawDigit::ADCvector_t> adc_vectors(256);
286  uint32_t slot = 0, fiber = 0;
287  for (size_t i = 0; i < n_frames; ++i)
288  {
289  auto frame = reinterpret_cast<WIBFrame*>(static_cast<uint8_t*>(frag.get_data()) + i*sizeof(WIBFrame));
290  for (size_t j = 0; j < adc_vectors.size(); ++j)
291  {
292  adc_vectors[j].push_back(frame->get_channel(j));
293  }
294 
295  if (i == 0)
296  {
297  slot = frame->get_wib_header()->slot_no;
298  fiber = frame->get_wib_header()->fiber_no;
299  }
300  }
301  //std::cout << "slot, fiber: " << slot << ", " << fiber << std::endl;
302  for (size_t iChan = 0; iChan < 256; ++iChan)
303  {
304  const raw::RawDigit::ADCvector_t & v_adc = adc_vectors[iChan];
305  //std::cout << "Channel: " << iChan << " N ticks: " << v_adc.size() << " Timestamp: " << frag.get_trigger_timestamp() << std::endl;
306 
307  int offline_chan = channelMap->getOfflChanFromSlotFiberChan(slot, fiber, iChan);
308  if (offline_chan < 0) continue;
309  if (offline_chan > maxchan) continue;
310  raw::RDTimeStamp rd_ts(frag.get_trigger_timestamp(), offline_chan);
311  timestamps.push_back(rd_ts);
312 
313  float median = 0., sigma = 0.;
314  getMedianSigma(v_adc, median, sigma);
315  raw::RawDigit rd(offline_chan, v_adc.size(), v_adc);
316  rd.SetPedestal(median, sigma);
317  raw_digits.push_back(rd);
318  }
319 
320  }
321  H5Gclose(linkGroup);
322  }
323 
324 }
Collection of charge vs time digitized from a single readout channel.
Definition: RawDigit.h:69
std::vector< short > ADCvector_t
Type representing a (compressed) vector of ADC counts.
Definition: RawDigit.h:73
std::list< std::string > getMidLevelGroupNames(hid_t grp)
void getMedianSigma(const raw::RawDigit::ADCvector_t &v_adc, float &median, float &sigma)
void SetPedestal(float ped, float sigma=1.)
Set pedestal and its RMS (the latter is 0 by default)
Definition: RawDigit.cxx:68
int getOfflChanFromSlotFiberChan(int slot, int fiber, int chan)
hid_t getGroupFromPath(HDFFileInfoPtr &hdfFileInfoPtr, const std::string &path)
double median(sqlite3 *db, std::string const &table_name, std::string const &column_name)
Definition: statistics.cc:26
QTextStream & endl(QTextStream &s)
void VDColdboxDataInterface::getMedianSigma ( const raw::RawDigit::ADCvector_t v_adc,
float &  median,
float &  sigma 
)
private

Definition at line 327 of file VDColdboxDataInterface_tool.cc.

329  {
330  size_t asiz = v_adc.size();
331  int imed=0;
332  if (asiz == 0) {
333  median = 0;
334  sigma = 0;
335  }
336  else {
337  // the RMS includes tails from bad samples and signals and may not be the best RMS calc.
338 
339  imed = TMath::Median(asiz,v_adc.data()) + 0.01; // add an offset to make sure the floor gets the right integer
340  median = imed;
341  sigma = TMath::RMS(asiz,v_adc.data());
342 
343  // add in a correction suggested by David Adams, May 6, 2019
344 
345  size_t s1 = 0;
346  size_t sm = 0;
347  for (size_t i = 0; i < asiz; ++i) {
348  if (v_adc.at(i) < imed) s1++;
349  if (v_adc.at(i) == imed) sm++;
350  }
351  if (sm > 0) {
352  float mcorr = (-0.5 + (0.5*(float) asiz - (float) s1)/ ((float) sm) );
353  //if (std::abs(mcorr)>1.0) std::cout << "mcorr: " << mcorr << std::endl;
354  median += mcorr;
355  }
356  }
357 
358 }
double median(sqlite3 *db, std::string const &table_name, std::string const &column_name)
Definition: statistics.cc:26
void VDColdboxDataInterface::readFragmentsForEvent ( art::Event evt)
int VDColdboxDataInterface::retrieveData ( art::Event evt,
std::string  inputlabel,
std::vector< raw::RawDigit > &  raw_digits,
std::vector< raw::RDTimeStamp > &  rd_timestamps,
std::vector< raw::RDStatus > &  rdstatuses 
)
virtual

Implements PDSPTPCDataInterfaceParent.

Definition at line 167 of file VDColdboxDataInterface_tool.cc.

171  {
172  return 0;
173 }
int VDColdboxDataInterface::retrieveDataAPAListWithLabels ( art::Event evt,
std::string  inputlabel,
std::vector< raw::RawDigit > &  raw_digits,
std::vector< raw::RDTimeStamp > &  rd_timestamps,
std::vector< raw::RDStatus > &  rdstatuses,
std::vector< int > &  apalist 
)
virtual

Implements PDSPTPCDataInterfaceParent.

Definition at line 230 of file VDColdboxDataInterface_tool.cc.

236  {
237  return 0;
238 }
int VDColdboxDataInterface::retrieveDataForSpecifiedAPAs ( art::Event evt,
std::vector< raw::RawDigit > &  raw_digits,
std::vector< raw::RDTimeStamp > &  rd_timestamps,
std::vector< raw::RDStatus > &  rdstatuses,
std::vector< int > &  apalist 
)
virtual

Implements PDSPTPCDataInterfaceParent.

Definition at line 176 of file VDColdboxDataInterface_tool.cc.

181 {
182  using namespace dune::HDF5Utils;
183  auto infoHandle = evt.getHandle<raw::DUNEHDF5FileInfo>(fFileInfoLabel);
184  const std::string & toplevel_groupname = infoHandle->GetEventGroupName();
185  const std::string & file_name = infoHandle->GetFileName();
186  hid_t file_id = infoHandle->GetHDF5FileHandle();
187 
188  std::cout << "HDF5 FileName: " << file_name << std::endl;
189  std::cout << "Top-Level Group Name: " << toplevel_groupname << std::endl;
190 
191  //If the fcl file said to force open the file
192  //(i.e. because one is just running DataPrep), then open
193  //but only if we are on a new file -- identified by if the handle
194  //stored in the event is different
195  if (fForceOpen && (file_id != fPrevStoredHandle))
196  {
197  std::cout << "Opening" << std::endl;
198  fHDFFile = H5Fopen(file_name.data(), H5F_ACC_RDONLY, H5P_DEFAULT);
199  }//If the handle is the same, fHDFFile won't change
200  else if (!fForceOpen)
201  {
202  fHDFFile = file_id;
203  }
204  fPrevStoredHandle = file_id;
205 
206  hid_t the_group = getGroupFromPath(fHDFFile, toplevel_groupname);
207 
208  std::cout << "Retrieving Data for " << apalist.size() << " APA " << std::endl;
209 
210  // NOTE: The "apalist" that DataPrep hands to the method is always of size 1.
211  // Also "apalist" should technically hand you the current APA No. we are looking at but there is exception.
212  // CAUTION: This is only and only for VDColdBox.The reason is VDColdBox has only one APA/CRU.
213  for (const int & i : apalist)
214  {
215  int apano = i;
216  std::cout << "apano: " << i << std::endl;
217 
218  getFragmentsForEvent(the_group, raw_digits, rd_timestamps, apano, fMaxChan);
219 
220  //Currently putting in dummy values for the RD Statuses
221  rdstatuses.clear();
222  rdstatuses.emplace_back(false, false, 0);
223  }
224 
225  return 0;
226 }
Handle< PROD > getHandle(SelectorBase const &) const
Definition: DataViewImpl.h:382
std::string string
Definition: nybbler.cc:12
QCString file_name
void getFragmentsForEvent(hid_t the_group, RawDigits &raw_digits, RDTimeStamps &timestamps, int apano, int maxchan)
hid_t getGroupFromPath(HDFFileInfoPtr &hdfFileInfoPtr, const std::string &path)
QTextStream & endl(QTextStream &s)

Member Data Documentation

std::map<int,std::vector<std::string> > VDColdboxDataInterface::_input_labels_by_apa
private

Definition at line 56 of file VDColdboxDataInterface.h.

std::string VDColdboxDataInterface::fFileInfoLabel
private

Definition at line 69 of file VDColdboxDataInterface.h.

bool VDColdboxDataInterface::fForceOpen
private

Definition at line 68 of file VDColdboxDataInterface.h.

hid_t VDColdboxDataInterface::fHDFFile = -1
private

Definition at line 67 of file VDColdboxDataInterface.h.

int VDColdboxDataInterface::fMaxChan = 1000000
private

Definition at line 71 of file VDColdboxDataInterface.h.

hid_t VDColdboxDataInterface::fPrevStoredHandle = -1
private

Definition at line 66 of file VDColdboxDataInterface.h.

std::string VDColdboxDataInterface::logname = "VDColdboxDataInterface"
private

Definition at line 65 of file VDColdboxDataInterface.h.


The documentation for this class was generated from the following files: