All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
opDetSBNDTriggerAlg.hh
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 //// File: opDetSBNDTriggerAlg.hh
3 ////
4 //// This algorithm emulates the behavior of the SBND trigger going
5 //// to the photon detection system. Created by Gray Putnam
6 //// <grayputnam@uchicago.edu>
7 ////
8 //// Changes:
9 ////
10 //// 2021.04.21 (wforeman)
11 //// - Implemented TriggerHoldoff functionality (default behavior unchanged)
12 //// - Fixed bug where PreTrigBeam was getting assigned wrong value
13 //// - Changed name of TriggerEnableWindow parameter to better reflect
14 //// the actual behavior of module
15 //// - Added drift period as fcl parameter, which sets how early the trigger
16 //// enable window starts
17 //// - Beam trigger is now added in FindTriggers instead of in ApplyTriggers
18 //// - Added option to toggle trigger overlap (default behavior unchanged)
19 //// - Modified the ApplyTriggers function to accommodate trigger overlap
20 //// - Restored distinction of beam vs non-beam readout window sizes
21 ////
22 //// 2021.03.24 (wforeman)
23 //// - Fixed wrong threshold used for ARAPUCAs
24 //// - TriggerEnableWindow now starts 1 drift period prior to TPC RO
25 //// - Set readout length to 12us (1us + 11us)
26 //// - Made changes to ApplyTriggers function that enforced no overlap of
27 //// triggers (ie, uniform readout size)
28 ////
29 ////
30 //////////////////////////////////////////////////////////////////////////
31 
32 #ifndef SBND_OPDETSIM_OPDETSBNDTRIGGERALG_HH
33 #define SBND_OPDETSIM_OPDETSBNDTRIGGERALG_HH
34 
35 #include "fhiclcpp/types/Atom.h"
36 #include "fhiclcpp/types/Sequence.h"
37 #include "fhiclcpp/types/OptionalSequence.h"
38 #include "messagefacility/MessageLogger/MessageLogger.h"
39 
40 #include <memory>
41 #include <vector>
42 #include <cmath>
43 #include <string>
44 #include <map>
45 #include <unordered_map>
46 #include <set>
47 #include <sstream>
48 #include <fstream>
49 
53 namespace detinfo {
54  class DetectorClocksData;
55  class DetectorPropertiesData;
56 }
57 
58 #include "TMath.h"
59 #include "TRandom3.h"
60 #include "TF1.h"
61 #include "TH1D.h"
62 
64 
65 namespace opdet {
66 
68 
69  public:
70 
71  struct Config {
72  using Name = fhicl::Name;
73  using Comment = fhicl::Comment;
74 
75  fhicl::Atom<int> PulsePolarityPMT {
76  Name("PulsePolarityPMT"),
77  Comment("Whether pulses go down (-1) or up (1)."),
78  -1
79  };
80 
81  fhicl::Atom<int> PulsePolarityArapuca {
82  Name("PulsePolarityArapuca"),
83  Comment("Whether pulses go down (-1) or up (1)."),
84  1
85  };
86 
87  fhicl::Atom<int> TriggerThresholdADCArapuca {
88  Name("TriggerThresholdADCArapuca"),
89  Comment("Threshold for channel to issue trigger across all arapuca channels. [ADC]")
90  };
91 
92  fhicl::Atom<int> TriggerThresholdADCPMT {
93  Name("TriggerThresholdADCPMT"),
94  Comment("Threshold for channel to issue trigger across all pmt channels. [ADC]")
95  };
96 
97  fhicl::OptionalSequence<unsigned> MaskedChannels {
98  Name("MaskedChannels"),
99  Comment("Channels which are ignored for issuing triggers.")
100  };
101 
102  fhicl::Atom<bool> MaskLightBars {
103  Name("MaskLightBars"),
104  Comment("Whether to mask all light bar readout channels from issuing triggers."),
105  false
106  };
107 
108  fhicl::Atom<bool> MaskPMTs {
109  Name("MaskPMTs"),
110  Comment("Whether to mask all PMT readout channels from issuing triggers."),
111  false
112  };
113 
114  fhicl::Atom<bool> MaskBarePMTs {
115  Name("MaskBarePMTs"),
116  Comment("Whether to mask all Bare PMT readout channels from issuing triggers."),
117  false
118  };
119 
120  fhicl::Atom<bool> MaskXArapucaPrimes {
121  Name("MaskXArapucaPrimes"),
122  Comment("Whether to mask all Arapuca Prime readout channels from issuing triggers."),
123  false
124  };
125 
126  fhicl::Atom<bool> MaskXArapucas {
127  Name("MaskXArapucas"),
128  Comment("Whether to mask all X-Arapuca readout channels from issuing triggers."),
129  false
130  };
131 
132  fhicl::Atom<bool> MaskArapucaT1s {
133  Name("MaskArapucaT1s"),
134  Comment("Whether to mask all Arapuca T1 readout channels from issuing triggers."),
135  false
136  };
137 
138  fhicl::Atom<bool> MaskArapucaT2s {
139  Name("MaskArapucaT1s"),
140  Comment("Whether to mask all Arapuca T2 readout channels from issuing triggers."),
141  false
142  };
143 
144  fhicl::Atom<bool> SelfTriggerPerChannel {
145  Name("SelfTriggerPerChannel"),
146  Comment("If true, each channel's OpDetWaveform will be 'triggered' when that channel goes above threshold (this setting ignores all mask settings). If false, a group of global triggers to be issued across all channel's will instead be collected."),
147  true
148  };
149 
150  fhicl::Atom<double> TriggerHoldoff {
151  Name("TriggerHoldoff"),
152  Comment("When collecting global triggers, sets an amount of time where one trigger will not be issued following the previous one. [us]"),
153  0.
154  };
155 
156  fhicl::Atom<double> TriggerCountWindow {
157  Name("TriggerCountWindow"),
158  Comment("Size of window to count triggers arriving from different readout channels. [us]"),
159  0.
160  };
161 
162  fhicl::Atom<unsigned> TriggerChannelCount {
163  Name("TriggerChannelCount"),
164  Comment("Number of trigger signals from different readout channels required to issue global trigger."),
165  1
166  };
167 
168  fhicl::Atom<bool> BeamTriggerEnable {
169  Name("BeamTriggerEnable"),
170  Comment("Whether to also send a beam trigger."),
171  false
172  };
173 
174  fhicl::Atom<double> BeamTriggerTime {
175  Name("BeamTriggerTime"),
176  Comment("Time at which the beam trigger will be issued. [us]"),
177  0.
178  };
179 
180  fhicl::Atom<double> BeamTriggerHoldoff {
181  Name("BeamTriggerHoldoff"),
182  Comment("Time to ignore other triggers after the beam trigger time. [us]"),
183  0.
184  };
185 
186  fhicl::Atom<bool> AllowTriggerOverlap {
187  Name("AllowTriggerOverlap"),
188  Comment("Extend the readout length if a second trigger occurs during a readout (CAEN board setting)."),
189  true
190  };
191 
193  Name("TriggerEnableWindowOneDriftBeforeTPCReadout"),
194  Comment("If true, triggers will be enabled starting one full drift period prior to the start of the TPC readout, continuing for the length of time that the TPC is read-out (DetectorClocks and DetectorProperties are used to determine this time span. Otherwise, the config options TriggerEnableWindowStart and TriggerEnableWindowLength are checked."),
195  true
196  };
197 
198  fhicl::Atom<double> DriftPeriod {
199  Name("DriftPeriod"),
200  Comment("Drift period, used to determine how soon before the TPC readout the trigger enable window begins. [us]"),
201  true
202  };
203 
204  fhicl::Atom<double> TriggerEnableWindowStart {
205  Name("TriggerEnableWindowStart"),
206  Comment("Start of the window of time for which triggers are enabled. Ignored if TriggerEnableWindowOneDriftBeforeTPCReadout is true. [us]"),
207  0.
208  };
209 
210  fhicl::Atom<double> TriggerEnableWindowLength {
211  Name("TriggerEnableWindowLength"),
212  Comment("Length of time for which trigger are enabled. Ignored if TriggerEnableWindowOneDriftBeforeTPCReadout is true. [us]"),
213  0.
214  };
215 
216  fhicl::Atom<double> ReadoutWindowPreTrigger {
217  Name("ReadoutWindowPreTrigger"),
218  Comment("Size of readout window prior to a non-beam trigger. [us]"),
219  0.
220  };
221 
222  fhicl::Atom<double> ReadoutWindowPostTrigger {
223  Name("ReadoutWindowPostTrigger"),
224  Comment("Size of readout window post a non-beam trigger. [us]"),
225  0.
226  };
227 
228  fhicl::Atom<double> ReadoutWindowPreTriggerBeam {
229  Name("ReadoutWindowPreTriggerBeam"),
230  Comment("Size of readout window pre a beam trigger. [us]"),
231  0.
232  };
233 
234  fhicl::Atom<double> ReadoutWindowPostTriggerBeam {
235  Name("ReadoutWindowPostTriggerBeam"),
236  Comment("Size of readout window post a beam trigger. [us]"),
237  0.
238  };
239  };
240 
241  // construct with config or with fhicl
242  opDetSBNDTriggerAlg(const Config &config);
243 
244  opDetSBNDTriggerAlg(const fhicl::ParameterSet &pset) :
245  opDetSBNDTriggerAlg(fhicl::Table<Config>(pset, {})())
246  {}
247 
248  // Clear out at the end of an event
249  void ClearTriggerLocations();
250 
251  // Add in a waveform to define trigger locations
254  const raw::OpDetWaveform &waveform,
256 
257  // Merge all of the triggers together
258  void MergeTriggerLocations();
259 
260  // Apply trigger locations to an input OpDetWaveform
261  std::vector<raw::OpDetWaveform> ApplyTriggerLocations(detinfo::DetectorClocksData const& clockData, const raw::OpDetWaveform &waveform) const;
262 
263  // Returns the time range over which triggers are enabled over a range [start, end]
264  std::array<double, 2> TriggerEnableWindow(detinfo::DetectorClocksData const& clockData,
266 
267  private:
268 
269  // internal functions
270  bool IsChannelMasked(raw::Channel_t channel) const;
271  bool IsTriggerEnabled(detinfo::DetectorClocksData const& clockData,
274  raw::TimeStamp_t Tick2Timestamp(raw::TimeStamp_t waveform_start, size_t waveform_index) const;
275  const std::vector<raw::TimeStamp_t> &GetTriggerTimes(raw::Channel_t channel) const;
276  double ReadoutWindowPreTrigger(raw::Channel_t channel) const;
277  double ReadoutWindowPostTrigger(raw::Channel_t channel) const;
278  double ReadoutWindowPreTriggerBeam(raw::Channel_t channel) const;
279  double ReadoutWindowPostTriggerBeam(raw::Channel_t channel) const;
280  double OpticalPeriod() const;
281 
282  // fhicl config
284 
285  // OpDet channel map
287 
288  // keeping track of triggers
289  std::map<raw::Channel_t, std::vector<std::array<raw::TimeStamp_t, 2>>> fTriggerRangesPerChannel;
290  std::map<raw::Channel_t, std::vector<raw::TimeStamp_t>> fTriggerLocationsPerChannel;
291  std::vector<raw::TimeStamp_t> fTriggerLocations;
292 
293  std::vector<unsigned> fMaskedChannels;
294 
295  };//class opDetSBNDTriggerAlg
296 
297 
298 } //namespace
299 
300 #endif // SBND_OPDETSIM_OPDETSBNDTRIGGERALG_HH
bool IsTriggerEnabled(detinfo::DetectorClocksData const &clockData, detinfo::DetectorPropertiesData const &detProp, raw::TimeStamp_t trigger_time) const
double ReadoutWindowPostTriggerBeam(raw::Channel_t channel) const
double ReadoutWindowPreTrigger(raw::Channel_t channel) const
fhicl::Atom< double > ReadoutWindowPostTriggerBeam
const std::vector< raw::TimeStamp_t > & GetTriggerTimes(raw::Channel_t channel) const
fhicl::Atom< double > TriggerEnableWindowStart
opDetSBNDTriggerAlg(const fhicl::ParameterSet &pset)
std::array< double, 2 > TriggerEnableWindow(detinfo::DetectorClocksData const &clockData, detinfo::DetectorPropertiesData const &detProp) const
fhicl::Atom< bool > TriggerEnableWindowOneDriftBeforeTPCReadout
std::map< raw::Channel_t, std::vector< std::array< raw::TimeStamp_t, 2 > > > fTriggerRangesPerChannel
void FindTriggerLocations(detinfo::DetectorClocksData const &clockData, detinfo::DetectorPropertiesData const &detProp, const raw::OpDetWaveform &waveform, raw::ADC_Count_t baseline)
double TimeStamp_t
us since 1970, based on TimeService
opDetSBNDTriggerAlg(const Config &config)
fhicl::Atom< unsigned > TriggerChannelCount
fhicl::Atom< double > ReadoutWindowPreTrigger
fhicl::Atom< double > ReadoutWindowPreTriggerBeam
bool IsChannelMasked(raw::Channel_t channel) const
timescale_traits< TriggerTimeCategory >::time_point_t trigger_time
A point in time on the trigger time scale.
double ReadoutWindowPreTriggerBeam(raw::Channel_t channel) const
double OpticalPeriod() const
BEGIN_PROLOG baseline
BEGIN_PROLOG vertical distance to the surface Name
fhicl::OptionalSequence< unsigned > MaskedChannels
std::vector< raw::OpDetWaveform > ApplyTriggerLocations(detinfo::DetectorClocksData const &clockData, const raw::OpDetWaveform &waveform) const
fhicl::Atom< double > ReadoutWindowPostTrigger
std::vector< unsigned > fMaskedChannels
Contains all timing reference information for the detector.
std::vector< raw::TimeStamp_t > fTriggerLocations
fhicl::Atom< double > TriggerEnableWindowLength
double ReadoutWindowPostTrigger(raw::Channel_t channel) const
std::map< raw::Channel_t, std::vector< raw::TimeStamp_t > > fTriggerLocationsPerChannel
raw::TimeStamp_t Tick2Timestamp(raw::TimeStamp_t waveform_start, size_t waveform_index) const
auto const detProp