DUNE-DAQ
DUNE Trigger and Data Acquisition software
Loading...
Searching...
No Matches
TCMakerMichelElectronAlgorithm.hpp
Go to the documentation of this file.
1
9#ifndef TRIGGERALGS_MICHELELECTRON_TRIGGERCANDIDATEMAKERMICHELELECTRON_HPP_
10#define TRIGGERALGS_MICHELELECTRON_TRIGGERCANDIDATEMAKERMICHELELECTRON_HPP_
11
13
14//#include "triggeralgs/triggercandidatemakerhorizontalmuon/Nljs.hpp"
15
16#include <fstream>
17#include <vector>
18
19namespace triggeralgs {
21{
22
23public:
25 void process(const TriggerActivity&, std::vector<TriggerCandidate>&);
26
27 void configure(const nlohmann::json& config);
28
29 // void flush(timestamp_t, std::vector<TriggerCandidate>& output_tc);
30
31private:
32 class Window
33 {
34 public:
35 bool is_empty() const { return inputs.empty(); };
36 void add(const TriggerActivity& input_ta)
37 {
38 // Add the input TA's contribution to the total ADC, increase the hit count
39 // of all of the channels which feature and add it to the TA list keeping the TA
40 // list time ordered by time_start. Preserving time order makes moving easier.
41 adc_integral += input_ta.adc_integral;
42 for (TriggerPrimitive tp : input_ta.inputs) {
44 }
45 // Perform binary search based on time_start.
46 uint16_t insert_at = 0;
47 for (auto ta : inputs) {
48 if (input_ta.time_start < ta.time_start)
49 break;
50 insert_at++;
51 }
52 inputs.insert(inputs.begin() + insert_at, input_ta);
53 };
54 void clear() { inputs.clear(); };
55 uint16_t n_channels_hit() { return channel_states.size(); };
56 void move(TriggerActivity const& input_ta, timestamp_t const& window_length)
57 {
58 // Find all of the TAs in the window that need to be removed
59 // if the input_ta is to be added and the size of the window
60 // is to be conserved.
61 // Subtract those TAs' contribution from the total window ADC and remove their
62 // contributions to the hit counts.
63 uint32_t n_tas_to_erase = 0;
64 for (auto ta : inputs) {
65 if (!(input_ta.time_start - ta.time_start < window_length)) {
66 n_tas_to_erase++;
67 adc_integral -= ta.adc_integral;
68 for (TriggerPrimitive tp : ta.inputs) {
70 // If a TA being removed from the window results in a channel no longer having
71 // any hits, remove from the states map so map.size() can be used for number
72 // channel hits.
73 if (channel_states[tp.channel] == 0)
74 channel_states.erase(tp.channel);
75 }
76 } else
77 break;
78 }
79 // Erase the TAs from the window.
80 inputs.erase(inputs.begin(), inputs.begin() + n_tas_to_erase);
81 // Make the window start time the start time of what is now the
82 // first TA.
83 if (inputs.size() != 0) {
84 time_start = inputs.front().time_start;
85 add(input_ta);
86 } else
87 add(input_ta);
88 }
89 void reset(TriggerActivity const& input_ta)
90 {
91 // Empty the channel and TA lists.
92 channel_states.clear();
93 inputs.clear();
94 // Set the start time of the window to be the start time of the
95 // input_ta.
96 time_start = input_ta.time_start;
97 // Start the total ADC integral.
98 adc_integral = input_ta.adc_integral;
99 // Start hit count for the hit channels.
100 for (TriggerPrimitive tp : input_ta.inputs) {
102 }
103 // Add the input TA to the TA list.
104 inputs.push_back(input_ta);
105 }
106 friend std::ostream& operator<<(std::ostream& os, const Window& window)
107 {
108 if (window.is_empty())
109 os << "Window is empty!\n";
110 // else{
111 // os << "Window start: " << window.time_start << ", end: " << window.inputs.back().time_start;
112 // os << ". Total of: " << window.adc_integral << " ADC counts with " << window.inputs.size() << " TPs.\n";
113 // os << window.channel_states.size() << " independent channels have hits.\n";
114 // }
115 return os;
116 };
117
119 uint64_t adc_integral;
120 std::unordered_map<channel_t, uint16_t> channel_states;
121 std::vector<TriggerActivity> inputs;
122 };
123
125 bool check_adjacency() const;
126
128 uint64_t m_activity_count = 0; // NOLINT(build/unsigned)
129
130 // Configurable parameters.
131 // triggercandidatemakerhorizontalmuon::ConfParams m_conf;
132 // If both m_trigger_on_adc and m_trigger_on_n_channels is false, nothing is done at
133 // the candidate level, candidates are made 1 for 1 with activities.
134 // Use any other combination of m_trigger_on_adc and m_trigger_on_n_channels with caution,
135 // they have not been tested.
136 bool m_trigger_on_adc = false;
138 uint32_t m_adc_threshold = 1200000;
139 uint16_t m_n_channels_threshold = 600; // 80ish for frames, O(200 - 600) for tpslink
141 int tc_number = 0;
142 // Might not be the best type for this map.
143 // std::unordered_map<std::pair<detid_t,channel_t>,channel_t> m_channel_map;
144
145 // For debugging purposes.
146 void add_window_to_record(Window window);
147 void dump_window_record();
148 std::vector<Window> m_window_record;
149};
150} // namespace triggeralgs
151
152#endif // TRIGGERALGS_MICHELELECTRON_TRIGGERCANDIDATEMAKERMICHELELECTRON_HPP_
void move(TriggerActivity const &input_ta, timestamp_t const &window_length)
friend std::ostream & operator<<(std::ostream &os, const Window &window)
void process(const TriggerActivity &, std::vector< TriggerCandidate > &)
The function that gets call when there is a new activity.
dunedaq::trgdataformats::timestamp_t timestamp_t
Definition Types.hpp:16
A single energy deposition on a TPC or PDS channel.
std::vector< TriggerPrimitive > inputs