All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
DaqDecoderICARUSPMT_module.cc
Go to the documentation of this file.
1 /**
2  * @file DaqDecoderICARUSPMT_module.cc
3  * @brief Produces `raw::OpDetWaveform` from V1730 artDAQ data fragments.
4  * @authors Andrea Scarpelli, Gianluca Petrillo (petrillo@slac.stanford.edu)
5  * @date May 21, 2021
6  *
7  */
8 
9 
10 // ICARUS/SBN libraries
16 #include "icarusalg/Utilities/BinaryDumpUtils.h" // icarus::ns::util::bin()
18 
19 #include "sbnobj/Common/PMT/Data/PMTconfiguration.h" // sbn::PMTconfiguration
20 #include "sbndaq-artdaq-core/Overlays/Common/CAENV1730Fragment.hh"
21 #include "sbndaq-artdaq-core/Overlays/FragmentType.hh" // sbndaq::FragmentType
22 
23 // LArSoft libraries
25 
28 #include "lardataalg/Utilities/quantities/spacetime.h" // nanoseconds
29 #include "lardataalg/Utilities/intervals_fhicl.h" // for nanoseconds in FHiCL
32 
35 
36 // artDAQ
37 #include "artdaq-core/Data/ContainerFragment.hh"
38 #include "artdaq-core/Data/Fragment.hh"
39 
40 // framework libraries
41 #include "art_root_io/TFileService.h"
42 #include "art/Framework/Services/Registry/ServiceHandle.h"
43 #include "art/Framework/Core/EDProducer.h"
44 #include "art/Framework/Core/ModuleMacros.h"
45 #include "art/Framework/Principal/Run.h"
46 #include "art/Framework/Principal/Event.h"
47 #include "art/Framework/Principal/Handle.h"
48 #include "canvas/Persistency/Provenance/EventID.h"
49 #include "canvas/Persistency/Provenance/Timestamp.h"
50 #include "canvas/Utilities/InputTag.h"
51 #include "messagefacility/MessageLogger/MessageLogger.h"
52 #include "fhiclcpp/types/TableAs.h"
53 #include "fhiclcpp/types/OptionalSequence.h"
54 #include "fhiclcpp/types/Sequence.h"
55 #include "fhiclcpp/types/OptionalAtom.h"
56 #include "fhiclcpp/types/Atom.h"
57 #include "cetlib_except/exception.h"
58 
59 // ROOT libraries
60 #include "TTree.h"
61 
62 // C/C++ standard libraries
63 #include <memory>
64 #include <ostream>
65 #include <unordered_map>
66 #include <map>
67 #include <set>
68 #include <vector>
69 #include <string>
70 #include <optional>
71 #include <cassert>
72 
73 
74 //------------------------------------------------------------------------------
75 using namespace util::quantities::time_literals;
77 
78 //------------------------------------------------------------------------------
79 namespace icarus { class DaqDecoderICARUSPMT; }
80 /**
81  * @brief Produces `raw::OpDetWaveform` from V1730 artDAQ data fragments.
82  *
83  * The module can read fragments from CAEN V1730 readout boards delivered by
84  * artDAQ. It produces optical detector waveforms in LArSoft standards.
85  *
86  * This decoder must support both a off-line mode (for storage and downstream
87  * processing) and a on-line mode (for monitoring).
88  * In particular, the on-line workflow is such that it may not be possible to
89  * access the FHiCL configuration of the job and therefore the PMT configuration
90  * data (see `icarus::PMTconfigurationExtraction` module).
91  *
92  *
93  * Configuration
94  * --------------
95  *
96  * The set of supported parameters can be seen on command line by running
97  * `lar --print-description DaqDecoderICARUSPMT`.
98  *
99  * Description of the configuration parameters:
100  * * `FragmentsLabels` (list of input tags): a list of possible input data
101  * products. A valid input data product contain a list of fragments from
102  * V1730. All input tags are tried. If only one is available, its content
103  * is used for decoding; otherwise, an exception is thrown.
104  * * `DiagnosticOutput` (flag, default: `false`): enables additional console
105  * output, including dumping of the fragments (that is huge output).
106  * * `PMTconfigTag` (data product tag, optional): if specified, the pre-trigger
107  * buffer duration is read from there; although optional, it is strongly
108  * recommended that this information be provided, since it is essential for
109  * the correct timing of the PMT waveforms (see
110  * @ref icarus_PMTDecoder_timestamps "the discussion on time stamps below").
111  * * `BoardSetup` (list of board setup information): each entry specifies some
112  * information about a specific readout board; the boards are identified by
113  * their name; if a board is found in input that has no setup information,
114  * some time corrections are not applied (see
115  * @ref icarus_PMTDecoder_timestamps "the discussion on time stamps below").
116  * Each entry is in the form of a table:
117  * * `Name` (string, mandatory): the name of the board
118  * (e.g. `"icaruspmtwwtop01"`); this is use to match the setup
119  * information to a fragment ID in the PMT configuration.
120  * * `FragmentID` (integral, optional): if specified, allows the corrections
121  * using setup information to be applied even when no PMT configuration is
122  * provided (if neither PMT configuration nor setup information including
123  * `FragmentID` is available, no time correction is applied).
124  * * `TriggerDelay` (nanoseconds, default: 0 ns): measured delay from the
125  * primitive trigger time to the execution of the PMT trigger; specify
126  * the unit! (e.g. `"43 ns"`).
127  * * `SpecialChannels` (list of configurations): each entry described
128  * special features and settings of a specific channel, identified by
129  * its V1730B board channel number (`0` to `15`). The supported keys are:
130  * * `ChannelIndex` (integer, mandatory) the index of the channel in
131  * the board, from `0` to `15`; note that for this module to decode
132  * this channel index, the channel itself must be enabled on the
133  * readout board.
134  * * `Skip` (flag, optional): if specified, tells whether to always skip
135  * this channel or whether always store it. If not specified, a
136  * channel will be always stored unless it is not associated to any
137  * channel number (see `Channel` below). Note that if the channel is
138  * not enabled in the readout board, this option has no effect.
139  * * `OnlyOnGlobalTrigger` (flag, default: `false`): if set, waveforms
140  * are saved from this channel only when it includes the global
141  * trigger time; if no trigger time is available, such channels will
142  * never be saved.
143  * * `MinSpan` (integer, default: `0`): waveforms which have a span
144  * (highest sample minus lowest sample) smaller than this limit are
145  * not saved (unless `Skip` is set to `false`).
146  * * `Channel` (integral, optional): if specified, it overrides the
147  * channel number used for the waveforms on this channel; if not
148  * specified, the channel number is obtained from the channel mapping
149  * database. Note that insisting to save a waveform without an
150  * assigned channel number (from database or from this configuration)
151  * will trigger an exception.
152  * * `InstanceName` (string, default: empty): the waveform will be
153  * stored in a collection of `raw::OpDetWaveform` with the specified
154  * instance name; the default is an empty name, which is the standard
155  * data product where all the normal waveforms are saved.
156  * * `RequireKnownBoards` (flag, default: `true`): if set, the readout boards
157  * in input must each have a setup configuration (`BoardSetup`) *and* must
158  * be present in the PMT DAQ configuration, or an exception is thrown;
159  * if not set, readout boards in input will be processed even when their
160  * setup or DAQ configuration is not known, at the cost of fewer corrections
161  * on the timestamps (which should then be considered unreliable).
162  * * `RequireBoardConfig` (flag, default: `true`): if set, the readout boards
163  * which have a setup (`BoardSetup`) are required to be included in the DAQ
164  * configuration of the input file, or an exception is thrown; if not set,
165  * missing readout boards are unnoticed.
166  * * `TriggerTag` (data product tag): tag for the information
167  * (`sbn::ExtraTriggerInfo`, produced by trigger decoding); if not
168  * specified, the _art_ event timestamp will be used as trigger time (*not*
169  * recommended).
170  * * `TTTresetEverySecond` (optional): if set, the decoder will take advantage
171  * of the assumption that the Trigger Time Tag of all PMT readout boards is
172  * synchronised with the global trigger time and reset at every change of
173  * second of the timescale of the latter; this is currently the only
174  * implementation supporting multiple PMT readout windows on the same board;
175  * if this option is set to `false`, all PMT readout boards are assumed to
176  * have been triggered at the time of the global trigger. By default, this
177  * option is set to `true` unless `TriggerTag` is specified empty.
178  * * `DataTrees` (list of strings, default: none): list of data trees to be
179  * produced; if none (default), then `TFileService` is not required.
180  * * `SkipWaveforms` (flag, default: `false`) if set, waveforms won't be
181  * produced; this is intended as a debugging option for jobs where only the
182  * `DataTrees` are desired.
183  * * `DropRawDataAfterUse` (flag, default: `true`): at the end of processing,
184  * the framework will be asked to remove the PMT data fragment from memory.
185  * Set this to `false` in the unlikely case where raw PMT fragments are
186  * still needed after decoding.
187  * * `LogCategory` (string, default: `DaqDecoderICARUSPMT`): name of the message
188  * facility category where the output is sent.
189  *
190  *
191  * Requirements
192  * -------------
193  *
194  * Services required include:
195  *
196  * * `IICARUSChannelMap` for the association of fragments to LArSoft channel ID;
197  * * `DetectorClocksService` for the correct decoding of the time stamps
198  * (always required, even when dumbed-down timestamp decoding is requested);
199  * * `TFileService` only if the production of trees or plots is requested.
200  *
201  *
202  * Waveform time stamp
203  * --------------------
204  *
205  * @anchor icarus_PMTDecoder_timestamps
206  *
207  * All waveforms on the same readout board fragment share the same timestamp.
208  * The same readout board can produce different fragments at different times
209  * within an event in which case each fragment will be independently assigned
210  * a timstamp .
211  *
212  * The time stamp of the waveform is defined as the time when the first sample
213  * of the waveform started (that is, if the sample represent the value of the
214  * signal in an interval of 2 ns, the time stamp is pointing at the beginning
215  * of those 2 ns). Whether we can honour that definition, though, is a different
216  * matter.
217  * The representation of the time stamp is in the
218  * @ref DetectorClocksElectronicsTime "electronics time scale".
219  *
220  * Waveforms all originate from a trigger primitive signal being sent to the
221  * readout boards by NI7820 FPGA to start the acquisition of the waveform.
222  * One of these primitives matches the global event trigger, but the
223  * corresponding fragment is not treated in any special way.
224  * Every delay between when that signal is emitted and when the PMT trigger is
225  * executed shifts the timestamp of the waveform backward.
226  *
227  * We assign the the time stamp of the waveforms as follow:
228  * 1. the base time is the global trigger time; the global trigger time is read
229  * from the trigger data, where it is stored as an absolute time in TAI scale;
230  * the global trigger time itself effectively defines the electronics time
231  * scale used within a _art_ event, so its representation is a fixed number
232  * that is configured in LArSoft and can be accessed with
233  * `DetectorClocksData::TriggerTime()`;
234  * 2. each V1730 data fragment comes with a time tag ("TTT") describing the time
235  * of the end of the readout buffer (unknown whether the start or the end of
236  * the tick of the last sample). This tag is a counter internal to the V1730
237  * board, incremented every 8 ns. The counter is reset at the beginning of
238  * every "new" TAI second, and it can thus be easily transformed into a time
239  * in the same TAI scale as the global trigger; the difference between TTT
240  * and the global trigger pins down the end of the readout buffer in the
241  * electronics time scale;
242  * 3. a delay is subtracted to the timestamp, which encompasses all fixed delays
243  * occurring in the trigger signal transportation; the most prominent is the
244  * delay occurring between the start of the "new" TAI second and when the
245  * TTT reset signal reaches and is honoured by the readout board.
246  * This value must be independently measured and provided to this decoder via
247  * configuration as setup information (`TTTresetDelay`); if not present in the
248  * setup, this delay is not added;
249  * 4. finally, the time of the beginning of the waveform, that is the target
250  * for `raw::OpDetWaveform` timestamp, is obtained from the time of its end
251  * by simply subtracting the readout buffer length.
252  *
253  * A special decoding mode, hoped to be just historical by now, does not rely
254  * on the global trigger time. This mode can be safely used only when the (only)
255  * data fragment was triggered by the global trigger (e.g. in the simplest
256  * minimum/zero bias trigger). It has the advantage that it does not use the TTT
257  * information. It can be enabled explicitly by setting the option
258  * `TTTresetEverySecond` to `true`, or by removing the specification of the
259  * trigger time data product tag (`TriggerTag`).
260  * 1. The trigger primitive time is assumed to be the global trigger time too,
261  * so that the trigger primitive time in electronics time also matches the
262  * global trigger time.
263  * 2. upon receiving the trigger primitive signal, the readout board will keep
264  * some of the samples already digitized, in what we call pre-trigger buffer;
265  * the size of this buffer is a fixed number of samples which is specified in
266  * DAQ as a fraction of the complete buffer that is _post-trigger_; this
267  * amount, converted in time, is subtracted to the trigger time to point back
268  * to the beginning of the waveform instead that to the trigger primitive
269  * time. The necessary information is read from the PMT configuration
270  * (`PMTconfigTag`); if no configuration is available, this offset is not
271  * subtracted; note that this is a major shift (typically, a few microseconds)
272  * that should be always included. This step is equivalent to the step (4) of
273  * the regular mode, where here the adjustment starts from the trigger
274  * primitive time instead than from the end-of-buffer time.
275  * 3. a delay is subtracted to the timestamp, which encompasses all fixed delays
276  * occurring in the trigger signal transportation; one component of it is the
277  * relative delay in the propagation of the trigger primitive signal from a
278  * board to the next (in fact, every three boards have the trigger signal in
279  * daisy chain). This value must be independently measured and provided to
280  * this decoder via configuration as setup information (`TriggerDelay`); if
281  * not present in the setup, this delay is not added.
282  * Note that the particular contribution of the daisy chain to the delay does
283  * not need to be explicitly taken into account in the main mode, because the
284  * TTT reset is independent and not daisy-chained, so that the TTT times
285  * are all synchronized and when the primitive trigger arrives (via daisy
286  * chain) the TTT value at that instant is already including the delay.
287  *
288  *
289  *
290  * Data trees
291  * -----------
292  *
293  * The module supports the following ROOT trees production on demand:
294  *
295  * * `PMTfragments`: data pertaining a single fragment; each entry is about a
296  * single fragment, and it includes full event ID, event time stamp (from
297  * _art_, i.e. the one assigned to the event by artDAQ), the ID of the
298  * fragment the entry is describing, and then the rest of the data, including:
299  * * `TTT` (64-bit integer): (Extended) Trigger Time Tag value, in readout
300  * board ticks (each worth 8 ns) from the last board reset;
301  * currently the value includes the 31-bit counter and in addition the
302  * overflow bit as MSB; the overflow bit is set by the readout board
303  * the first time the counter passes its limit (2^31) and wraps, and never
304  * cleared until the next board reset.
305  * While the tree has room for the 48-bit version (ETTT), the rest of the
306  * decoder does not yet support it.
307  * * `trigger` (64-bit signed integer): global trigger time (from the
308  * trigger decoder), in nanoseconds from The Epoch.
309  * * `triggerSec`, `triggerNS` (32-bit integer each): same time as `trigger`
310  * branch, split into second and nanosecond components.
311  * * `relBeamGateNS`, (32-bit integer): beam gate time opening relative to
312  * the trigger time, in nanoseconds; it may be affected by rounding.
313  * branch, split into second and nanosecond components.
314  * * `fragTime` (64-bit signed integer), `fragTimeSec` (32-bit signed
315  * integer): the timestamp of the PMT fragment, assigned by the board
316  * reader constructing the fragment.
317  * * `fragCount` (unsigned integer): the counter of the fragment within a
318  * readout board (its "event counter").
319  * * `waveformTime` (double): time assigned to the waveforms from this
320  * fragment (electronics time scale).
321  * * `waveformSize` (unsigned integer): number of ticks for the waveforms
322  * from this fragment.
323  * * `triggerBits` (unsigned integer): bits from the `raw::Trigger`.
324  * * `gateCount` (unsigned integer): number of this gate from run start.
325  * * `onGlobalTrigger` (boolean): whether the waveform covers the nominal
326  * trigger time (which should be equivalent to whether the fragment was
327  * triggered by the global trigger).
328  *
329  *
330  * Technical notes
331  * ----------------
332  *
333  * In order to correctly reconstruct the time stamp, this module needs several
334  * pieces of information.
335  * These include the size of the pre-trigger buffer, which is set by the readout
336  * board configuration, and the delay either between the start of the "new"
337  * second and the execution of the TTT reset, or between the global trigger and
338  * the time that trigger is received and acted upon in the readout board,
339  * depending of the time stamping mode; in both cases, the delays need to be
340  * measured.
341  * The first category of information, from readout board configuration, are read
342  * from the input file (`sbn::PMTconfiguration`), while the second category
343  * needs to be specified in the module FHiCL configuration.
344  *
345  * PMT configuration is optional, in the sense that it can be omitted; in that
346  * case, some standard values will be used for it. This kind of setup may be
347  * good for monitoring, but it does not meet the requirements for physics
348  * analyses.
349  * For a board to be served, an entry of that board must be present in the
350  * module configuration (`BoardSetup`). It is an error for a fragment in input
351  * not to have an entry for the corresponding board setup.
352  *
353  * The module extracts the needed information and matches it into a
354  * sort-of-database keyed by fragment ID, so that it can be quickly applied
355  * when decoding a fragment. The matching is performed by board name.
356  *
357  *
358  * ### Proto-waveforms
359  *
360  * The module has grown complex enough that some decisions need to be taken much
361  * later after a waveform is created, but using information available only
362  * during the creation. For example, a waveform can be ultimately routed to a
363  * different data product depending on which readout board it comes from,
364  * but that information is lost by the time the waveform needs to be written.
365  * Rather than attempting to recreate that information, it is saved when still
366  * available and travels together with the waveform itself in a data structure
367  * called "proto-waveform". The extra-information is eventually discarded when
368  * putting the actual waveform into the _art_ event.
369  *
370  *
371  * ### Processing pipeline
372  *
373  * Processing happens according to the following structure (in `produce()`):
374  * 1. pre-processing: currently nothing
375  * 2. processing of each board data independently: at this level, all the
376  * buffers from the 16 channels of a single board are processed together
377  * (`processBoardFragments()`)
378  * 1. the configuration and parameters specific to this board are fetched
379  * 2. each data fragment is processed independently: at this level, data
380  * from all 16 channels _at a given time_ are processed together,
381  * producing up to 16 proto-waveforms
382  * 3. merging of contiguous waveforms is performed
383  * 3. post-processing of proto-waveforms:
384  * * sorting by time (as opposed as roughly by channel, as they come)
385  * 4. conversion to data products and output
386  *
387  *
388  *
389  * Glossary
390  * ---------
391  *
392  * * **setup**, **[PMT] configuration**: this is jargon specific to this module.
393  * Information about a readout board can come from two sources: the "setup"
394  * is information included in the `BoardSetup` configuration list of this
395  * module; the "PMT configuration" is information included in the DAQ
396  * configuration that is delivered via `PMTconfigTag`.
397  * * **TAI** (International Atomic Time): a time standard defining a universal
398  * time with a precision higher than it will ever matter for ICARUS.
399  * * **TTT**: trigger time tag, from the V1730 event record (31 bits); may be:
400  * * **ETTT**: extended trigger time tag, from the V1730 event record (48 bits).
401  * * **trigger delay**: time point when a V1730 board processes a (PMT) trigger
402  * signal (and increments the TTT register) with respect to the time of the
403  * time stamp of the (SPEXi) global trigger that acquired the event.
404  *
405  */
406 class icarus::DaqDecoderICARUSPMT: public art::EDProducer {
407 
408  // --- BEGIN -- some debugging tree declarations -----------------------------
409 
410  /// Enumerate the supported data trees.
411  enum class DataTrees: std::size_t {
412  Fragments, ///< Information about fragments
413  N ///< Counter.
414  };
415  using TreeNameList_t
416  = std::array<std::string, static_cast<std::size_t>(DataTrees::N)>;
417  static TreeNameList_t const TreeNames;
418 
419  /// Returns a string with all supported tree names.
420  static std::string listTreeNames(std::string const& sep = "\n");
421 
422  // --- END ---- some debugging tree declarations -----------------------------
423 
424  public:
425 
426  // --- BEGIN -- public data types --------------------------------------------
427 
429 
430  /// Data structure for trigger time.
432  static_assert(sizeof(int) >= 4U);
433  struct Split_t {
434  int seconds = std::numeric_limits<int>::min(); ///< The ongoing second.
435  /// Nanoseconds from the start of current second.
436  unsigned int nanoseconds = std::numeric_limits<unsigned int>::max();
437  }; // Split_t
438 
439  /// Trigger time in nanoseconds from The Epoch.
440  long long int time = std::numeric_limits<long long int>::min();
441  /// Trigger time in nanoseconds from The Epoch (in components).
443 
444  constexpr SplitTimestamp_t() = default;
445  constexpr SplitTimestamp_t(int sec, unsigned int ns);
446  constexpr SplitTimestamp_t(long long int triggerTime);
447  }; // SplitTimestamp_t
448 
449  // --- END ---- public data types --------------------------------------------
450 
451 
452  // --- BEGIN -- FHiCL configuration ------------------------------------------
453 
454  /// Configuration of the V1730 readout board setup.
456 
458 
459  fhicl::Atom<unsigned short int> ChannelIndex {
460  fhicl::Name("ChannelIndex"),
461  fhicl::Comment
462  ("index of the channel on the board these settings pertain")
463  // mandatory
464  };
465 
466  fhicl::OptionalAtom<bool> Skip {
467  fhicl::Name("Skip"),
468  fhicl::Comment(
469  "set to true to force skipping this channel, false to force saving it"
470  )
471  };
472 
473  fhicl::Atom<bool> OnlyOnGlobalTrigger {
474  fhicl::Name("OnlyOnGlobalTrigger"),
475  fhicl::Comment
476  ("save this channel only if its data includes global trigger time"),
477  false // default
478  };
479 
480  fhicl::Atom<std::uint16_t> MinSpan {
481  fhicl::Name("MinSpan"),
482  fhicl::Comment(
483  "discard this channel if its span (maximum minus minimum)"
484  " is smallerthan this"
485  ),
486  0 // default
487  };
488 
489  fhicl::Atom<raw::Channel_t> Channel {
490  fhicl::Name("Channel"),
491  fhicl::Comment("Off-line channel ID associated to this board channel"),
493  };
494 
495  fhicl::Atom<std::string> InstanceName {
496  fhicl::Name("InstanceName"),
497  fhicl::Comment
498  ("name of the data product instance where to add this channel"),
499  "" // default
500  };
501 
502  }; // struct ChannelSetupConfig
503 
504 
505  fhicl::Atom<std::string> Name {
506  fhicl::Name("Name"),
507  fhicl::Comment("board name, as specified in the DAQ configuration")
508  };
509 
510  fhicl::OptionalAtom<unsigned int> FragmentID {
511  fhicl::Name("FragmentID"),
512  fhicl::Comment("ID of the fragments associated with the board")
513  };
514 
515  fhicl::Atom<nanoseconds> TriggerDelay {
516  fhicl::Name("TriggerDelay"),
517  fhicl::Comment
518  ("from delay from the trigger timestamp to the PMT trigger [ns]"),
519  0_ns
520  };
521 
522  fhicl::Atom<nanoseconds> TTTresetDelay {
523  fhicl::Name("TTTresetDelay"),
524  fhicl::Comment
525  ("assume that V1730 counter (Trigger Time Tag) is reset every second"),
526  0_ns
527  };
528 
529  fhicl::OptionalSequence<fhicl::Table<ChannelSetupConfig>> SpecialChannels {
530  fhicl::Name("SpecialChannels"),
531  fhicl::Comment("special settings for selected channels on the board")
532  };
533 
534  }; // struct BoardSetupConfig
535 
536 
537  /// Main module configuration.
538  struct Config {
539 
540  using Name = fhicl::Name;
541  using Comment = fhicl::Comment;
542 
543  fhicl::Sequence<art::InputTag> FragmentsLabels {
544  Name("FragmentsLabels"),
545  Comment("data product candidates with the PMT fragments from DAQ"),
546  std::vector<art::InputTag>{ "daq:CAENV1730", "daq:ContainerCAENV1730" }
547  };
548 
549  fhicl::Atom<bool> SurviveExceptions {
550  Name("SurviveExceptions"),
551  Comment
552  ("when the decoding module throws an exception, print a message and move on"),
553  true // default
554  };
555 
556  fhicl::Atom<bool> DiagnosticOutput {
557  Name("DiagnosticOutput"),
558  Comment("enable additional console output"),
559  false // default
560  };
561 
562  fhicl::Atom<bool> PacketDump {
563  Name("PacketDump"),
564  Comment("enable dump of the whole V1730 data (huge)"),
565  false // default
566  };
567 
568  fhicl::Atom<bool> RequireKnownBoards {
569  Name("RequireKnownBoards"),
570  Comment
571  ("all readout boards in input must be known (setup+PMT configuration)"),
572  true
573  };
574 
575  fhicl::Atom<bool> RequireBoardConfig {
576  Name("RequireBoardConfig"),
577  Comment
578  ("all readout boards in setup must have a matching PMT configuration"),
579  true
580  };
581 
582  fhicl::OptionalAtom<art::InputTag> PMTconfigTag {
583  Name("PMTconfigTag"),
584  Comment("input tag for the PMT readout board configuration information")
585  };
586 
587  fhicl::Sequence
588  <fhicl::TableAs<daq::details::BoardSetup_t, BoardSetupConfig>>
589  BoardSetup {
590  Name("BoardSetup"),
591  Comment("list of the setup settings for all relevant V1730 boards")
592  };
593 
594  fhicl::OptionalAtom<art::InputTag> TriggerTag {
595  Name("TriggerTag"),
596  Comment("input tag for the global trigger object (sbn::ExtraTriggerInfo)")
597  };
598 
599  fhicl::OptionalAtom<bool> TTTresetEverySecond {
600  Name("TTTresetEverySecond"),
601  Comment
602  ("assume that V1730 counter (Trigger Time Tag) is reset every second")
603  };
604 
605  fhicl::Sequence<std::string> DataTrees {
606  fhicl::Name("DataTrees"),
607  fhicl::Comment
608  ("produces the specified ROOT trees (" + listTreeNames(",") + ")"),
609  std::vector<std::string>{} // default
610  };
611 
612  fhicl::Atom<bool> SkipWaveforms {
613  Name("SkipWaveforms"),
614  Comment("do not decode and produce waveforms"),
615  false // default
616  };
617 
618  fhicl::Atom<bool> DropRawDataAfterUse {
619  Name("DropRawDataAfterUse"),
620  Comment("drop PMT data fragments from memory after use"),
621  true // default
622  };
623 
624  fhicl::Atom<std::string> LogCategory {
625  Name("LogCategory"),
626  Comment("name of the category for message stream"),
627  "PMTDecoder" // default
628  };
629 
630  }; // Config
631 
632  using Parameters = art::EDProducer::Table<Config>;
633 
634 
635  static constexpr electronics_time NoTimestamp
636  = std::numeric_limits<electronics_time>::min();
637 
638 
639  // --- END ---- FHiCL configuration ------------------------------------------
640 
641 
642  /// Constructor.
643  explicit DaqDecoderICARUSPMT(Parameters const& params);
644 
645  /// On a new run: cache PMT configuration information.
646  void beginRun(art::Run& run) override;
647 
648  /// Processes the event.
649  void produce(art::Event& event) override;
650 
651  /// Prints a end-of-job message.
652  void endJob() override;
653 
654 
655  private:
656 
657  /// Type of setup of all channels in a readout board.
659 
660  /// Collection of useful information from fragment data.
661  struct FragmentInfo_t {
662  artdaq::Fragment::fragment_id_t fragmentID
663  = std::numeric_limits<artdaq::Fragment::fragment_id_t>::max();
664  artdaq::Fragment::timestamp_t fragmentTimestamp;
665  unsigned int eventCounter = 0U;
666  std::uint32_t TTT = 0U;
667  std::uint16_t enabledChannels = 0U;
668  std::size_t nSamplesPerChannel = 0U;
669  std::uint16_t const* data = nullptr;
670  }; // FragmentInfo_t
671 
672  /// Information used in decoding from a board.
674  static AllChannelSetup_t const DefaultChannelSetup; // default-initialized
675 
676  std::string const name;
681  AllChannelSetup_t const* specialChannelSetup = nullptr;
682 
684  { return specialChannelSetup? *specialChannelSetup: DefaultChannelSetup; }
685 
686  }; // NeededBoardInfo_t
687 
688  /// Collection of information about the global trigger.
689  struct TriggerInfo_t {
690  SplitTimestamp_t time; ///< Time of the trigger (absolute).
691  long int relBeamGateTime; ///< Time of beam gate relative to trigger [ns].
692  sbn::triggerSourceMask bits; ///< Trigger bits.
693  unsigned int gateCount = 0U; ///< Gate number from the beginning of run.
694  }; // TriggerInfo_t
695 
696  /// All the information collected about a waveform (with the waveform itself).
698 
699  raw::OpDetWaveform waveform; ///< The complete waveform.
700 
701  /// Pointer to the settings for the channel this waveform belongs to.
702  daq::details::BoardSetup_t::ChannelSetup_t const* channelSetup = nullptr;
703 
704  /// Whether the waveform includes the global trigger time.
705  bool onGlobal = false;
706 
707  ///< Lowest sample in the original waveform.
708  std::uint16_t minSample = std::numeric_limits<std::uint16_t>::max();
709  ///< Highest sample in the original waveform.
710  std::uint16_t maxSample = std::numeric_limits<std::uint16_t>::max();
711 
712  ///< Returns the span of the waveform (cached, not computed anew!).
713  std::uint16_t span() const
714  { return (maxSample > minSample)? (maxSample - minSample): 0; }
715 
716  /// Ordering: the same as the contained waveform.
717  bool operator< (ProtoWaveform_t const& than) const
718  { return waveform < than.waveform; }
719 
720  }; // struct ProtoWaveform_t
721 
722  // --- BEGIN -- Configuration parameters -------------------------------------
723 
724  ///< List of candidate data products with artDAQ data fragments.
725  std::vector<art::InputTag> const fInputTags;
726 
727  bool const fSurviveExceptions; ///< Whether to "ignore" errors.
728 
729  /// If true will spew endless messages to output.
730  bool const fDiagnosticOutput;
731 
732  bool const fPacketDump; ///< Dump V1730 data.
733 
734  /// Whether info on all input boards is required.
736 
737  /// Whether setup info on all boards is required.
739 
740  /// Input tag of the PMT configuration.
741  std::optional<art::InputTag> const fPMTconfigTag;
742 
743  /// Input tag of the global trigger.
744  std::optional<art::InputTag> const fTriggerTag;
745 
746  bool const fTTTresetEverySecond; ///< Whether V1730 TTT is reset every second.
747 
748  /// All board setup settings.
749  std::vector<daq::details::BoardSetup_t> const fBoardSetup;
750 
751  bool const fSkipWaveforms; ///< Whether to skip waveform decoding.
752 
753  /// Clear fragment data product cache after use.
755 
756  std::string const fLogCategory; ///< Message facility category.
757 
758  // --- END ---- Configuration parameters -------------------------------------
759 
760 
761  // --- BEGIN -- Services -----------------------------------------------------
762 
763  /// Interface to LArSoft configuration for detector timing.
765 
766  /// Fragment/channel mapping database.
768 
769  // --- END ---- Services -----------------------------------------------------
770 
771 
772  // --- BEGIN -- Cached values ------------------------------------------------
773 
774  /// Duration of the optical detector readout sampling tick (i.e. 2 ns; hush!).
776 
777  /// Trigger time as reported by `detinfo::DetectorClocks` service.
779 
780  // --- END ---- Cached values ------------------------------------------------
781 
782 
783  // --- BEGIN -- Per-run data cache -------------------------------------------
784 
785  /// Find the information on a readout boards by fragment ID.
786  std::optional<daq::details::BoardInfoLookup> fBoardInfoLookup;
787 
788  // --- END -- Per-run data cache ---------------------------------------------
789 
790 
791  unsigned int fNFailures = 0U; ///< Number of event failures encountered.
792 
793 
794  // --- BEGIN -- PMT readout configuration ------------------------------------
795 
796  /// Returns whether PMT configuration information is expected to be available.
797  bool hasPMTconfiguration() const { return fPMTconfigTag.has_value(); }
798 
799  /// Updates the PMT configuration cache. How? Dunno. Placeholder.
800  bool UpdatePMTConfiguration(sbn::PMTconfiguration const* PMTconfig);
801 
802 
803  /**
804  * @brief Returns a lookup object with board setup and configuration info.
805  * @param PMTconfig the PMT configuration, if available
806  * @return an object working like lookup table for all fragment information
807  *
808  * This method merges the setup information from the module configuration with
809  * the PMT configuration specified in the argument, and returns an object
810  * that can look up all the information as a single record, with the
811  * fragment ID as key. In addition, a few intermediate quantities ("facts",
812  * see `BoardFacts_t`) are computed and stored in this object.
813  *
814  * If a fragment ID is missing, it means that no PMT configuration was
815  * provided and that the setup information did not include a fragment ID.
816  * If some information (configuration or setup) is missing, the "facts"
817  * depending on the missing information will have default values.
818  */
819  daq::details::BoardInfoLookup matchBoardConfigurationAndSetup
820  (sbn::PMTconfiguration const* PMTconfig) const;
821 
822  /// Puts together all the needed information for a board.
823  NeededBoardInfo_t fetchNeededBoardInfo(
825  unsigned int fragmentID
826  ) const;
827 
828  /// Extracts the Trigger Time Tag (31+1 bits) value from the fragment
829  static unsigned int extractTriggerTimeTag(artdaq::Fragment const& fragment);
830 
831 
832  // --- END -- PMT readout configuration --------------------------------------
833 
834 
835  /**
836  * @brief Returns the count of set bits for each set bit.
837  * @tparam NBits the number of bits to test
838  * @param value the bit mask to be analyzed
839  * @return a pair: `first` an array with the `NBits` values,
840  * `second` the number of bits set
841  *
842  * Much better to go with an example.
843  * Let's process `setBitIndices<16U>(0x6701)`: `value` is `0110 0111 0000 0001`,
844  * and we request the `16` least significant bits (`NBits`).
845  * There are six bit that are set, in position 0, 8, 9, 10, 13 and 14.
846  * The returned value includes the number of set bits (6) as the `second`
847  * element of the pair, an as `first` an array of 16 indices, one for each bit
848  * position.
849  * Its values are `0` at `[0]`, `1` at `[8]`, `2` at `[9]`, `3` at `[10]`,
850  * `4` at `[13]` and `5` at `[14]`. That is, each value represents the number
851  * of that bit in the list of set bits. All the other indices, associated to
852  * the bits which are not set, are assigned a value that is equal or larger than
853  * the number of bits set (i.e. `6` or larger).
854  */
855  template <std::size_t NBits, typename T>
856  static constexpr std::pair<std::array<std::size_t, NBits>, std::size_t>
857  setBitIndices(T value) noexcept;
858 
859 
860  // --- BEGIN -- Trees and their data ---------------------------------------
861 
862  /// Data structure for basic event information in simple ROOT trees.
864  unsigned int run; ///< Run number.
865  unsigned int subrun; ///< Subrun number.
866  unsigned int event; ///< Event number.
867  SplitTimestamp_t timestamp; ///< Event timestamp (seconds from the epoch).
868  }; // TreeData_EventID_t
869 
870  /// Structure collecting all data for a fragment ROOT tree.
871  struct TreeFragment_t {
872  struct Data_t: public TreeData_EventID_t {
873 
874  int fragmentID = 0; ///< ID of the fragment of this entry.
875 
876  unsigned int fragCount = 0U; ///< Counter of the fragment in the board.
877 
878  ///< Trigger time tag from the fragment.
879  unsigned long int TriggerTimeTag = 0;
880 
881  SplitTimestamp_t trigger; ///< Global trigger time.
882 
883  SplitTimestamp_t fragTime; ///< PMT fragment time stamp.
884 
885  long int relBeamGate; ///< Beam gate start relative to trigger [ns].
886 
887  /// Time assigned to the waveforms.
888  double waveformTime = std::numeric_limits<double>::lowest();
889 
890  unsigned int waveformSize = 0U; ///< Ticks in the waveforms.
891 
892  unsigned int triggerBits = 0x0; ///< Trigger bits, from `raw::Trigger`.
893 
894  unsigned int gateCount = 0U; ///< The number of gate from run start.
895 
896  ///< Whether waveforms cover nominal trigger time.
897  bool onGlobalTrigger = false;
898 
899 
900  }; // Data_t
901 
903  TTree* tree = nullptr;
904  }; // TreeFragment_t
905 
906 
907  std::unique_ptr<TreeData_EventID_t> fEventInfo; ///< Event ID for trees.
908 
909  ///< Tree with fragment information.
910  std::unique_ptr<TreeFragment_t> fTreeFragment;
911 
912  // --- END ---- Trees and their data -----------------------------------------
913 
914 
915  // --- BEGIN -- Timestamps ---------------------------------------------------
916 
917  /// Retrieves the global trigger time stamp from the event.
918  TriggerInfo_t fetchTriggerTimestamp(art::Event const& event) const;
919 
920  /// Returns the timestamp for the waveforms in the specified fragment.
921  electronics_time fragmentWaveformTimestamp(
922  FragmentInfo_t const& fragInfo,
923  NeededBoardInfo_t const& boardInfo,
924  SplitTimestamp_t triggerTime
925  ) const;
926 
927  /**
928  * @brief Returns the timestamp for the waveforms in the specified fragment.
929  *
930  * This method assumes that all PMT readout board triggers happened at the
931  * same time as the global trigger.
932  * Unsurprisingly, this does not work well with multiple PMT windows.
933  */
934  electronics_time fragmentWaveformTimestampOnTrigger(
935  FragmentInfo_t const& fragInfo,
936  NeededBoardInfo_t const& boardInfo,
937  SplitTimestamp_t triggerTime
938  ) const;
939 
940  /**
941  * @brief Returns the timestamp for the waveforms in the specified fragment.
942  *
943  * This method assumes that the Trigger Time Tag counters are synchronised
944  * with the global trigger and their value is reset on each new second of
945  * the global trigger timescale (International Atomic Time).
946  *
947  * This assumptions enables timestamping of waveforms from the same readout
948  * boards at different times ("multi-window PMT readout").
949  *
950  * See `TTTresetEverySecond` configuration option.
951  */
952  electronics_time fragmentWaveformTimestampFromTTT(
953  FragmentInfo_t const& fragInfo,
954  NeededBoardInfo_t const& boardInfo,
955  SplitTimestamp_t triggerTime
956  ) const;
957 
958 
959  /// Returns the timestamp difference `a - b`.
960  static long long int timestampDiff(std::uint64_t a, std::uint64_t b)
961  { return static_cast<long long int>(a) - static_cast<long long int>(b); }
962 
963 
964  /// Returns the fragment ID to be used with databases.
965  static constexpr std::size_t effectivePMTboardFragmentID
966  (artdaq::Fragment::fragment_id_t id)
967  { return id & 0x0fff; }
968 
969  // --- END ---- Timestamps ---------------------------------------------------
970 
971 
972  // --- BEGIN -- Input data management ----------------------------------------
973 
974  /// Tracks _art_ data products and removes their cached data on demand.
975  // As often happens, caches need to be mutable. Also, not thread-safe.
977 
978  /// Reads the fragments to be processed.
979  artdaq::Fragments const& readInputFragments(art::Event const& event) const;
980 
981  /// Throws an exception if `artdaqFragment` is not of type `CAEN1730`.
982  void checkFragmentType(artdaq::Fragment const& artdaqFragment) const;
983 
984  /// Converts a fragment into a fragment collection
985  /// (dispatcher based on fragment type).
986  artdaq::FragmentPtrs makeFragmentCollection
987  (artdaq::Fragment const& sourceFragment) const;
988 
989  /// Converts a plain fragment into a fragment collection.
990  artdaq::FragmentPtrs makeFragmentCollectionFromFragment
991  (artdaq::Fragment const& sourceFragment) const;
992 
993  /// Converts a container fragment into a fragment collection.
994  artdaq::FragmentPtrs makeFragmentCollectionFromContainerFragment
995  (artdaq::Fragment const& sourceFragment) const;
996 
997  /// Extracts waveforms from the specified fragments from a board.
998  std::vector<ProtoWaveform_t> processBoardFragments(
999  artdaq::FragmentPtrs const& artdaqFragment,
1000  TriggerInfo_t const& triggerInfo
1001  );
1002 
1003  // --- END ---- Input data management ----------------------------------------
1004 
1005 
1006  // --- BEGIN -- Output waveforms ---------------------------------------------
1007 
1008  /// Merges "in place" `waveforms` based on timestamp.
1009  /// @return the number of original waveforms merged into another
1010  unsigned int mergeWaveforms(std::vector<ProtoWaveform_t>& waveforms) const;
1011 
1012  /// Sorts in place the specified waveforms in channel order, then in time.
1013  void sortWaveforms(std::vector<ProtoWaveform_t>& waveforms) const;
1014 
1015  /// Returns pointers to all waveforms including the nominal trigger time.
1016  std::vector<ProtoWaveform_t const*> findWaveformsWithNominalTrigger
1017  (std::vector<ProtoWaveform_t> const& waveforms) const;
1018 
1019  /// Returns whether `waveform` includes the tick of the nominal trigger time.
1020  bool containsGlobalTrigger(raw::OpDetWaveform const& waveform) const;
1021 
1022  /// Returns whether nominal trigger time is within `nTicks` from `time`.
1023  bool containsGlobalTrigger(electronics_time time, std::size_t nTicks) const;
1024 
1025  /// Returns a waveform merging of the `indices` ones from `allWaveforms`.
1026  /// The merged waveforms are emptied of their content.
1027  ProtoWaveform_t mergeWaveformGroup(
1028  std::vector<ProtoWaveform_t>& allWaveforms,
1029  std::vector<std::size_t> const& indices
1030  ) const;
1031 
1033  { return electronics_time{ wf.TimeStamp() }; }
1034 
1036  { return waveformStartTime(wf.waveform); }
1037 
1039  { return waveformStartTime(wf) + fOpticalTick * wf.size(); }
1040 
1042  { return waveformEndTime(wf.waveform); }
1043 
1044  // --- END ---- Output waveforms ---------------------------------------------
1045 
1046  using BoardID_t = short int; ///< Type used internally to represent board ID.
1047 
1048 
1049  /**
1050  * @brief Create waveforms and fills trees for the specified artDAQ fragment.
1051  * @param artdaqFragment the fragment to process
1052  * @param boardInfo board information needed, from configuration/setup
1053  * @param triggerTime absolute time of the trigger
1054  * @return collection of PMT waveforms from the fragment
1055  *
1056  * This method fills the information for the PMT fragment tree
1057  * (`fillPMTfragmentTree()`) and creates PMT waveforms from the fragment data
1058  * (`createFragmentWaveforms()`).
1059  */
1060  std::vector<ProtoWaveform_t> processFragment(
1061  artdaq::Fragment const& artdaqFragment,
1062  NeededBoardInfo_t const& boardInfo,
1063  TriggerInfo_t const& triggerInfo
1064  );
1065 
1066 
1067  /**
1068  * @brief Creates `raw::OpDetWaveform` objects from the fragment data.
1069  * @param fragInfo information extracted from the fragment
1070  * @param channelSetup settings channel by channel
1071  * @param timeStamp timestamp of the waveforms in the fragment
1072  * @return collection of newly created `raw::OpDetWaveform`
1073  *
1074  * All fragment information needed is enclosed in `fragInfo`
1075  * (`extractFragmentInfo()`). The timestamp can be obtained with a call to
1076  * `fragmentWaveformTimestamp()`.
1077  */
1078  std::vector<ProtoWaveform_t> createFragmentWaveforms(
1079  FragmentInfo_t const& fragInfo,
1080  AllChannelSetup_t const& channelSetup,
1081  electronics_time const timeStamp
1082  ) const;
1083 
1084  /// Extracts useful information from fragment data.
1085  FragmentInfo_t extractFragmentInfo
1086  (artdaq::Fragment const& artdaqFragment) const;
1087 
1088  /// Extracts the fragment ID (i.e. board ID) from the specified `fragment`.
1089  static BoardID_t extractFragmentBoardID(artdaq::Fragment const& fragment);
1090 
1091  /// Returns the board information for this fragment.
1092  NeededBoardInfo_t neededBoardInfo
1093  (artdaq::Fragment::fragment_id_t fragment_id) const;
1094 
1095  /// Returns all the instance names we will produce.
1096  std::set<std::string> getAllInstanceNames() const;
1097 
1098  /// Throws an exception if the configuration of boards shows errors.
1099  void checkBoardSetup
1100  (std::vector<daq::details::BoardSetup_t> const& allBoardSetup) const;
1101 
1102 
1103  // --- BEGIN -- Tree-related methods -----------------------------------------
1104 
1105  /// Declares the use of event information.
1106  void usesEventInfo();
1107 
1108  /// Initializes all requested data trees.
1109  void initTrees(std::vector<std::string> const& treeNames);
1110 
1111  /// Initializes the event ID part of a tree.
1112  void initEventIDtree(TTree& tree, TreeData_EventID_t& data);
1113 
1114  /// Initializes the fragment data tree (`fTreeFragment`).
1115  void initFragmentsTree();
1116 
1117  /// Fills the base information of a tree data entry from an _art_ event.
1118  void fillTreeEventID
1119  (art::Event const& event, TreeData_EventID_t& treeData) const;
1120 
1121  /// Assigns the cached event information to the specified tree data.
1122  void assignEventInfo(TreeData_EventID_t& treeData) const;
1123 
1124  /// Fills the PMT fragment tree with the specified information
1125  /// (additional information needs to have been set already).
1126  void fillPMTfragmentTree(
1127  FragmentInfo_t const& fragInfo,
1128  TriggerInfo_t const& triggerInfo,
1129  electronics_time waveformTimestamp
1130  );
1131 
1132 
1133  /// Returns the name of the specified tree.
1134  static std::string const& treeName(DataTrees treeID);
1135 
1136  /// Static initialization.
1137  static TreeNameList_t initTreeNames();
1138 
1139  // --- END ---- Tree-related methods -----------------------------------------
1140 
1141  friend struct dumpChannel;
1142  friend std::ostream& operator<< (std::ostream&, ProtoWaveform_t const&);
1143 
1144 }; // icarus::DaqDecoderICARUSPMT
1145 
1146 
1147 //------------------------------------------------------------------------------
1148 // --- implementation
1149 //------------------------------------------------------------------------------
1150 namespace {
1151 
1152  /// Moves the contend of `src` into the end of `dest`.
1153  template <typename T>
1154  std::vector<T>& appendTo(std::vector<T>& dest, std::vector<T>&& src) {
1155  if (dest.empty()) dest = std::move(src);
1156  else {
1157  dest.reserve(dest.size() + src.size());
1158  std::move(src.begin(), src.end(), std::back_inserter(dest));
1159  }
1160  src.clear();
1161  return dest;
1162  } // appendTo()
1163 
1164  /// Returns whether an element `value` is found in `coll`.
1165  template <typename Coll, typename T>
1166  bool contains(Coll const& coll, T const& value) {
1167  auto const cend = end(coll);
1168  return std::find(begin(coll), cend, value) != cend;
1169  } // contains()
1170 
1171 } // local namespace
1172 
1173 
1174 //------------------------------------------------------------------------------
1175 namespace icarus {
1176 
1177  /// Special function `fhicl::TableAs` uses to convert BoardSetupConfig.
1180  {
1181 
1183  config.Name() // name
1184  , config.FragmentID()
1185  .value_or(daq::details::BoardSetup_t::NoFragmentID) // fragmentID
1186  , config.TriggerDelay() // triggerDelay
1187  , config.TTTresetDelay() // TTTresetDelay
1188  };
1189 
1190  // set the special configuration for the board channels that have one;
1191  // the others will stay with the default one (which implies that a channel
1192  // is saved in the standard data product, and only if it has a channel ID
1193  // entry from the database)
1194  for (
1196  : config.SpecialChannels().value_or(
1197  std::vector<DaqDecoderICARUSPMT::BoardSetupConfig::ChannelSetupConfig>{}
1198  )
1199  ) {
1200  try {
1201  bs.channelSettings.at(chConfig.ChannelIndex()) = {
1202  chConfig.Channel() // channelID
1203  , chConfig.Skip() // forcedSkip
1204  , chConfig.OnlyOnGlobalTrigger() // onGlobalOnly
1205  , chConfig.MinSpan() // minSpan
1206  , chConfig.InstanceName() // category
1207  };
1208  }
1209  catch (std::out_of_range const&) {
1210  throw art::Exception(art::errors::Configuration)
1211  << "Configuration requested for invalid channel index "
1212  << chConfig.ChannelIndex() << " of the board '"
1213  << config.Name() << "'\n";
1214  }
1215  } // for
1216 
1217  return bs;
1218  } // convert(BoardSetupConfig)
1219 
1220 
1221  struct dumpChannel {
1224  };
1225 
1226  std::ostream& operator<< (std::ostream& out, dumpChannel const& d) {
1227  return out << (d.wf.channelSetup->category.empty()? std::dec: std::hex)
1228  << d.wf.waveform.ChannelNumber() << std::dec;
1229  }
1230 
1231 } // namespace icarus
1232 
1233 
1234 //------------------------------------------------------------------------------
1235 // --- icarus::DaqDecoderICARUSPMT::SplitTimestamp_t
1236 //------------------------------------------------------------------------------
1238  (int sec, unsigned int ns)
1239  : time { static_cast<long long int>(sec) * 1'000'000'000LL + ns }
1240  , split { sec, ns }
1241 {}
1242 
1243 
1244 //------------------------------------------------------------------------------
1245 constexpr icarus::DaqDecoderICARUSPMT::SplitTimestamp_t::SplitTimestamp_t
1246  (long long int triggerTime)
1247  : time { triggerTime }
1248  , split {
1249  static_cast<int>(time / 1'000'000'000), // seconds
1250  static_cast<unsigned int>(time % 1'000'000'000) // nanoseconds
1251  }
1252 {}
1253 
1254 
1255 //------------------------------------------------------------------------------
1256 namespace icarus {
1257 
1258  std::ostream& operator<<
1259  (std::ostream& out, DaqDecoderICARUSPMT::SplitTimestamp_t const& time)
1260  {
1261  out << time.split.seconds << '.'
1262  << std::setfill('0') << std::setw(9) << time.split.nanoseconds;
1263  return out;
1264  } // operator<< (std::ostream&, PMTDecoder::SplitTimestamp_t)
1265 
1266 } // namespace icarus
1267 
1268 
1269 //------------------------------------------------------------------------------
1270 // --- icarus::DaqDecoderICARUSPMT
1271 //------------------------------------------------------------------------------
1272 // --- template implementation
1273 //------------------------------------------------------------------------------
1274 template <std::size_t NBits, typename T>
1275 constexpr std::pair<std::array<std::size_t, NBits>, std::size_t>
1276 icarus::DaqDecoderICARUSPMT::setBitIndices(T value) noexcept {
1277 
1278  std::pair<std::array<std::size_t, NBits>, std::size_t> res;
1279  auto& [ indices, nSetBits ] = res;
1280  for (std::size_t& index: indices) {
1281  index = (value & 1)? nSetBits++: NBits;
1282  value >>= 1;
1283  } // for
1284  return res;
1285 
1286 } // icarus::DaqDecoderICARUSPMT::setBitIndices()
1287 
1288 
1289 //------------------------------------------------------------------------------
1290 // --- static definitions
1291 //------------------------------------------------------------------------------
1292 icarus::DaqDecoderICARUSPMT::AllChannelSetup_t const
1293 icarus::DaqDecoderICARUSPMT::NeededBoardInfo_t::DefaultChannelSetup;
1294 
1295 //------------------------------------------------------------------------------
1296 icarus::DaqDecoderICARUSPMT::TreeNameList_t const
1297 icarus::DaqDecoderICARUSPMT::TreeNames
1298  = icarus::DaqDecoderICARUSPMT::initTreeNames();
1299 
1300 auto icarus::DaqDecoderICARUSPMT::initTreeNames() -> TreeNameList_t {
1301  TreeNameList_t names;
1302  names[static_cast<std::size_t>(DataTrees::Fragments)] = "PMTfragments";
1303  return names;
1304 } // icarus::DaqDecoderICARUSPMT::initTreeNames()
1305 
1306 
1307 //------------------------------------------------------------------------------
1308 std::string const& icarus::DaqDecoderICARUSPMT::treeName(DataTrees treeID)
1309  { return TreeNames[static_cast<std::size_t>(treeID)]; }
1310 
1311 
1312 //------------------------------------------------------------------------------
1313 std::string icarus::DaqDecoderICARUSPMT::listTreeNames
1314  (std::string const& sep /* = " " */)
1315 {
1316  std::string l;
1317  for (std::string const& name: TreeNames) {
1318  if (!l.empty()) l += sep;
1319  l += '\'';
1320  l += name;
1321  l += '\'';
1322  } // for
1323  return l;
1324 } // icarus::DaqDecoderICARUSPMT::listTreeNames()
1325 
1326 
1327 //------------------------------------------------------------------------------
1328 // --- implementation
1329 //------------------------------------------------------------------------------
1331  : art::EDProducer(params)
1332  , fInputTags{ params().FragmentsLabels() }
1333  , fSurviveExceptions{ params().SurviveExceptions() }
1334  , fDiagnosticOutput{ params().DiagnosticOutput() }
1335  , fPacketDump{ params().PacketDump() }
1336  , fRequireKnownBoards{ params().RequireKnownBoards() }
1337  , fRequireBoardConfig{ params().RequireBoardConfig() }
1338  , fPMTconfigTag{ params().PMTconfigTag() }
1339  , fTriggerTag{ params().TriggerTag() }
1340  , fTTTresetEverySecond
1341  { params().TTTresetEverySecond().value_or(fTriggerTag.has_value()) }
1342  , fBoardSetup{ params().BoardSetup() }
1343  , fSkipWaveforms{ params().SkipWaveforms() }
1344  , fDropRawDataAfterUse{ params().DropRawDataAfterUse() }
1345  , fLogCategory{ params().LogCategory() }
1346  , fDetTimings
1347  { art::ServiceHandle<detinfo::DetectorClocksService const>()->DataForJob() }
1348  , fChannelMap{ *(art::ServiceHandle<icarusDB::IICARUSChannelMap const>{}) }
1349  , fOpticalTick{ fDetTimings.OpticalClockPeriod() }
1350  , fNominalTriggerTime{ fDetTimings.TriggerTime() }
1351 {
1352  //
1353  // configuration check
1354  //
1355  checkBoardSetup(fBoardSetup); // throws on error
1356 
1357  //
1358  // consumed data products declaration
1359  //
1360  for (art::InputTag const& inputTag: fInputTags)
1361  consumes<artdaq::Fragments>(inputTag);
1362  if (fPMTconfigTag) consumes<sbn::PMTconfiguration>(*fPMTconfigTag);
1363  if (fTriggerTag) {
1364  consumes<std::vector<sbn::ExtraTriggerInfo>>(*fTriggerTag);
1365  if (contains(params().DataTrees(), treeName(DataTrees::Fragments)))
1366  consumes<std::vector<raw::Trigger>>(*fTriggerTag);
1367  } // if trigger
1368 
1369  //
1370  // produced data products declaration
1371  //
1372  if (!fSkipWaveforms) {
1373  for (std::string const& instanceName: getAllInstanceNames())
1374  produces<std::vector<raw::OpDetWaveform>>(instanceName);
1375  }
1376 
1377  //
1378  // additional initialization
1379  //
1380  initTrees(params().DataTrees());
1381 
1382 
1383  //
1384  // configuration dump
1385  //
1386  mf::LogInfo log(fLogCategory);
1387  log << "Configuration:"
1388  << "\n * data from one of " << fInputTags.size() << " data products:";
1389  for (art::InputTag const& inputTag: fInputTags)
1390  log << " '" << inputTag.encode() << "'";
1391  log
1392  << "\n * boards with setup: " << fBoardSetup.size();
1393  if (fPMTconfigTag)
1394  log << "\n * PMT configuration from '" << fPMTconfigTag->encode() << "'";
1395  else
1396  log << "\n * PMT configuration not used (and some corrections will be skipped)";
1397  if (fTriggerTag)
1398  log << "\n * trigger information from: '" << fTriggerTag->encode() << '\'';
1399  else
1400  log << "\n * trigger time from event timestamp [fallback]";
1401  if (fRequireKnownBoards) {
1402  log << "\n * all readout boards in input must be known (from `"
1403  << params().BoardSetup.name() << "` or `"
1404  << params().PMTconfigTag.name() << "`)"
1405  ;
1406  }
1407  else {
1408  log << "\n * readout boards with no information (from neither `"
1409  << params().BoardSetup.name() << "` or `"
1410  << params().PMTconfigTag.name()
1411  << "`) are processed at the best we can (skipping corrections)"
1412  ;
1413  }
1414  if (fRequireBoardConfig) {
1415  log << "\n * all readout boards in `"
1416  << params().BoardSetup.name()
1417  << "` must appear in the PMT configuration from `"
1418  << params().PMTconfigTag.name() << "`"
1419  ;
1420  }
1421  else {
1422  log << "\n * all readout boards in `"
1423  << params().BoardSetup.name()
1424  << "` may lack a matching PMT configuration from `"
1425  << params().PMTconfigTag.name() << "`"
1426  ;
1427  }
1428  if (fSkipWaveforms) {
1429  log << "\n * PMT WAVEFORMS WILL NOT BE DECODED AND STORED";
1430  }
1431 
1432  //
1433  // sanity checks
1434  //
1435 
1436  if (fChannelMap.nPMTfragmentIDs() == 0) {
1437  throw cet::exception("DaqDecoderICARUSPMT")
1438  << "Channel mapping database does not report any PMT fragment ID!\n";
1439  } // if no PMT in database
1440 
1441 
1442 } // icarus::DaqDecoderICARUSPMT::DaqDecoderICARUSPMT()
1443 
1444 
1445 //------------------------------------------------------------------------------
1446 void icarus::DaqDecoderICARUSPMT::beginRun(art::Run& run) {
1447 
1448  //sbn::PMTconfiguration const* PMTconfig = fPMTconfigTag
1449  // ? run.getPointerByLabel<sbn::PMTconfiguration>(*fPMTconfigTag): nullptr;
1450  sbn::PMTconfiguration const* PMTconfig = fPMTconfigTag
1451  ? run.getHandle<sbn::PMTconfiguration>(*fPMTconfigTag).product(): nullptr;
1452 
1453  UpdatePMTConfiguration(PMTconfig);
1454 
1455 } // icarus::DaqDecoderICARUSPMT::beginRun()
1456 
1457 
1458 //------------------------------------------------------------------------------
1459 void icarus::DaqDecoderICARUSPMT::produce(art::Event& event) {
1460 
1461  // ---------------------------------------------------------------------------
1462  // preparation
1463  //
1464 
1465  fDataCacheRemover.useEvent(event);
1466 
1467  //
1468  // global trigger
1469  //
1470  TriggerInfo_t const triggerInfo = fetchTriggerTimestamp(event);
1471  {
1472  mf::LogDebug log { fLogCategory };
1473  if (fTriggerTag)
1474  log << "Trigger time ('" << fTriggerTag->encode() << "'): ";
1475  else
1476  log << "Trigger from event timestamp: ";
1477  log << triggerInfo.time << " s, bits: "
1478  << icarus::ns::util::bin(triggerInfo.bits.bits);
1479  if (triggerInfo.bits) {
1480  log << " {";
1481  for (std::string const& name: names(triggerInfo.bits)) log << ' ' << name;
1482  log << " }";
1483  } // if
1484  if (fTriggerTag) log << ", spill count: " << triggerInfo.gateCount;
1485  } // local block
1486 
1487  //
1488  // event ID
1489  //
1490 
1491  // if needed, fill the record with the basic information of the event
1492  if (fEventInfo) fillTreeEventID(event, *fEventInfo);
1493 
1494  //
1495  // output data product initialization
1496  //
1497  std::vector<ProtoWaveform_t> protoWaveforms; // empty if `fSkipWaveforms`
1498 
1499 
1500  // ---------------------------------------------------------------------------
1501  // pre-processing
1502  //
1503 
1504 
1505  // ---------------------------------------------------------------------------
1506  // processing
1507  //
1508 
1509  std::unordered_map<BoardID_t, unsigned int> boardCounts;
1510  bool duplicateBoards = false;
1511  try { // catch-all
1512  auto const& fragments = readInputFragments(event);
1513 
1514  for (artdaq::Fragment const& fragment: fragments) {
1515 
1516  artdaq::FragmentPtrs const& fragmentCollection
1517  = makeFragmentCollection(fragment);
1518 
1519  if (empty(fragmentCollection)) {
1520  mf::LogWarning("DaqDecoderICARUSPMT")
1521  << "Found a data fragment (ID=" << extractFragmentBoardID(fragment)
1522  << ") containing no data.";
1523  continue;
1524  } // if no data
1525 
1526  BoardID_t const boardID
1527  = extractFragmentBoardID(*(fragmentCollection.front()));
1528  if (++boardCounts[boardID] > 1U) duplicateBoards = true;
1529 
1530  appendTo(
1531  protoWaveforms,
1532  processBoardFragments(fragmentCollection, triggerInfo)
1533  );
1534 
1535  } // for all input fragments
1536 
1537  }
1538  catch (cet::exception const& e) {
1539  if (!fSurviveExceptions) throw;
1540  mf::LogError("DaqDecoderICARUSPMT")
1541  << "Error while attempting to decode PMT data:\n" << e.what() << '\n';
1542  protoWaveforms.clear();
1543  ++fNFailures;
1544  }
1545  catch (...) {
1546  if (!fSurviveExceptions) throw;
1547  mf::LogError("DaqDecoderICARUSPMT")
1548  << "Error while attempting to decode PMT data.\n";
1549  protoWaveforms.clear();
1550  ++fNFailures;
1551  }
1552 
1553  if (duplicateBoards) {
1554  mf::LogWarning log { "DaqDecoderICARUSPMT" };
1555  log << "found multiple data product entries for the same board:";
1556  for (auto [ boardID, count ]: boardCounts) {
1557  if (count < 2U) continue;
1558  log << " " << std::hex << boardID << std::dec << " (x" << count << ");";
1559  } // for
1560  } // if duplicate board fragments
1561 
1562 
1563  // we are done with the input: drop the caches
1564  // (if we were asked not to, no data is registered)
1565  fDataCacheRemover.removeCachedProducts();
1566 
1567  //
1568  // post-processing
1569  //
1570  sortWaveforms(protoWaveforms);
1571 
1572  if (!fSkipWaveforms) {
1573  std::vector<ProtoWaveform_t const*> const waveformsWithTrigger
1574  = findWaveformsWithNominalTrigger(protoWaveforms);
1575  mf::LogTrace(fLogCategory) << waveformsWithTrigger.size() << "/"
1576  << protoWaveforms.size() << " decoded waveforms include trigger time ("
1577  << fNominalTriggerTime << ").";
1578  } // if !fSkipWaveforms
1579 
1580  // ---------------------------------------------------------------------------
1581  // output
1582  //
1583 
1584  if (!fSkipWaveforms) {
1585  // split the waveforms by destination
1586  std::map<std::string, std::vector<raw::OpDetWaveform>> waveformProducts;
1587  for (std::string const& instanceName: getAllInstanceNames())
1588  waveformProducts.emplace(instanceName, std::vector<raw::OpDetWaveform>{});
1589  for (ProtoWaveform_t& waveform: protoWaveforms) {
1590 
1591  // on-global and span requirements override even `mustSave()` requirement;
1592  // if this is not good, user should not set `mustSave()`!
1593  bool const keep =
1594  (waveform.onGlobal || !waveform.channelSetup->onGlobalOnly)
1595  && (waveform.span() >= waveform.channelSetup->minSpan)
1596  ;
1597 
1598  if (!keep) continue;
1599  waveformProducts.at(waveform.channelSetup->category).push_back
1600  (std::move(waveform.waveform));
1601  } // for
1602 
1603  // put all the categories
1604  for (auto&& [ category, waveforms ]: waveformProducts) {
1605  mf::LogTrace(fLogCategory)
1606  << waveforms.size() << " PMT waveforms saved for "
1607  << (category.empty()? "standard": category) << " instance.";
1608  event.put(
1609  std::make_unique<std::vector<raw::OpDetWaveform>>(std::move(waveforms)),
1610  category // the instance name is the category the waveforms belong to
1611  );
1612  }
1613  } // if !fSkipWaveforms
1614 
1615 } // icarus::DaqDecoderICARUSPMT::produce()
1616 
1617 
1618 //------------------------------------------------------------------------------
1620 
1621  if (fNFailures > 0U) {
1622  mf::LogError(fLogCategory) << "Encountered errors on " << fNFailures
1623  << " events. Errors were ignored.";
1624  }
1625 
1626 } // icarus::DaqDecoderICARUSPMT::endJob()
1627 
1628 
1629 //------------------------------------------------------------------------------
1630 artdaq::Fragments const& icarus::DaqDecoderICARUSPMT::readInputFragments
1631  (art::Event const& event) const
1632 {
1633  art::Handle<artdaq::Fragments> handle;
1634  art::InputTag selectedInputTag; // empty
1635  for (art::InputTag const& inputTag: fInputTags) {
1636 
1637  mf::LogTrace(fLogCategory)
1638  << "DaqDecoderICARUSPMT trying data product: '" << inputTag.encode()
1639  << "'";
1640 
1641  auto const& thisHandle = event.getHandle<artdaq::Fragments>(inputTag);
1642  if (!thisHandle.isValid() || thisHandle->empty()) continue;
1643 
1644  mf::LogTrace(fLogCategory)
1645  << " => data product: '" << inputTag.encode() << "' is present and has "
1646  << thisHandle->size() << " entries";
1647 
1648  if (!selectedInputTag.empty()) {
1649  throw cet::exception("DaqDecoderICARUSPMT")
1650  << "Found multiple suitable input candidates: '"
1651  << inputTag.encode() << "' and '" << selectedInputTag << "'\n";
1652  }
1653  selectedInputTag = inputTag;
1654  handle = thisHandle;
1655  } // for
1656 
1657  if (!handle.isValid()) {
1658  cet::exception e { "DaqDecoderICARUSPMT" };
1659  e << "No suitable input data product found among:";
1660  for (art::InputTag const& inputTag: fInputTags)
1661  e << " '" << inputTag.encode() << "'";
1662  throw e << "\n";
1663  }
1664 
1665  if (fDropRawDataAfterUse) fDataCacheRemover.registerHandle(handle);
1666 
1667  return *handle;
1668 } // icarus::DaqDecoderICARUSPMT::readInputFragments()
1669 
1670 
1671 //------------------------------------------------------------------------------
1673  (std::vector<ProtoWaveform_t> const& waveforms) const
1674  -> std::vector<ProtoWaveform_t const*>
1675 {
1676  std::vector<ProtoWaveform_t const*> matchedWaveforms;
1677  for (ProtoWaveform_t const& waveform: waveforms) {
1678  if (waveform.onGlobal) matchedWaveforms.push_back(&waveform);
1679  } // for
1680  return matchedWaveforms;
1681 } // icarus::DaqDecoderICARUSPMT::findWaveformsWithNominalTrigger()
1682 
1683 
1684 //------------------------------------------------------------------------------
1686  (electronics_time time, std::size_t nTicks) const
1687 {
1688  return (fNominalTriggerTime >= time)
1689  && (fNominalTriggerTime < time + nTicks * fOpticalTick);
1690 } // icarus::DaqDecoderICARUSPMT::containsGlobalTrigger(time, ticks)
1691 
1692 
1693 //------------------------------------------------------------------------------
1695  (raw::OpDetWaveform const& waveform) const
1696 {
1697  return containsGlobalTrigger
1698  (electronics_time{ waveform.TimeStamp() }, waveform.size());
1699 } // icarus::DaqDecoderICARUSPMT::containsGlobalTrigger(waveform)
1700 
1701 
1702 //------------------------------------------------------------------------------
1704  (sbn::PMTconfiguration const* PMTconfig)
1705 {
1706  fBoardInfoLookup.emplace(matchBoardConfigurationAndSetup(PMTconfig));
1707 
1708  mf::LogDebug(fLogCategory)
1709  << "Board information as cached:\n" << *fBoardInfoLookup;
1710 
1711  return true;
1712 } // icarus::DaqDecoderICARUSPMT::UpdatePMTConfiguration()
1713 
1714 
1716  (sbn::PMTconfiguration const* PMTconfig) const
1718 {
1719  /*
1720  * We need to support the case where no PMT configuration is known
1721  * (that is the standard situation in the online monitor).
1722  * The "strategy" is that in such cases we give up the correct time stamp
1723  * decoding; if the setup information contains a fragment ID, it may be
1724  * possible to do a little better, that is to use the setup information
1725  * (this is not possible without knowing the fragment ID that each bit of
1726  * setup information pertains).
1727  *
1728  * So the cases for a board are:
1729  * * setup information is not present: encountering such a board will cause
1730  * an exception to be thrown (implemented elsewhere)
1731  * * PMT configuration and setup present: full configuration
1732  * * exception thrown if setup fragment ID is present and inconsistent
1733  * * PMT configuration not present: a general warning is printed;
1734  * * boards with setup fragment ID information: add setup information
1735  * to the "database" for the board: it will be used for partial
1736  * timestamp reconstruction
1737  * * boards without setup fragment ID information: board will not be
1738  * added into the database; no specific correction will be performed;
1739  * a warning is printed for each board
1740  *
1741  */
1742 
1743  // dictionary of board configurations (if any)
1744  std::vector<std::pair<std::string, sbn::V1730Configuration const*>>
1745  configByName;
1746  if (PMTconfig) {
1747  if (!PMTconfig->boards.empty())
1748  configByName.reserve(PMTconfig->boards.size());
1749  for (sbn::V1730Configuration const& boardConfig: PMTconfig->boards)
1750  configByName.emplace_back(boardConfig.boardName, &boardConfig);
1751  std::sort(configByName.begin(), configByName.end()); // sorted by board name
1752  } // if we have configuration
1753 
1754 
1755  auto findPMTconfig = [this, &configByName]
1756  (std::string const& name) -> sbn::V1730Configuration const*
1757  {
1758  if (!hasPMTconfiguration()) return nullptr;
1759  auto const* ppBoardConfig
1760  = daq::details::binarySearch(configByName, name);
1761  if (!ppBoardConfig) {
1762  if (!fRequireBoardConfig) return nullptr;
1763  throw cet::exception("DaqDecoderICARUSPMT")
1764  << "No DAQ configuration found for PMT readout board '"
1765  << name << "'\n"
1766  "If this is expected, you may skip this check by setting "
1767  "DaqDecoderICARUSPMT module configuration `RequireBoardConfig`"
1768  " to `false`.\n"
1769  ;
1770  }
1771  return ppBoardConfig->second;
1772  }; // findPMTconfig()
1773 
1774  // the filling is driven by boards configured in the module
1775  // (which is the reason a setup entry is mandatory)
1776  daq::details::BoardInfoLookup::Database_t boardInfoByFragment;
1777 
1778  for (daq::details::BoardSetup_t const& boardSetup: fBoardSetup) {
1779 
1780  std::string const& boardName = boardSetup.name;
1781 
1782  sbn::V1730Configuration const* pBoardConfig = findPMTconfig(boardName);
1783 
1784  if (pBoardConfig) {
1785  // fragment ID from configuration and setup must match if both present
1786  if (boardSetup.hasFragmentID()
1787  && (boardSetup.fragmentID != pBoardConfig->fragmentID)
1788  ) {
1789  throw cet::exception("DaqDecoderICARUSPMT")
1790  << "Board '" << boardName << "' has fragment ID "
1791  << std::hex << pBoardConfig->fragmentID << std::dec
1792  << " but it is set up as "
1793  << std::hex << boardSetup.fragmentID << std::dec
1794  << "!\n";
1795  } // if fragment ID in setup
1796  }
1797  else {
1798  if (boardSetup.hasFragmentID()) {
1799  mf::LogPrint(fLogCategory)
1800  << "Board '" << boardName
1801  << "' has no configuration information;"
1802  " some time stamp corrections will be skipped.";
1803  // to avoid this, make a PMT configuration available from input file
1804  }
1805  else {
1806  mf::LogPrint(fLogCategory)
1807  << "Board '" << boardName
1808  << "' can't be associated to a fragment ID;"
1809  " its time stamp corrections will be skipped.";
1810  // to avoid this, add a `BoardSetup.FragmentID` entry for it in the
1811  // configuration of this module, or make a PMT configuration available
1812  continue; // no entry for this board at all
1813  }
1814  }
1815 
1816  unsigned int const fragmentID
1817  = pBoardConfig? pBoardConfig->fragmentID: boardSetup.fragmentID;
1818  assert(fragmentID != daq::details::BoardSetup_t::NoFragmentID);
1819 
1820  nanoseconds const preTriggerTime
1821  = pBoardConfig
1822  ? (pBoardConfig->bufferLength * (1.0f - pBoardConfig->postTriggerFrac))
1823  * fOpticalTick
1824  : nanoseconds{ 0.0 }
1825  ;
1826 
1827  daq::details::BoardFacts_t boardFacts {
1828  preTriggerTime // ditto
1829  };
1830 
1831  boardInfoByFragment.push_back({
1832  fragmentID, // fragmentID
1833  &boardSetup, // setup
1834  pBoardConfig, // config
1835  std::move(boardFacts) // facts
1836  });
1837  } // for
1838 
1839  return daq::details::BoardInfoLookup{ std::move(boardInfoByFragment) };
1840 
1841 } // icarus::DaqDecoderICARUSPMT::matchBoardConfigurationAndSetup()
1842 
1843 
1844 //------------------------------------------------------------------------------
1847  unsigned int fragmentID
1848 ) const -> NeededBoardInfo_t {
1849 
1851 
1852  return NeededBoardInfo_t{
1853  // name
1854  ((boardInfo && boardInfo->config)
1855  ? boardInfo->config->boardName: ("<ID=" + std::to_string(fragmentID)))
1856  // bufferLength
1857  , ((boardInfo && boardInfo->config)
1858  ? boardInfo->config->bufferLength * fOpticalTick: nanoseconds{ 0.0 }
1859  )
1860  // preTriggerTime
1861  , (boardInfo? boardInfo->facts.preTriggerTime: nanoseconds{ 0.0 })
1862  // PMTtriggerDelay
1863  , ((boardInfo && boardInfo->setup)
1864  ? boardInfo->setup->triggerDelay: nanoseconds{ 0.0 })
1865  // TTTresetDelay
1866  , ((boardInfo && boardInfo->setup)
1867  ? boardInfo->setup->TTTresetDelay: nanoseconds{ 0.0 })
1868  , ((boardInfo && boardInfo->setup)?
1869  &(boardInfo->setup->channelSettings): nullptr)
1870  };
1871 
1872 } // icarus::DaqDecoderICARUSPMT::fetchNeededBoardInfo()
1873 
1874 
1875 //------------------------------------------------------------------------------
1877  (art::Event const& event) const -> TriggerInfo_t
1878 {
1879 
1880  if (!fTriggerTag)
1881  return { SplitTimestamp_t(event.time().value()), 0U };
1882 
1883 
1884  auto const& extraTrigger
1885  = event.getProduct<sbn::ExtraTriggerInfo>(*fTriggerTag);
1886  if (!extraTrigger.isValid()) {
1887  // this means there is some problem from trigger decoder;
1888  // while we might recover additional information from other data products,
1889  // we expect those to be as problematic.
1890  throw cet::exception("DaqDecoderICARUSPMT")
1891  << "Extra trigger information from '"
1892  << fTriggerTag->encode() << "' is marked as invalid!\n";
1893  }
1894 
1895  if (fDiagnosticOutput) {
1896  mf::LogVerbatim(fLogCategory)
1897  << "Extended trigger information:\n" << extraTrigger;
1898  }
1899 
1900  auto const& triggers
1901  = event.getProduct<std::vector<raw::Trigger>>(*fTriggerTag);
1902  if (triggers.size() != 1U) {
1903  // if this is hit, the decoder needs some development to correctly deal
1904  // with input with no trigger, or more than one
1905  throw cet::exception("DaqDecoderICARUSPMT")
1906  << "Found " << triggers.size() << " raw::Trigger from '"
1907  << fTriggerTag->encode() << "', can deal only with 1.\n";
1908  }
1909  raw::Trigger const& trigger = triggers.front();
1910 
1911  long long int const relBeamGate = timestampDiff
1912  (extraTrigger.beamGateTimestamp, extraTrigger.triggerTimestamp);
1913 
1914  unsigned int const gateCount = extraTrigger.gateID;
1915 
1916  return {
1917  SplitTimestamp_t
1918  { static_cast<long long int>(extraTrigger.triggerTimestamp) }
1919  , relBeamGate
1920  , {trigger.TriggerBits()}
1921  , gateCount
1922  };
1923 
1924 } // icarus::DaqDecoderICARUSPMT::fetchTriggerTimestamp()
1925 
1926 
1927 //------------------------------------------------------------------------------
1929  (artdaq::Fragment const& sourceFragment) const
1930 {
1931  switch (sourceFragment.type()) {
1932  case sbndaq::FragmentType::CAENV1730:
1933  return makeFragmentCollectionFromFragment(sourceFragment);
1934  case artdaq::Fragment::ContainerFragmentType:
1935  return makeFragmentCollectionFromContainerFragment(sourceFragment);
1936  default:
1937  throw cet::exception("DaqDecoderICARUSPMT")
1938  << "Unexpected PMT data product fragment type: "
1939  << static_cast<int>(sourceFragment.type()) << " ('"
1940  << sbndaq::fragmentTypeToString
1941  (static_cast<sbndaq::FragmentType>(sourceFragment.type()))
1942  << "')\n";
1943  } // switch
1944 } // icarus::DaqDecoderICARUSPMT::makeFragmentCollection()
1945 
1946 
1947 //------------------------------------------------------------------------------
1948 artdaq::FragmentPtrs
1950  (artdaq::Fragment const& sourceFragment) const
1951 {
1952  assert(sourceFragment.type() == sbndaq::FragmentType::CAENV1730);
1953  artdaq::FragmentPtrs fragColl;
1954  fragColl.push_back(std::make_unique<artdaq::Fragment>(sourceFragment));
1955  return fragColl;
1956 } // icarus::DaqDecoderICARUSPMT::makeFragmentCollectionFromFragment()
1957 
1958 
1959 //------------------------------------------------------------------------------
1960 artdaq::FragmentPtrs
1962  (artdaq::Fragment const& sourceFragment) const
1963 {
1964  assert(sourceFragment.type() == artdaq::Fragment::ContainerFragmentType);
1965  artdaq::ContainerFragment const containerFragment{ sourceFragment };
1966 
1967  if (containerFragment.block_count() == 0) return {};
1968 
1969  artdaq::FragmentPtrs fragColl;
1970  for (auto const iFrag: util::counter(containerFragment.block_count()))
1971  fragColl.push_back(containerFragment.at(iFrag));
1972 
1973  return fragColl;
1974 } // icarus::DaqDecoderICARUSPMT::makeFragmentCollectionFromContainerFragment()
1975 
1976 
1977 //------------------------------------------------------------------------------
1979  (artdaq::Fragment const& artdaqFragment) const
1980 {
1981  if (artdaqFragment.type() == sbndaq::FragmentType::CAENV1730) return;
1982 
1983  throw cet::exception("DaqDecoderICARUSPMT")
1984  << "Unexpected PMT fragment data type: '"
1985  << sbndaq::fragmentTypeToString
1986  (static_cast<sbndaq::FragmentType>(artdaqFragment.type()))
1987  << "'\n";
1988 
1989 } // icarus::DaqDecoderICARUSPMT::checkFragmentType
1990 
1991 
1992 //------------------------------------------------------------------------------
1994  artdaq::FragmentPtrs const& artdaqFragments,
1995  TriggerInfo_t const& triggerInfo
1996 ) -> std::vector<ProtoWaveform_t> {
1997 
1998  if (artdaqFragments.empty()) return {};
1999 
2000  artdaq::Fragment const& referenceFragment = *(artdaqFragments.front());
2001  assert(&referenceFragment);
2002 
2003  checkFragmentType(referenceFragment);
2004 
2005  NeededBoardInfo_t const boardInfo
2006  = neededBoardInfo(artdaqFragments.front()->fragmentID());
2007 
2008  mf::LogTrace(fLogCategory)
2009  << " - " << boardInfo.name << ": " << artdaqFragments.size()
2010  << " fragments";
2011 
2012  std::vector<ProtoWaveform_t> waveforms;
2013  for (artdaq::FragmentPtr const& fragment: artdaqFragments)
2014  appendTo(waveforms, processFragment(*fragment, boardInfo, triggerInfo));
2015 
2016  mergeWaveforms(waveforms);
2017 
2018  return { waveforms };
2019 
2020 } // icarus::DaqDecoderICARUSPMT::processBoardFragments()
2021 
2022 
2023 //------------------------------------------------------------------------------
2025  artdaq::Fragment const& artdaqFragment,
2026  NeededBoardInfo_t const& boardInfo,
2027  TriggerInfo_t const& triggerInfo
2028 ) -> std::vector<ProtoWaveform_t> {
2029 
2030  checkFragmentType(artdaqFragment);
2031 
2032  if (fPacketDump) {
2033  mf::LogVerbatim{ fLogCategory } << "PMT packet:"
2034  << "\n" << std::string(80, '-')
2035  << "\n" << sbndaq::dumpFragment(artdaqFragment)
2036  << "\n" << std::string(80, '-')
2037  ;
2038  } // if diagnostics
2039 
2040  FragmentInfo_t const fragInfo = extractFragmentInfo(artdaqFragment);
2041 
2042  auto const timeStamp
2043  = fragmentWaveformTimestamp(fragInfo, boardInfo, triggerInfo.time);
2044 
2045  if (fTreeFragment) fillPMTfragmentTree(fragInfo, triggerInfo, timeStamp);
2046 
2047  return ((timeStamp != NoTimestamp) && !fSkipWaveforms)
2048  ? createFragmentWaveforms(fragInfo, boardInfo.channelSetup(), timeStamp)
2049  : std::vector<ProtoWaveform_t>{}
2050  ;
2051 
2052 } // icarus::DaqDecoderICARUSPMT::processFragment()
2053 
2054 
2055 //------------------------------------------------------------------------------
2057  FragmentInfo_t const& fragInfo, AllChannelSetup_t const& channelSetup,
2058  electronics_time const timeStamp
2059 ) const -> std::vector<ProtoWaveform_t>
2060 {
2061 
2062  assert(timeStamp != NoTimestamp);
2063 
2064  std::vector<ProtoWaveform_t> protoWaveforms; // output collection
2065 
2066  std::optional<mf::LogVerbatim> diagOut;
2067  if (fDiagnosticOutput) diagOut.emplace(fLogCategory);
2068 
2069  icarusDB::DigitizerChannelChannelIDPairVec const& digitizerChannelVec
2070  = fChannelMap.getChannelIDPairVec
2071  (effectivePMTboardFragmentID(fragInfo.fragmentID))
2072  ;
2073 
2074  // allocate the vector outside the loop since we'll reuse it over and over
2075  std::vector<std::uint16_t> wvfm(fragInfo.nSamplesPerChannel);
2076 
2077  // all waveforms share the same timestamp,
2078  // so either all contain the global trigger, or they all do not
2079  bool const onGlobal = containsGlobalTrigger(timeStamp, wvfm.size());
2080 
2081  auto channelNumberToChannel
2082  = [&digitizerChannelVec](unsigned short int channelNumber) -> raw::Channel_t
2083  {
2084  for (auto const [ chNo, chID ]: digitizerChannelVec)
2085  if (chNo == channelNumber) return chID;
2087  };
2088 
2089  if (diagOut)
2090  (*diagOut) << " " << digitizerChannelVec.size() << " channels:";
2091 
2092  std::size_t iNextChunk = 0;
2093  for (unsigned short int const channelNumber: util::counter(16U)) {
2094 
2095  if ((fragInfo.enabledChannels & (1 << channelNumber)) == 0) {
2096  if (diagOut)
2097  (*diagOut) << " " << channelNumber << " [disabled];";
2098  continue;
2099  }
2100 
2101  std::size_t const iChunk = iNextChunk++;
2102 
2103  daq::details::BoardSetup_t::ChannelSetup_t const& thisChannelSetup
2104  = channelSetup.at(channelNumber); // this setup must be available
2105 
2106  if (thisChannelSetup.mustSkip()) {
2107  mf::LogTrace(fLogCategory)
2108  << "Channel number " << channelNumber << " of board 0x"
2109  << std::hex << fragInfo.fragmentID << std::dec
2110  << " is enabled but is requested to be skipped.";
2111  continue;
2112  }
2113 
2114  //
2115  // assign the offline channel ID
2116  //
2117  raw::Channel_t channel = thisChannelSetup.hasChannel()
2118  ? thisChannelSetup.channelID: channelNumberToChannel(channelNumber);
2119 
2120  if (!thisChannelSetup.isChannel(channel)) {
2121  if (thisChannelSetup.mustSave()) {
2122  /*
2123  * Assuming that there are no errors in the database and in the logic
2124  * of this function, this is a "special" channel not connected to a PMT,
2125  * and the configuration of the module is explicitly asking for it to
2126  * be saved but it does not specify a channel ID for it
2127  * (the module configuration is the only source of ID for such channels)
2128  * Together with the request for saving the channel is,
2129  * the configuration should also provide the ID to save it with.
2130  */
2131  throw cet::exception("DaqDecoderICARUSPMT")
2132  << "Channel number " << channelNumber << " of board 0x"
2133  << std::hex << fragInfo.fragmentID << std::dec
2134  << " is not associated to any channel ID"
2135  << " but was demanded to be saved.\n"
2136  ;
2137  }
2138  mf::LogTrace(fLogCategory)
2139  << "Channel number " << channelNumber << " of board 0x"
2140  << std::hex << fragInfo.fragmentID << std::dec
2141  << " is enabled but not associated to any channel ID:"
2142  " it will be skipped.";
2143  continue;
2144  } // if no channel ID
2145  assert(thisChannelSetup.isChannel(channel));
2146 
2147  //
2148  // global trigger
2149  //
2150 
2151  // we never skip waveforms not on global trigger here, because after merging
2152  // they may become part of a waveform that is on it; same about minimum span
2153 
2154 
2155  //
2156  // fill the waveform data
2157  //
2158  std::size_t const ch_offset = iChunk * fragInfo.nSamplesPerChannel;
2159  std::copy_n(fragInfo.data + ch_offset, wvfm.size(), wvfm.begin());
2160 
2161  //
2162  // create the proto-waveform
2163  //
2164  auto const [ itMin, itMax ] = std::minmax_element(wvfm.begin(), wvfm.end());
2165  protoWaveforms.push_back({ // create the waveform and its ancillary info
2166  raw::OpDetWaveform{ timeStamp.value(), channel, wvfm } // waveform
2167  , &thisChannelSetup // channelSetup
2168  , onGlobal // onGlobal
2169  , *itMin // minSample
2170  , *itMax // maxSample
2171  });
2172 
2173  mf::LogTrace log(fLogCategory);
2174  log << "PMT channel " << dumpChannel(protoWaveforms.back())
2175  << " has " << wvfm.size() << " samples (read from entry #" << iChunk
2176  << " in fragment data) starting at electronics time " << timeStamp;
2177  if (protoWaveforms.back().onGlobal) log << ", on global trigger";
2178  if (!protoWaveforms.back().channelSetup->category.empty()) {
2179  log << "; category: '" << protoWaveforms.back().channelSetup->category
2180  << "'";
2181  }
2182 
2183  } // for all channels in the board
2184 
2185 
2186  if (diagOut) {
2187  (*diagOut)
2188  << " - number of waveforms decoded: " << protoWaveforms.size();
2189  if (!protoWaveforms.empty() && protoWaveforms.back().onGlobal)
2190  (*diagOut) << " - matches global trigger!";
2191  } // if diagnostics
2192 
2193  return protoWaveforms;
2194 
2195 } // icarus::DaqDecoderICARUSPMT::createFragmentWaveforms()
2196 
2197 
2198 //------------------------------------------------------------------------------
2200  (std::vector<ProtoWaveform_t>& waveforms) const
2201 {
2202  std::size_t const nWaveforms = waveforms.size();
2203  if (nWaveforms < 2) return 0U;
2204 
2205  // matching criteria: allow rounding error (10 ps)
2206  // and require the second time closer than half a tick from the first one
2207  auto const matchTimes = [this](electronics_time a, electronics_time b)
2208  { return ((a - 0.01_ns) < b) && (b < a + fOpticalTick / 2.0); };
2209 
2210  sortWaveforms(waveforms);
2211 
2212  // define the merging groups
2213  std::vector<std::vector<std::size_t>> waveformGroups;
2214  std::size_t iWave = 0U;
2215  do {
2216  // start a new group with the next available waveform
2217  std::vector<std::size_t> group{ iWave }; // one element: `iWave`
2218 
2219  // merge all following waveforms contiguous to this group
2220  electronics_time currentEnd = waveformEndTime(waveforms[iWave]);
2221  raw::Channel_t const currentChannel
2222  = waveforms[iWave].waveform.ChannelNumber();
2223  while (++iWave < nWaveforms) {
2224  raw::OpDetWaveform const& waveform = waveforms[iWave].waveform;
2225  if (waveform.ChannelNumber() != currentChannel) break;
2226  if (!matchTimes(currentEnd, waveformStartTime(waveform))) break;
2227  group.push_back(iWave);
2228  currentEnd = waveformEndTime(waveform);
2229  } // while matching times
2230 
2231  waveformGroups.push_back(std::move(group));
2232  } while (iWave < nWaveforms);
2233 
2234  std::vector<ProtoWaveform_t> mergedWaveforms;
2235  mergedWaveforms.reserve(waveformGroups.size());
2236  for (auto const& group: waveformGroups)
2237  mergedWaveforms.push_back(mergeWaveformGroup(waveforms, group));
2238  waveforms = std::move(mergedWaveforms);
2239  return nWaveforms - waveforms.size();
2240 } // icarus::DaqDecoderICARUSPMT::mergeWaveforms()
2241 
2242 
2243 //------------------------------------------------------------------------------
2245  std::vector<ProtoWaveform_t>& allWaveforms,
2246  std::vector<std::size_t> const& indices
2247 ) const -> ProtoWaveform_t {
2248 
2249  if (indices.empty()) {
2250  throw std::logic_error
2251  { "DaqDecoderICARUSPMT::mergeWaveformGroup(): empty waveform group." };
2252  }
2253 
2254  auto itIndex = indices.begin();
2255  auto const iend = indices.end();
2256  ProtoWaveform_t mergedWaveform{ std::move(allWaveforms.at(*itIndex)) };
2257  /*
2258  mf::LogTrace("DaqDecoderICARUSPMT")
2259  << "Extending waveform [#" << (*itIndex) << "] channel="
2260  << dumpChannel(mergedWaveform) << " time="
2261  << waveformStartTime(mergedWaveform) << " -- "
2262  << waveformEndTime(mergedWaveform) << " (" << mergedWaveform.waveform.size()
2263  << " samples)"
2264  ;
2265  */
2266  while (++itIndex != iend) {
2267  ProtoWaveform_t& wf = allWaveforms.at(*itIndex);
2268  mf::LogTrace("DaqDecoderICARUSPMT")
2269  << " - extending waveform channel=" << dumpChannel(mergedWaveform)
2270  << " time=" << waveformStartTime(mergedWaveform)
2271  << " -- " << waveformEndTime(mergedWaveform)
2272  << " (" << mergedWaveform.waveform.size()
2273  << " samples) with waveform [#" << (*itIndex) << "] channel="
2274  << dumpChannel(wf) << " at time=" << waveformStartTime(wf.waveform)
2275  << " (" << wf.waveform.size() << " samples)"
2276  ;
2277  if (mergedWaveform.waveform.ChannelNumber() != wf.waveform.ChannelNumber())
2278  {
2279  throw std::logic_error{
2280  "DaqDecoderICARUSPMT::mergeWaveformGroup(): "
2281  "attempt to merge waveforms from channels "
2282  + std::to_string(mergedWaveform.waveform.ChannelNumber())
2283  + " and " + std::to_string(wf.waveform.ChannelNumber())
2284  };
2285  }
2286  if (wf.waveform.empty()) {
2287  throw std::logic_error{
2288  "DaqDecoderICARUSPMT::mergeWaveformGroup(): "
2289  "attempt to merge a waveform (channel "
2290  + std::to_string(mergedWaveform.waveform.ChannelNumber())
2291  + ", timestamp " + std::to_string(wf.waveform.TimeStamp()) + " again"
2292  };
2293  }
2294  if (wf.channelSetup->category != mergedWaveform.channelSetup->category) {
2295  throw std::logic_error{
2296  "DaqDecoderICARUSPMT::mergeWaveformGroup(): "
2297  "attempt to merge a waveform (channel "
2298  + std::to_string(mergedWaveform.waveform.ChannelNumber())
2299  + " of category '" + mergedWaveform.channelSetup->category
2300  + "' with a waveform (channel "
2301  + std::to_string(wf.waveform.ChannelNumber())
2302  + " of the different category '" + wf.channelSetup->category + "'"
2303  };
2304  }
2305  if (wf.channelSetup != mergedWaveform.channelSetup) {
2306  throw std::logic_error{
2307  "DaqDecoderICARUSPMT::mergeWaveformGroup(): "
2308  "attempt to merge a waveform (channel "
2309  + std::to_string(mergedWaveform.waveform.ChannelNumber())
2310  + " with channel settings different than the other waveform (channel "
2311  + std::to_string(wf.waveform.ChannelNumber())
2312  };
2313  }
2314  // raw::OpDetWaveform happen to be `std::vector`, so this will work:
2315  std::size_t const expectedSize [[maybe_unused]]
2316  = mergedWaveform.waveform.size() + wf.waveform.size();
2317  appendTo(mergedWaveform.waveform, std::move(wf.waveform));
2318  mergedWaveform.onGlobal |= wf.onGlobal;
2319  if (wf.minSample < mergedWaveform.minSample)
2320  mergedWaveform.minSample = wf.minSample;
2321  if (wf.maxSample > mergedWaveform.maxSample)
2322  mergedWaveform.maxSample = wf.maxSample;
2323  assert(wf.waveform.empty());
2324  assert(mergedWaveform.waveform.size() == expectedSize);
2325  } // while
2326  return mergedWaveform;
2327 } // icarus::DaqDecoderICARUSPMT::mergeWaveformGroup()
2328 
2329 
2330 //------------------------------------------------------------------------------
2332  FragmentInfo_t const& fragInfo,
2333  TriggerInfo_t const& triggerInfo,
2334  electronics_time waveformTimestamp
2335 ) {
2336 
2337  if (!fTreeFragment) return;
2338 
2339  fTreeFragment->data.fragmentID = fragInfo.fragmentID;
2340  fTreeFragment->data.fragCount = fragInfo.eventCounter;
2341  fTreeFragment->data.TriggerTimeTag = fragInfo.TTT;
2342  fTreeFragment->data.trigger = triggerInfo.time;
2343  fTreeFragment->data.relBeamGate = triggerInfo.relBeamGateTime;
2344  fTreeFragment->data.fragTime
2345  = { static_cast<long long int>(fragInfo.fragmentTimestamp) };
2346  fTreeFragment->data.waveformTime = waveformTimestamp.value();
2347  fTreeFragment->data.waveformSize = fragInfo.nSamplesPerChannel;
2348  fTreeFragment->data.triggerBits = triggerInfo.bits;
2349  fTreeFragment->data.gateCount = triggerInfo.gateCount;
2350  fTreeFragment->data.onGlobalTrigger
2351  = containsGlobalTrigger(waveformTimestamp, fragInfo.nSamplesPerChannel);
2352  assignEventInfo(fTreeFragment->data);
2353  fTreeFragment->tree->Fill();
2354 
2355 } // icarus::DaqDecoderICARUSPMT::fillPMTfragmentTree()
2356 
2357 
2358 //------------------------------------------------------------------------------
2360  (artdaq::Fragment const& fragment) -> BoardID_t
2361 {
2362  return static_cast<BoardID_t>(fragment.fragmentID());
2363 } // icarus::DaqDecoderICARUSPMT::extractFragmentBoardID()
2364 
2365 
2366 //------------------------------------------------------------------------------
2368  (artdaq::Fragment const& artdaqFragment) const -> FragmentInfo_t
2369 {
2370  //
2371  // fragment ID, timestamp and data begin
2372  //
2373  artdaq::Fragment::fragment_id_t const fragment_id
2374  = artdaqFragment.fragmentID();
2375  artdaq::Fragment::timestamp_t const fragmentTimestamp
2376  = artdaqFragment.timestamp();
2377  std::uint16_t const* data_begin = reinterpret_cast<std::uint16_t const*>
2378  (artdaqFragment.dataBeginBytes() + sizeof(sbndaq::CAENV1730EventHeader));
2379 
2380  //
2381  // event counter, trigger time tag, enabled channels
2382  //
2383  sbndaq::CAENV1730Fragment const fragment { artdaqFragment };
2384  sbndaq::CAENV1730FragmentMetadata const& metafrag = *(fragment.Metadata());
2385  sbndaq::CAENV1730EventHeader const& header = fragment.Event()->Header;
2386 
2387  unsigned int const eventCounter = header.eventCounter;
2388 
2389  unsigned int const TTT = header.triggerTimeTag;
2390 
2391  std::uint16_t const enabledChannels = header.ChannelMask();
2392 
2393  //
2394  // samples per channel
2395  //
2396  // artDAQ size is measured in 4-byte words
2397  std::size_t const eventSize = header.eventSize * sizeof(std::uint32_t);
2398 
2399  constexpr std::size_t headerSize = sizeof(sbndaq::CAENV1730EventHeader);
2400 
2401  std::size_t const sampleDataSize = eventSize - headerSize;
2402 
2403  std::size_t const samplesInFragment = sampleDataSize / sizeof(std::uint16_t);
2404 
2405  std::size_t const nChannelsPerBoard = metafrag.nChannels;
2406  std::size_t const nSamplesPerChannel = samplesInFragment / nChannelsPerBoard;
2407 
2408  //
2409  // diagnostics dump
2410  //
2411  if (fDiagnosticOutput) {
2412 
2413  mf::LogVerbatim{ fLogCategory }
2414  << "----> PMT Fragment ID: " << std::hex << fragment_id << std::dec
2415  << ", size: " << eventSize << " B"
2416  << ", data size: " << samplesInFragment << " samples"
2417  << "\n "
2418  << " channels/board: " << nChannelsPerBoard
2419  << ", enabled: " << icarus::ns::util::bin(enabledChannels)
2420  << ", samples/channel: " << nSamplesPerChannel
2421  << "\n "
2422  << " event counter: " << eventCounter
2423  << ", trigger time tag: " << TTT
2424  << ", time stamp: " << (fragmentTimestamp / 1'000'000'000UL)
2425  << "." << (fragmentTimestamp % 1'000'000'000UL) << " s"
2426  ;
2427 
2428  } // if diagnostics
2429 
2430  //
2431  // all done
2432  //
2433  return { // C++20: write the member names explicitly
2434  fragment_id,
2435  fragmentTimestamp,
2436  eventCounter,
2437  TTT,
2438  enabledChannels,
2439  nSamplesPerChannel,
2440  data_begin
2441  };
2442 
2443 } // icarus::DaqDecoderICARUSPMT::extractFragmentInfo()
2444 
2445 
2446 //------------------------------------------------------------------------------
2448  FragmentInfo_t const& fragInfo,
2449  NeededBoardInfo_t const& boardInfo,
2450  SplitTimestamp_t triggerTime
2451 ) const {
2452 
2453  // check availability of mapping for this board, otherwise can't do anything
2454  std::size_t const fragment_id = fragInfo.fragmentID;
2455  std::size_t const eff_fragment_id = effectivePMTboardFragmentID(fragment_id);
2456 
2457  if (!fChannelMap.hasPMTDigitizerID(eff_fragment_id)) {
2458  mf::LogError(fLogCategory)
2459  << "*** PMT could not find channel information for fragment: "
2460  << fragment_id
2461  ;
2462  return NoTimestamp;
2463  }
2464 
2465  if (fTTTresetEverySecond)
2466  return fragmentWaveformTimestampFromTTT(fragInfo, boardInfo, triggerTime);
2467  else
2468  return fragmentWaveformTimestampOnTrigger(fragInfo, boardInfo, triggerTime);
2469 
2470 } // icarus::DaqDecoderICARUSPMT::fragmentWaveformTimestamp()
2471 
2472 
2473 //------------------------------------------------------------------------------
2475  FragmentInfo_t const& /* fragInfo */,
2476  NeededBoardInfo_t const& boardInfo,
2477  SplitTimestamp_t /* triggerTime */
2478 ) const -> electronics_time {
2479 
2480  nanoseconds const preTriggerTime = boardInfo.preTriggerTime;
2481  nanoseconds const PMTtriggerDelay = boardInfo.PMTtriggerDelay;
2482 
2483  auto const timestamp = fNominalTriggerTime - PMTtriggerDelay - preTriggerTime;
2484  mf::LogTrace(fLogCategory) << "V1730 board '" << boardInfo.name
2485  << "' has data starting at electronics time " << timestamp
2486  << " = " << fNominalTriggerTime
2487  << " - " << PMTtriggerDelay << " - " << preTriggerTime
2488  ;
2489  return timestamp;
2490 
2491 } // icarus::DaqDecoderICARUSPMT::fragmentWaveformTimestampOnTrigger()
2492 
2493 
2494 //------------------------------------------------------------------------------
2496  FragmentInfo_t const& fragInfo,
2497  NeededBoardInfo_t const& boardInfo,
2498  SplitTimestamp_t triggerTime
2499 ) const -> electronics_time {
2500 
2501  /*
2502  * 1. goal is a timestamp in electronics time
2503  * 2. we have the global trigger time in electronics time
2504  * (from raw::Trigger data product or from DetectorClocks service)
2505  * 3. board TTT and global trigger time are on the same timescale:
2506  * their difference directly describes the time of the board trigger
2507  * relative to the global trigger
2508  * 4. TTT tags the last (or after the last?) sample of the collected waveform;
2509  * the time of the first sample precedes that tag by the full buffer length
2510  * 5. the PMT trigger itself is subject to a sequence of delays compared to
2511  * the (local or global) trigger from SPEXi; here we quantify these delays
2512  * from calibration offsets collectively passed via job configuration.
2513  */
2514 
2515  using namespace util::quantities::time_literals;
2516 
2517  //
2518  // 2. global trigger time
2519  //
2520  electronics_time waveformTime = fNominalTriggerTime;
2521 
2522  //
2523  // 3. PMT readout board trigger relative to global trigger time
2524  //
2525  unsigned int const triggerTimeNS = triggerTime.split.nanoseconds;
2526 
2527  // converted into nanoseconds (each TTT tick is 8 nanoseconds):
2528  unsigned int const TTT = fragInfo.TTT * 8;
2529 
2530  /*
2531  * The trigger time tag (TTT) on the PMT readout board is incremented every 8
2532  * nanoseconds, and the board is sent a reset signal every second, matching
2533  * the time of the change of second of the global trigger time scale
2534  * (International Atomic Time from White Rabbit). If the global trigger and
2535  * the trigger of the PMT readout happen at the same instant (which would be
2536  * ideally true for one fragment for each board and each event), the TTT will
2537  * represent exactly the number of nanoseconds of the global trigger passed
2538  * since the last crossing of a second boundary.
2539  * (in practice, this is biassed by the signal creation, propagation and
2540  * interpretation delays, and smeared by the clocks of the SPEXi board sending
2541  * the signal and the CAEN V1730 readout board, which have a period of 25 ns
2542  * and 8 ns, respectively).
2543  *
2544  * Multiple fragments can be collected from one PMT readout board for the same
2545  * event by sending the board multiple "local" triggers; these triggers are
2546  * all supposed to be in the neighbourhood of the global trigger time;
2547  * in fact, they should not be more than 1 ms away from it, since PMT readout
2548  * enable window is +/- 1 ms around the expected beam arrival time.
2549  *
2550  * In the most common case when global trigger and fragment tag are not
2551  * separated by a boundary of the second
2552  * (e.g. global trigger at 1620'284'028 seconds + 799'800'000 nanoseconds
2553  * and fragment tag 0.5 milliseconds earlier, at 1620'284'028 seconds
2554  * + 799'300'000 nanoseconds, i.e. with a trigger tag around 799'300'000);
2555  * in this example, the fragment is 0.5 ms earlier:
2556  * -500'000 nanoseconds = 799'300'000 (TTT) - 799'800'000 (glob. trigger)
2557  */
2558  int fragmentRelTime = static_cast<int>(TTT) - triggerTimeNS;
2559  if (TTT > triggerTimeNS + 500'000'000U) { // 0.5 seconds
2560  /*
2561  * case when global trigger arrives just after the boundary of the second
2562  * (e.g. global trigger at 1620'284'029 seconds + 300'000 nanoseconds)
2563  * and this fragment was tagged before that crossing (e.g. 0.5 milliseconds
2564  * earlier, at 1620'284'028 seconds + 999'800'000 nanoseconds, i.e. with a
2565  * trigger tag around 999'800'000);
2566  * in this example, the fragment is 0.5 ms earlier: -500'000 nanoseconds =
2567  * 999'800'000 (TTT) - 1'000'000'000 (second step) - 300'000 (glob. trigger)
2568  * and the plain difference,
2569  * 999'800'000 (TTT) - 300'000 (glob. trigger) = 999'500'000,
2570  * must be corrected by removing a whole second:
2571  */
2572  fragmentRelTime -= 1'000'000'000;
2573  }
2574  else if (TTT + 500'000'000U < triggerTimeNS) { // 0.5 seconds
2575  /*
2576  * case when global trigger arrives just before the boundary of the second
2577  * (e.g. global trigger at 1620'284'028 seconds + 999'800'000 nanoseconds)
2578  * and this fragment was tagged after that crossing (e.g. 0.5 milliseconds
2579  * later, at 1620'284'029 seconds + 300'000 nanoseconds, i.e. with a
2580  * trigger tag around 300'000 because of the pulse-per-second reset);
2581  * in this example, the fragment is 0.5 ms late: 500'000 nanoseconds =
2582  * 300'000 (TTT) + 1'000'000'000 (second step) - 999'800'000 (glob. trigger)
2583  * and the plain difference,
2584  * 300'000 (TTT) - 999'800'000 (glob. trigger) = -999'500'000,
2585  * must be corrected by adding a whole second:
2586  */
2587  fragmentRelTime += 1'000'000'000;
2588  }
2589 
2590  waveformTime += nanoseconds::castFrom(fragmentRelTime);
2591 
2592  //
2593  // 4. correction for relative position of sample tagged by TTT in the buffer
2594  //
2595  waveformTime -= fragInfo.nSamplesPerChannel * fOpticalTick;
2596 
2597  //
2598  // 5. correction for calibrated delays
2599  //
2600  /*
2601  * Waveform time has been expressed based on the "absolute" trigger time plus
2602  * an offset based on the Trigger Time Tag, which is synchronous with the
2603  * global trigger and reset every second.
2604  * We are missing a possible delay between the time of the trigger time scale
2605  * stepping into a new second and the the time TTT reset is effective.
2606  *
2607  *
2608  */
2609 
2610  waveformTime += boardInfo.TTTresetDelay;
2611 
2612  mf::LogTrace(fLogCategory) << "V1730 board '" << boardInfo.name
2613  << "' has data starting at electronics time " << waveformTime
2614  << " = " << fNominalTriggerTime << " (global trigger)"
2615  << " + " << nanoseconds(fragmentRelTime) << " (TTT - global trigger)"
2616  << " - " << (fragInfo.nSamplesPerChannel * fOpticalTick) << " (buffer size)"
2617  << " + " << boardInfo.TTTresetDelay << " (reset delay)"
2618  ;
2619 
2620  return waveformTime;
2621 
2622 
2623 } // icarus::DaqDecoderICARUSPMT::fragmentWaveformTimestampFromTTT()
2624 
2625 
2626 //------------------------------------------------------------------------------
2628  (artdaq::Fragment::fragment_id_t fragment_id) const -> NeededBoardInfo_t
2629 {
2630 
2631  assert(fBoardInfoLookup);
2632 
2633  /*
2634  * The trigger time is always the nominal one, because that is the
2635  * reference time of the whole DAQ (PMT, TPC...).
2636  * We only need to know how sooner than the trigger the V1730 buffer
2637  * starts. Oh, and the delay from the global trigger time to when
2638  * the readout board receives and processes the trigger signal.
2639  */
2641  = fBoardInfoLookup->findBoardInfo(fragment_id);
2642  if (!boardInfo) {
2643  if (fRequireKnownBoards) {
2644  cet::exception e("DaqDecoderICARUSPMT");
2645  e << "Input fragment has ID " << fragment_id
2646  << " which has no associated board information (`BoardSetup`";
2647  if (!hasPMTconfiguration()) e << " + `.FragmentID`";
2648  throw e << ").\n";
2649  }
2650  }
2651  else {
2652  assert(boardInfo->fragmentID == fragment_id);
2653  assert(boardInfo->setup);
2654  }
2655 
2656  return fetchNeededBoardInfo(boardInfo, fragment_id);
2657 } // icarus::DaqDecoderICARUSPMT::neededBoardInfo()
2658 
2659 
2660 //------------------------------------------------------------------------------
2661 std::set<std::string> icarus::DaqDecoderICARUSPMT::getAllInstanceNames() const {
2662  std::set<std::string> names;
2663  for (daq::details::BoardSetup_t const& setup: fBoardSetup) {
2664  for (daq::details::BoardSetup_t::ChannelSetup_t const& chSetup
2665  : setup.channelSettings
2666  ) {
2667  if (chSetup.mustSkip()) continue;
2668  names.insert(chSetup.category);
2669  } // for all channels
2670  } // for all boards
2671  return names;
2672 } // icarus::DaqDecoderICARUSPMT::getAllInstanceNames()
2673 
2674 
2675 //------------------------------------------------------------------------------
2677  (std::vector<daq::details::BoardSetup_t> const& allBoardSetup) const
2678 {
2679  //
2680  // duplicate special channel ID check
2681  //
2682  std::unordered_map<raw::Channel_t, unsigned int> assignedChannels;
2683  bool hasDuplicates = false;
2684  for (daq::details::BoardSetup_t const& setup: allBoardSetup) {
2685  for (daq::details::BoardSetup_t::ChannelSetup_t const& chSetup
2686  : setup.channelSettings
2687  ) {
2688  if (chSetup.mustSkip()) continue; // skipped channels do not matter anyway
2689  // special channels by definition have a non-empty category
2690  if (!chSetup.hasChannel() || chSetup.category.empty()) continue;
2691  if (++assignedChannels[chSetup.channelID] > 1) hasDuplicates = true;
2692  } // for all channels
2693  } // for all boards
2694  if (!hasDuplicates) return;
2695 
2696  art::Exception e { art::errors::Configuration };
2697  e << "Some channel ID are specified in multiple board special channel setup:";
2698  for (auto [ channelID, entries ]: assignedChannels) {
2699  if (entries <= 1) continue; // not a duplicate
2700  e << "\n - channel ID 0x" << std::hex << channelID << std::dec << ": "
2701  << entries << " times";
2702  } // for
2703 
2704  throw e;
2705 } // icarus::DaqDecoderICARUSPMT::checkBoardSetup()
2706 
2707 
2708 //------------------------------------------------------------------------------
2710  (artdaq::Fragment const& fragment)
2711 {
2712  sbndaq::CAENV1730Fragment const V1730fragment { fragment };
2713  sbndaq::CAENV1730EventHeader const header = V1730fragment.Event()->Header;
2714 
2715  unsigned int TTT { header.triggerTimeTag }; // prevent narrowing
2716  return TTT;
2717 
2718 } // icarus::DaqDecoderICARUSPMT::extractTriggerTimeTag()
2719 
2720 
2721 //------------------------------------------------------------------------------
2723  (std::vector<ProtoWaveform_t>& waveforms) const
2724 {
2725  auto byChannelThenTime = []
2726  (ProtoWaveform_t const& left, ProtoWaveform_t const& right)
2727  {
2728  return (left.waveform.ChannelNumber() != right.waveform.ChannelNumber())
2729  ? left.waveform.ChannelNumber() < right.waveform.ChannelNumber()
2730  : left.waveform.TimeStamp() < right.waveform.TimeStamp();
2731  };
2732 
2733  std::sort(waveforms.begin(), waveforms.end(), byChannelThenTime);
2734 
2735 } // icarus::DaqDecoderICARUSPMT::sortWaveforms()
2736 
2737 
2738 //------------------------------------------------------------------------------
2740  (std::vector<std::string> const& treeNames)
2741 {
2742 
2743  auto findTree = [](std::string const& name)
2744  {
2745  return static_cast<DataTrees>(
2746  std::distance(TreeNames.begin(),
2747  std::find(TreeNames.begin(), TreeNames.end(), name))
2748  );
2749  };
2750 
2751  for (std::string const& name: treeNames) {
2752  switch (findTree(name)) {
2753  case DataTrees::Fragments: initFragmentsTree(); break;
2754  case DataTrees::N:
2755  default:
2756  throw cet::exception("DaqDecoderICARUSPMT")
2757  << "initTrees(): no data tree supported with name '" << name
2758  << "'.\n";
2759  } // switch
2760  } // for names
2761 
2762 } // icarus::DaqDecoderICARUSPMT::initTrees()
2763 
2764 
2765 //------------------------------------------------------------------------------
2767  (TTree& tree, TreeData_EventID_t& data)
2768 {
2769 
2770  usesEventInfo(); // this tree includes event information
2771 
2772  tree.Branch("run", &data.run);
2773  tree.Branch("subrun", &data.subrun);
2774  tree.Branch("event", &data.event);
2775  tree.Branch("timestamp", &data.timestamp.time, "timestamp/L");
2776 
2777 } // icarus::DaqDecoderICARUSPMT::initEventIDtree()
2778 
2779 
2780 //------------------------------------------------------------------------------
2782 
2783  if (fTreeFragment) return;
2784 
2785  TTree* tree = art::ServiceHandle<art::TFileService>()
2786  ->make<TTree>("PMTfragments", "PMT fragment data");
2787 
2788  fTreeFragment = std::make_unique<TreeFragment_t>();
2789  fTreeFragment->tree = tree;
2790  auto& data = fTreeFragment->data;
2791 
2792  initEventIDtree(*tree, data);
2793 
2794  tree->Branch("fragmentID", &data.fragmentID);
2795  tree->Branch("fragCount", &data.fragCount);
2796  tree->Branch("fragTime", &data.fragTime.time, "fragTime/L"); // ROOT 6.24 can't detect 64-bit
2797  tree->Branch("fragTimeSec", &data.fragTime.split.seconds);
2798  tree->Branch("TTT", &data.TriggerTimeTag, "TTT/l"); // ROOT 6.24 can't detect 64-bit
2799  tree->Branch("trigger", &data.trigger.time, "trigger/L"); // ROOT 6.24 can't detect 64-bit
2800  tree->Branch("triggerSec", &data.trigger.split.seconds);
2801  tree->Branch("triggerNS", &data.trigger.split.nanoseconds);
2802  tree->Branch("relBeamGateNS", &data.relBeamGate, "relBeamGateNS/I"); // ROOT 6.24 can't handle `long` neither
2803  tree->Branch("waveformTime", &data.waveformTime);
2804  tree->Branch("waveformSize", &data.waveformSize);
2805  tree->Branch("triggerBits", &data.triggerBits);
2806  tree->Branch("gateCount", &data.gateCount);
2807  tree->Branch("onGlobal", &data.onGlobalTrigger);
2808 
2809 } // icarus::DaqDecoderICARUSPMT::initFragmentsTree()
2810 
2811 
2812 //------------------------------------------------------------------------------
2814 
2815  // the allocation of fEventInfo is the flag for event information usage
2816  if (!fEventInfo) fEventInfo = std::make_unique<TreeData_EventID_t>();
2817 
2818 } // icarus::DaqDecoderICARUSPMT::usesEventInfo()
2819 
2820 
2822  (TreeData_EventID_t& treeData) const
2823 {
2824 
2825  assert(fEventInfo);
2826  treeData = *fEventInfo; // nice slicing
2827 
2828 } // icarus::DaqDecoderICARUSPMT::assignEventInfo()
2829 
2830 
2831 //------------------------------------------------------------------------------
2833  (art::Event const& event, TreeData_EventID_t& treeData) const
2834 {
2835  art::EventID const& id = event.id();
2836  treeData.run = id.run();
2837  treeData.subrun = id.subRun();
2838  treeData.event = id.event();
2839 
2840  art::Timestamp const& timestamp = event.time();
2841  treeData.timestamp
2842  = { static_cast<int>(timestamp.timeHigh()), timestamp.timeLow() };
2843 
2844 } // icarus::DaqDecoderICARUSPMT::fillTreeEventID()
2845 
2846 
2847 //------------------------------------------------------------------------------
2848 DEFINE_ART_MODULE(icarus::DaqDecoderICARUSPMT)
2849 
2850 
2851 //------------------------------------------------------------------------------
NeededBoardInfo_t fetchNeededBoardInfo(daq::details::BoardInfoLookup::BoardInfo_t const *boardInfo, unsigned int fragmentID) const
Puts together all the needed information for a board.
Definitions of the trigger bits for SBN.
detinfo::DetectorTimings const fDetTimings
Interface to LArSoft configuration for detector timing.
void beginRun(art::Run &run) override
On a new run: cache PMT configuration information.
Information used in decoding from a board.
electronics_time waveformStartTime(ProtoWaveform_t const &wf) const
std::vector< daq::details::BoardSetup_t > const fBoardSetup
All board setup settings.
std::string boardName
Name (mnemonic) of the board.
void usesEventInfo()
Declares the use of event information.
void endJob() override
Prints a end-of-job message.
then source grid fermiapp products dune setup_dune_fermiapp sh exit else echo No setup file found exit fi setup
std::array< std::string, static_cast< std::size_t >(DataTrees::N)> TreeNameList_t
void fillPMTfragmentTree(FragmentInfo_t const &fragInfo, TriggerInfo_t const &triggerInfo, electronics_time waveformTimestamp)
std::vector< DigitizerChannelChannelIDPair > DigitizerChannelChannelIDPairVec
bool const fSurviveExceptions
Whether to &quot;ignore&quot; errors.
NeededBoardInfo_t neededBoardInfo(artdaq::Fragment::fragment_id_t fragment_id) const
Returns the board information for this fragment.
void checkBoardSetup(std::vector< daq::details::BoardSetup_t > const &allBoardSetup) const
Throws an exception if the configuration of boards shows errors.
decltype(auto) constexpr cend(T &&obj)
ADL-aware version of std::cend.
Definition: StdUtils.h:87
void produce(art::Event &event) override
Processes the event.
Utility class for fast lookup of board data by fragment ID.
Derivative information of a V1730 readout board.
bool mustSave() const
Whether this channel is requested to be saved.
Data structure for basic event information in simple ROOT trees.
Definition of util::enumerate().
static unsigned int extractTriggerTimeTag(artdaq::Fragment const &fragment)
Extracts the Trigger Time Tag (31+1 bits) value from the fragment.
DataTrees
Enumerate the supported data trees.
microseconds_as<> microseconds
Type of time interval stored in microseconds, in double precision.
Definition: spacetime.h:259
std::array< ChannelSetup_t, NBoardChannels > AllChannelSetup_t
artdaq::FragmentPtrs makeFragmentCollectionFromFragment(artdaq::Fragment const &sourceFragment) const
Converts a plain fragment into a fragment collection.
walls no right
Definition: selectors.fcl:105
std::vector< ProtoWaveform_t > processFragment(artdaq::Fragment const &artdaqFragment, NeededBoardInfo_t const &boardInfo, TriggerInfo_t const &triggerInfo)
Create waveforms and fills trees for the specified artDAQ fragment.
All the information collected about a waveform (with the waveform itself).
bool const fRequireKnownBoards
Whether info on all input boards is required.
#define the
long int relBeamGate
Beam gate start relative to trigger [ns].
bool mustSkip() const
Whether this channel is requested to be skipped.
bool operator<(const VectorMap< _Key, _Tp, _Compare > &__x, const VectorMap< _Key, _Tp, _Compare > &__y)
Definition: VectorMap.h:501
Collection of useful information from fragment data.
electronics_time const fNominalTriggerTime
Trigger time as reported by detinfo::DetectorClocks service.
Some helpers for PMT decoder tool.
void sortWaveforms(std::vector< ProtoWaveform_t > &waveforms) const
Sorts in place the specified waveforms in channel order, then in time.
electronics_time waveformEndTime(raw::OpDetWaveform const &wf) const
TriggerInfo_t fetchTriggerTimestamp(art::Event const &event) const
Retrieves the global trigger time stamp from the event.
void assignEventInfo(TreeData_EventID_t &treeData) const
Assigns the cached event information to the specified tree data.
std::vector< BoardInfo_t > Database_t
pure virtual base interface for detector clocks
std::set< std::string > getAllInstanceNames() const
Returns all the instance names we will produce.
void checkFragmentType(artdaq::Fragment const &artdaqFragment) const
Throws an exception if artdaqFragment is not of type CAEN1730.
ProtoWaveform_t mergeWaveformGroup(std::vector< ProtoWaveform_t > &allWaveforms, std::vector< std::size_t > const &indices) const
util::ArtHandleTrackerManager< art::Event > fDataCacheRemover
Tracks art data products and removes their cached data on demand.
Structure collecting all data for a fragment ROOT tree.
static long long int timestampDiff(std::uint64_t a, std::uint64_t b)
Returns the timestamp difference a - b.
constexpr details::BinObj< T > bin(T value)
Returns a wrapper to print the specified data in binary format.
unsigned int fragmentID
DAQ fragment ID.
Interface to detinfo::DetectorClocks.
process_name gaushit a
artdaq::FragmentPtrs makeFragmentCollectionFromContainerFragment(artdaq::Fragment const &sourceFragment) const
Converts a container fragment into a fragment collection.
unsigned int bufferLength
Ticks in each buffer (recordLength).
fDetProps &fDetProps fDetProps &fDetProps fLogCategory
decltype(auto) constexpr end(T &&obj)
ADL-aware version of std::end.
Definition: StdUtils.h:77
raw::Channel_t channelID
Associated off-line (LArSoft) channel ID.
daq::details::BoardInfoLookup matchBoardConfigurationAndSetup(sbn::PMTconfiguration const *PMTconfig) const
Returns a lookup object with board setup and configuration info.
SplitTimestamp_t time
Time of the trigger (absolute).
electronics_time fragmentWaveformTimestampFromTTT(FragmentInfo_t const &fragInfo, NeededBoardInfo_t const &boardInfo, SplitTimestamp_t triggerTime) const
Returns the timestamp for the waveforms in the specified fragment.
void initEventIDtree(TTree &tree, TreeData_EventID_t &data)
Initializes the event ID part of a tree.
auto counter(T begin, T end)
Returns an object to iterate values from begin to end in a range-for loop.
Definition: counter.h:285
std::vector< ProtoWaveform_t > createFragmentWaveforms(FragmentInfo_t const &fragInfo, AllChannelSetup_t const &channelSetup, electronics_time const timeStamp) const
Creates raw::OpDetWaveform objects from the fragment data.
daq::details::BoardSetup_t::ChannelSetup_t const * channelSetup
Pointer to the settings for the channel this waveform belongs to.
double distance(geo::Point_t const &point, CathodeDesc_t const &cathode)
Returns the distance of a point from the cathode.
Configuration of the V1730 readout board setup.
Additional information on trigger.
second seconds
Alias for common language habits.
Definition: spacetime.h:88
long int relBeamGateTime
Time of beam gate relative to trigger [ns].
std::vector< ProtoWaveform_t const * > findWaveformsWithNominalTrigger(std::vector< ProtoWaveform_t > const &waveforms) const
Returns pointers to all waveforms including the nominal trigger time.
bool const fDropRawDataAfterUse
Clear fragment data product cache after use.
fhicl::OptionalSequence< fhicl::Table< ChannelSetupConfig > > SpecialChannels
electronics_time fragmentWaveformTimestamp(FragmentInfo_t const &fragInfo, NeededBoardInfo_t const &boardInfo, SplitTimestamp_t triggerTime) const
Returns the timestamp for the waveforms in the specified fragment.
daq::details::BoardSetup_t::AllChannelSetup_t AllChannelSetup_t
Type of setup of all channels in a readout board.
FragmentInfo_t extractFragmentInfo(artdaq::Fragment const &artdaqFragment) const
Extracts useful information from fragment data.
constexpr std::array< std::size_t, geo::vect::dimension< Vector >)> indices()
Returns a sequence of indices valid for a vector of the specified type.
std::vector< art::InputTag > const fInputTags
&lt; List of candidate data products with artDAQ data fragments.
process_name PMTconfig
std::optional< daq::details::BoardInfoLookup > fBoardInfoLookup
Find the information on a readout boards by fragment ID.
Special settings for one channel on the board.
std::vector< ProtoWaveform_t > processBoardFragments(artdaq::FragmentPtrs const &artdaqFragment, TriggerInfo_t const &triggerInfo)
Extracts waveforms from the specified fragments from a board.
art::EDProducer::Table< Config > Parameters
return match has_match and(match.match_pdg==11 or match.match_pdg==-11)
bool const fDiagnosticOutput
If true will spew endless messages to output.
BEGIN_PROLOG vertical distance to the surface Name
void initFragmentsTree()
Initializes the fragment data tree (fTreeFragment).
unsigned int mergeWaveforms(std::vector< ProtoWaveform_t > &waveforms) const
walls no left
Definition: selectors.fcl:105
Information from the configuration of PMT readout.
Test of util::counter and support utilities.
std::string const fLogCategory
Message facility category.
nanoseconds const fOpticalTick
Duration of the optical detector readout sampling tick (i.e. 2 ns; hush!).
static constexpr unsigned int NoFragmentID
Special value to mark the absence of fragment ID information.
static bool isChannel(raw::Channel_t channel)
Returns if channel is a valid channel ID.
An interval (duration, length, distance) between two quantity points.
Definition: intervals.h:114
Utility to dump a artDAQ fragment on screen.
Coll::value_type const * binarySearch(Coll const &coll, Key const &key, KeyExtractor &&extractKey)
artdaq::Fragments const & readInputFragments(art::Event const &event) const
Reads the fragments to be processed.
icarusDB::IICARUSChannelMap const & fChannelMap
Fragment/channel mapping database.
std::optional< art::InputTag > const fTriggerTag
Input tag of the global trigger.
electronics_time waveformStartTime(raw::OpDetWaveform const &wf) const
bool containsGlobalTrigger(raw::OpDetWaveform const &waveform) const
Returns whether waveform includes the tick of the nominal trigger time.
Split_t split
Trigger time in nanoseconds from The Epoch (in components).
Record of information about a readout board.
nanoseconds_as<> nanoseconds
Type of time interval stored in nanoseconds, in double precision.
Definition: spacetime.h:270
static const std::vector< std::string > names
Utilities to read interval and point quantity FHiCL configuration.
bool UpdatePMTConfiguration(sbn::PMTconfiguration const *PMTconfig)
Updates the PMT configuration cache. How? Dunno. Placeholder.
electronics_time waveformEndTime(ProtoWaveform_t const &wf) const
raw::OpDetWaveform waveform
The complete waveform.
bool const fRequireBoardConfig
Whether setup info on all boards is required.
Produces raw::OpDetWaveform from V1730 artDAQ data fragments.
static BoardID_t extractFragmentBoardID(artdaq::Fragment const &fragment)
Extracts the fragment ID (i.e. board ID) from the specified fragment.
details::DumpFragWrap dumpFragment(artdaq::Fragment const &frag)
Dump a artDAQ fragment into an output stream.
std::string category
Category of this channel (will become instance name).
std::string to_string(WindowPattern const &pattern)
DaqDecoderICARUSPMT(Parameters const &params)
Constructor.
bool const fTTTresetEverySecond
Whether V1730 TTT is reset every second.
std::ostream & operator<<(std::ostream &, Binner< T > const &)
Definition: Binner.h:218
process_name sequence::icarus_stage0_EastHits_TPC physics sequence::icarus_stage0_EastHits_TPC physics pathNUMI physics streamNUMI outputs outBNB drop *_ *_ *_DAQ drop *_daqTPCROI_ *_ drop *_decon1droi_ *_ drop *_decon1DroiTPC *_ *_ *outputs outNUMI drop *_ *_ *_DAQ drop *_daqTPCROI_ *_ drop *_decon1droi_ *_ drop *_decon1DroiTPC *_ *_ *physics producers daqPMT DecoderTool RequireBoardConfig
Dimensioned variables representing space or time quantities.
void fillTreeEventID(art::Event const &event, TreeData_EventID_t &treeData) const
Fills the base information of a tree data entry from an art event.
dumpChannel(DaqDecoderICARUSPMT::ProtoWaveform_t const &wf)
A class exposing an upgraded interface of detinfo::DetectorClocksData.
process_name largeant stream1 can override from command line with o or output physics producers generator N
do i e
mask_t< triggerSource > triggerSourceMask
Type of mask with triggerSource bits.
Definition: BeamBits.h:109
bool const fPacketDump
Dump V1730 data.
Collection of information about the global trigger.
then echo fcl name
SplitTimestamp_t timestamp
Event timestamp (seconds from the epoch).
Data product holding additional trigger information.
decltype(auto) constexpr begin(T &&obj)
ADL-aware version of std::begin.
Definition: StdUtils.h:72
float postTriggerFrac
Fraction of the waveform after the trigger signal (postPercent).
temporary value
bool hasChannel() const
Returns if this channel is associated to a off-line channel number.
DaqDecoderICARUSPMT::ProtoWaveform_t const & wf
static constexpr auto NoChannelID
Special value for unassigned channel ID.
void initTrees(std::vector< std::string > const &treeNames)
Initializes all requested data trees.
unsigned int TriggerBits() const
Trigger Bits.
Definition: TriggerData.h:54
std::optional< art::InputTag > const fPMTconfigTag
Input tag of the PMT configuration.
std::size_t count(Cont const &cont)
artdaq::FragmentPtrs makeFragmentCollection(artdaq::Fragment const &sourceFragment) const
bool hasPMTconfiguration() const
Returns whether PMT configuration information is expected to be available.
Tracks handles for cache deletion.
Information of the setup of a V1730 readout board.
std::unique_ptr< TreeFragment_t > fTreeFragment
bool empty(FixedBins< T, C > const &) noexcept
Definition: FixedBins.h:555
Class containing configuration for PMT readout.
timescale_traits< ElectronicsTimeCategory >::time_point_t electronics_time
A point in time on the electronics time scale.
electronics_time fragmentWaveformTimestampOnTrigger(FragmentInfo_t const &fragInfo, NeededBoardInfo_t const &boardInfo, SplitTimestamp_t triggerTime) const
Returns the timestamp for the waveforms in the specified fragment.
std::unique_ptr< TreeData_EventID_t > fEventInfo
Event ID for trees.
TimeTrackTreeStorage::TriggerInputSpec_t convert(TimeTrackTreeStorage::Config::TriggerSpecConfig const &config)
Functions to dump the content of binary data chunks to console.
short int BoardID_t
Type used internally to represent board ID.
static std::string name()
Name of this time scale.
std::vector< sbn::V1730Configuration > boards
Configuration of all PMT readout boards.
nanosecond nanoseconds
Alias for common language habits.
Definition: spacetime.h:139
static TreeNameList_t const TreeNames
Class containing configuration for a V1730 board.
bool const fSkipWaveforms
Whether to skip waveform decoding.