LArRawInputDriverJP250L.cxx
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 /// \file LArRawInputDriverJP250L.cxx
3 /// \brief Source to convert JP250L files to LArSoft files
4 ///
5 /// \author eito@post.kek.jp, brebel@fnal.gov
6 ////////////////////////////////////////////////////////////////////////
7 
12 
21 
22 #include "TFile.h"
23 #include "TObject.h"
24 #include "TTree.h"
25 
26 #include <memory>
27 #include <vector>
28 
29 // ======================================================================
30 // JP250L conversion of raw DAQ file to ART format
31 
32 namespace lris {
33  // ======================================================================
34  // class c'tor/d'tor:
37  art::SourceHelper const &pm)
38  : principalMaker_(pm)
39  , m_current(0)
40  , m_data(0)
41  {
43  helper.reconstitutes<std::vector<raw::RawDigit>, art::InEvent>("daq");
44  helper.reconstitutes<sumdata::RunData, art::InRun> ("daq");
45  }
46 
47  // ======================================================================
49  {
50  delete [] m_data;
51  }
52 
53  // ======================================================================
55  art::FileBlock* &fb)
56  {
57  TFile m_f(name.c_str());
58  TTree* runTree = dynamic_cast<TTree*>(m_f.Get("runTree"));
59 
60  m_eventTree = dynamic_cast<TTree*>(m_f.Get("eventTree"));
61  m_eventTree->SetDirectory(0);
62  m_nEvent = m_eventTree->GetEntries();
63 
64  // run information
65  runTree->SetBranchAddress("runID",&m_runID);
66  runTree->SetBranchAddress("unixtime",&m_unixtime);
67  runTree->SetBranchAddress("nChannels",&m_nChannels);
68  runTree->SetBranchAddress("nSamples",&m_nSamples);
69  runTree->GetEntry(0);
70 
71  // have to add 1 to the runID because it can't be zero
72  // adding 1 to every runID ensures that all runIDs are bumped in the same way
73  m_runID += 1;
74 
75  const int nLength=(m_nSamples+4)*m_nChannels;
76  m_data = new unsigned short[nLength];
77  m_eventTree->SetBranchAddress("data",m_data);
78 
79  // Fill and return a new Fileblock.
80  // The string tells you what the version of the LArRawInputDriver is
81  fb = new art::FileBlock(art::FileFormatVersion(1, "LArRawInputJP250L 2013_01"),
82  name);
83 
84 
85  return;
86  }
87 
88  // ======================================================================
90  art::SubRunPrincipal* const & /* inSR */,
91  art::RunPrincipal* &outR,
92  art::SubRunPrincipal* &outSR,
93  art::EventPrincipal* &outE)
94  {
95 
96  if(m_current > m_nEvent) return false;
97 
98  m_eventTree->GetEntry(m_current);
99 
100  raw::DAQHeader daqHeader;
101  daqHeader.SetRun(m_runID);
102  daqHeader.SetTimeStamp(m_unixtime);
103  daqHeader.SetNChannels(m_nChannels);
104  // set the event number to start at 1 - ART likes to start numbering
105  // things from 1
106  daqHeader.SetEvent(m_current+1);
107 
108  // the following DAQHeader fields are not used by JP250L at this time
109  //daqHeader.SetStatus(1);
110  //daqHeader.SetFixedWord(h1.fixed);
111  //daqHeader.SetFileFormat(h1.format);
112  //daqHeader.SetSoftwareVersion(h1.software);
113  //daqHeader.SetSpareWord(h1.spare);
114 
115  // make unique_ptrs for the data products to store in the event
116  std::unique_ptr<raw::DAQHeader> daqcol( new raw::DAQHeader(daqHeader) );
117  std::unique_ptr<std::vector<raw::RawDigit> > rdcol( new std::vector<raw::RawDigit> );
118 
119  // loop over the signals and break them into
120  // one RawDigit for each channel
121  std::vector<short> adcVec(m_nSamples,0);
122  for(unsigned int n = 0; n < m_nChannels; ++n){
123  adcVec.clear();
124  for(unsigned int i = 0; i < m_nSamples; ++i){
125  adcVec.push_back(m_data[(m_nSamples+4)*n+(4+i)]);
126  }
127  rdcol->push_back(raw::RawDigit(n,m_nSamples,adcVec));
128  }
129 
130  art::RunNumber_t rn = daqHeader.GetRun();
132  art::EventNumber_t en = daqHeader.GetEvent();
133  art::Timestamp tstamp = daqHeader.GetTimeStamp();
134 
135  // Make the Run and SubRun principals
136  // this step is done once per run.
137  if (m_current < 1){
138  std::unique_ptr<sumdata::RunData> rundata(new sumdata::RunData("jpl250l") );
139  outR = principalMaker_.makeRunPrincipal(rn, tstamp);
141  sn,
142  tstamp);
143  art::put_product_in_principal(std::move(rundata), *outR, "daq");
144  }
145 
146 
147  // make the Event principal
149  sn,
150  en,
151  tstamp);
152 
153  // Put products in the event.
154  // first argument places the desired data product in the file
155  // second argument is event record to associate it to
156  // third argument is the label of the module storing the
157  // information. "daq" is standard in LArSoft
159  *outE,
160  "daq"); // Module label
162  *outE,
163  "daq"); // Module label
164 
165  // increment the entry to get out of the TTree for the next time through the loop
166  ++m_current;
167 
168  return true;
169  }
170 
171 }
static QCString name
Definition: declinfo.cpp:673
unsigned int m_unixtime
unix timestamp of the start of the run
unsigned int m_nEvent
number of triggers in the TTree
Collection of charge vs time digitized from a single readout channel.
Definition: RawDigit.h:69
void SetRun(unsigned short i)
Definition: DAQHeader.h:93
std::string string
Definition: nybbler.cc:12
void readFile(std::string const &name, art::FileBlock *&fb)
unsigned short m_runID
run ID, has to start from 1
void SetTimeStamp(time_t t)
Definition: DAQHeader.h:96
SubRunPrincipal * makeSubRunPrincipal(SubRunAuxiliary const &subRunAux) const
LArRawInputDriverJP250L(fhicl::ParameterSet const &pset, art::ProductRegistryHelper &helper, art::SourceHelper const &pm)
TTree * m_eventTree
TTree containing information from each trigger.
RunPrincipal * makeRunPrincipal(RunAuxiliary const &runAux) const
Definition: SourceHelper.cc:92
Source to convert JP250L files to LArSoft files.
TypeLabel const & reconstitutes(std::string const &modLabel, std::string const &instanceName={})
std::void_t< T > n
def move(depos, offset)
Definition: depos.py:107
IDNumber_t< Level::SubRun > SubRunNumber_t
Definition: IDNumber.h:119
Definition: SNSlice.h:7
void SetNChannels(uint32_t i)
Definition: DAQHeader.h:98
unsigned short * m_data
the ADC of each time sample for each channel
bool readNext(art::RunPrincipal *const &inR, art::SubRunPrincipal *const &inSR, art::RunPrincipal *&outR, art::SubRunPrincipal *&outSR, art::EventPrincipal *&outE)
unsigned short m_nChannels
number of channels in the detector
std::enable_if_t<!detail::range_sets_supported(P::branch_type)> put_product_in_principal(std::unique_ptr< T > &&product, P &principal, std::string const &module_label, std::string const &instance_name={})
unsigned short m_nSamples
number of time samples per channel
Conversion of binary data to root files.
unsigned short GetRun() const
Definition: DAQHeader.h:103
art::SourceHelper const & principalMaker_
IDNumber_t< Level::Event > EventNumber_t
Definition: IDNumber.h:118
unsigned short GetEvent() const
Definition: DAQHeader.h:105
unsigned int m_current
current entry in the TTree
EventPrincipal * makeEventPrincipal(EventAuxiliary const &eventAux, std::unique_ptr< History > &&history) const
void SetEvent(unsigned short i)
Definition: DAQHeader.h:95
time_t GetTimeStamp() const
Definition: DAQHeader.h:106
IDNumber_t< Level::Run > RunNumber_t
Definition: IDNumber.h:120